Loading...
Searching...
No Matches
patchCloudSet.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "patchCloudSet.H"
30#include "polyMesh.H"
32#include "treeBoundBox.H"
33#include "treeDataFace.H"
34#include "Time.H"
35#include "meshTools.H"
36// For 'nearInfo' helper class only
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void Foam::patchCloudSet::calcSamples
51(
52 DynamicList<point>& samplingPts,
53 DynamicList<label>& samplingCells,
54 DynamicList<label>& samplingFaces,
55 DynamicList<label>& samplingSegments,
56 DynamicList<scalar>& samplingCurveDist
57) const
58{
59 DebugInfo << "patchCloudSet : sampling on patches :" << endl;
60
61 // Construct search tree for all patch faces.
62 label sz = 0;
63 for (const label patchi : patchSet_)
64 {
65 const polyPatch& pp = mesh().boundaryMesh()[patchi];
66
67 sz += pp.size();
68
69 DebugInfo << " " << pp.name() << " size " << pp.size() << endl;
70 }
71
72 labelList patchFaces(sz);
73 sz = 0;
74 treeBoundBox bb;
75 for (const label patchi : patchSet_)
76 {
77 const polyPatch& pp = mesh().boundaryMesh()[patchi];
78
79 forAll(pp, i)
80 {
81 patchFaces[sz++] = pp.start()+i;
82 }
83
84 // Without reduction.
85 bb.add(pp.points(), pp.meshPoints());
86 }
87
88 // Not very random
89 Random rndGen(123456);
90 // Make bb asymmetric just to avoid problems on symmetric meshes
91 // Make sure bb is 3D.
92 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
93
94
95 indexedOctree<treeDataFace> patchTree
96 (
97 treeDataFace(mesh(), patchFaces), // boundary faces only
98
99 bb, // overall search domain
100 8, // maxLevel
101 10, // leafsize
102 3.0 // duplicity
103 );
104
105 // Force calculation of face-diagonal decomposition
106 (void)mesh().tetBasePtIs();
107
108
109 // All the info for nearest. Construct to miss
110 List<mappedPatchBase::nearInfo> nearest(sampleCoords_.size());
111
112 forAll(sampleCoords_, sampleI)
113 {
114 const auto& treeData = patchTree.shapes();
115 const point& sample = sampleCoords_[sampleI];
116
117 pointIndexHit& nearInfo = nearest[sampleI].first();
118 auto& distSqrProc = nearest[sampleI].second();
119
120 // Find the nearest locally
121 if (patchFaces.size())
122 {
123 nearInfo = patchTree.findNearest(sample, sqr(searchDist_));
124 }
125 else
126 {
127 nearInfo.setMiss();
128 }
129
130
131 // Fill in the distance field and the processor field
132 if (!nearInfo.hit())
133 {
134 distSqrProc.first() = Foam::sqr(GREAT);
135 distSqrProc.second() = Pstream::myProcNo();
136 }
137 else
138 {
139 // Set nearest to mesh face label
140 const label objectIndex = treeData.objectIndex(nearInfo.index());
141
142 nearInfo.setIndex(objectIndex);
143
144 distSqrProc.first() = sample.distSqr(nearInfo.point());
145 distSqrProc.second() = Pstream::myProcNo();
146 }
147 }
148
149
150 // Find nearest - globally consistent
151 Pstream::listCombineReduce(nearest, mappedPatchBase::nearestEqOp());
152
153
154 if (debug && Pstream::master())
155 {
156 OFstream str
157 (
158 mesh().time().path()
159 / name()
160 + "_nearest.obj"
161 );
162 Info<< "Dumping mapping as lines from supplied points to"
163 << " nearest patch face to file " << str.name() << endl;
164
165 label vertI = 0;
166
167 forAll(nearest, i)
168 {
169 if (nearest[i].first().hit())
170 {
171 meshTools::writeOBJ(str, sampleCoords_[i]);
172 ++vertI;
173 meshTools::writeOBJ(str, nearest[i].first().point());
174 ++vertI;
175 str << "l " << vertI-1 << ' ' << vertI << nl;
176 }
177 }
178 }
179
180
181 // Store the sampling locations on the nearest processor
182 forAll(nearest, sampleI)
183 {
184 const pointIndexHit& nearInfo = nearest[sampleI].first();
185
186 if (nearInfo.hit())
187 {
188 if (nearest[sampleI].second().second() == Pstream::myProcNo())
189 {
190 label facei = nearInfo.index();
191
192 samplingPts.append(nearInfo.point());
193 samplingCells.append(mesh().faceOwner()[facei]);
194 samplingFaces.append(facei);
195 samplingSegments.append(0);
196 samplingCurveDist.append(1.0 * sampleI);
197 }
198 }
199 else
200 {
201 // No processor found point near enough. Mark with special value
202 // which is intercepted when interpolating
203 if (Pstream::master())
204 {
205 samplingPts.append(sampleCoords_[sampleI]);
206 samplingCells.append(-1);
207 samplingFaces.append(-1);
208 samplingSegments.append(0);
209 samplingCurveDist.append(1.0 * sampleI);
210 }
211 }
212 }
213}
214
215
216void Foam::patchCloudSet::genSamples()
217{
218 // Storage for sample points
219 DynamicList<point> samplingPts;
220 DynamicList<label> samplingCells;
221 DynamicList<label> samplingFaces;
222 DynamicList<label> samplingSegments;
223 DynamicList<scalar> samplingCurveDist;
224
225 calcSamples
226 (
227 samplingPts,
228 samplingCells,
229 samplingFaces,
230 samplingSegments,
231 samplingCurveDist
232 );
233
234 samplingPts.shrink();
235 samplingCells.shrink();
236 samplingFaces.shrink();
237 samplingSegments.shrink();
238 samplingCurveDist.shrink();
239
240 setSamples
241 (
242 samplingPts,
243 samplingCells,
244 samplingFaces,
245 samplingSegments,
246 samplingCurveDist
247 );
248
249 if (debug)
250 {
252 }
253}
254
255
256// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
257
259(
260 const word& name,
261 const polyMesh& mesh,
262 const meshSearch& searchEngine,
263 const word& axis,
264 const List<point>& sampleCoords,
265 const labelHashSet& patchSet,
266 const scalar searchDist
267)
268:
269 sampledSet(name, mesh, searchEngine, axis),
270 sampleCoords_(sampleCoords),
271 patchSet_(patchSet),
272 searchDist_(searchDist)
273{
274 genSamples();
275}
276
277
279(
280 const word& name,
281 const polyMesh& mesh,
282 const meshSearch& searchEngine,
283 const dictionary& dict
284)
285:
286 sampledSet(name, mesh, searchEngine, dict),
287 sampleCoords_(dict.get<pointField>("points")),
288 patchSet_
289 (
290 mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
291 ),
292 searchDist_(dict.get<scalar>("maxDistance"))
293{
294 genSamples();
295}
296
297
298// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static void listCombineReduce(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listReduce with an in-place cop.
vector & first()
Definition UList.H:957
bool get(const label i) const
Definition UList.H:868
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const word & axis() const
The sort axis name.
Definition coordSet.H:160
const word & name() const noexcept
The coord-set name.
Definition coordSet.H:152
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
Like Foam::cloudSet but samples nearest patch face.
patchCloudSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const List< point > &sampleCoords, const labelHashSet &patchSet, const scalar searchDist)
Construct from components.
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition polyMesh.C:884
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition sampledSet.H:82
const meshSearch & searchEngine() const noexcept
Definition sampledSet.H:378
sampledSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const coordSet::coordFormat axisType)
Construct from components.
Definition sampledSet.C:405
const polyMesh & mesh() const noexcept
Definition sampledSet.H:373
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
dynamicFvMesh & mesh
auto & name
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
runTime write()
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen