patchSeedSet.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) 2012-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2021 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 "patchSeedSet.H"
30#include "polyMesh.H"
32#include "treeBoundBox.H"
33#include "treeDataFace.H"
34#include "Time.H"
35#include "meshTools.H"
36#include "mappedPatchBase.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void Foam::patchSeedSet::calcSamples
51(
52 DynamicList<point>& samplingPts,
53 DynamicList<label>& samplingCells,
54 DynamicList<label>& samplingFaces,
55 DynamicList<label>& samplingSegments,
56 DynamicList<scalar>& samplingCurveDist
57)
58{
59 DebugInfo << "patchSeedSet : 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 for (const label patchi : patchSet_)
75 {
76 const polyPatch& pp = mesh().boundaryMesh()[patchi];
77 forAll(pp, i)
78 {
79 patchFaces[sz++] = pp.start()+i;
80 }
81 }
82
83
84 if (!rndGenPtr_)
85 {
86 rndGenPtr_.reset(new Random(0));
87 }
88 Random& rndGen = *rndGenPtr_;
89
90
91 if (selectedLocations_.size())
92 {
93 DynamicList<label> newPatchFaces(patchFaces.size());
94
95 // Find the nearest patch face
96 {
97 // 1. All processors find nearest local patch face for all
98 // selectedLocations
99
100 // All the info for nearest. Construct to miss
101 List<mappedPatchBase::nearInfo> nearest(selectedLocations_.size());
102
104 (
105 IndirectList<face>(mesh().faces(), patchFaces),
106 mesh().points()
107 );
108
109 treeBoundBox patchBb
110 (
111 treeBoundBox(pp.points(), pp.meshPoints()).extend
112 (
113 rndGen,
114 1e-4
115 )
116 );
117 patchBb.min() -= point::uniform(ROOTVSMALL);
118 patchBb.max() += point::uniform(ROOTVSMALL);
119
120 indexedOctree<treeDataFace> boundaryTree
121 (
122 treeDataFace // all information needed to search faces
123 (
124 false, // do not cache bb
125 mesh(),
126 patchFaces // boundary faces only
127 ),
128 patchBb, // overall search domain
129 8, // maxLevel
130 10, // leafsize
131 3.0 // duplicity
132 );
133
134 // Get some global dimension so all points are equally likely
135 // to be found
136 const scalar globalDistSqr
137 (
138 //magSqr
139 //(
140 // boundBox
141 // (
142 // pp.points(),
143 // pp.meshPoints(),
144 // true
145 // ).span()
146 //)
147 GREAT
148 );
149
150 forAll(selectedLocations_, sampleI)
151 {
152 const point& sample = selectedLocations_[sampleI];
153
154 pointIndexHit& nearInfo = nearest[sampleI].first();
155 nearInfo = boundaryTree.findNearest
156 (
157 sample,
158 globalDistSqr
159 );
160
161 if (!nearInfo.hit())
162 {
163 nearest[sampleI].second().first() = Foam::sqr(GREAT);
164 nearest[sampleI].second().second() =
166 }
167 else
168 {
169 point fc(pp[nearInfo.index()].centre(pp.points()));
170 nearInfo.setPoint(fc);
171 nearest[sampleI].second().first() = magSqr(fc-sample);
172 nearest[sampleI].second().second() =
174 }
175 }
176
177
178 // 2. Reduce on master. Select nearest processor.
179
180 // Find nearest - globally consistent
182 (
183 nearest,
184 mappedPatchBase::nearestEqOp()
185 );
186
187 // 3. Pick up my local faces that have won
188
189 forAll(nearest, sampleI)
190 {
191 if (nearest[sampleI].first().hit())
192 {
193 label procI = nearest[sampleI].second().second();
194 label index = nearest[sampleI].first().index();
195
196 if (procI == Pstream::myProcNo())
197 {
198 newPatchFaces.append(pp.addressing()[index]);
199 }
200 }
201 }
202 }
203
204 if (debug)
205 {
206 Pout<< "Found " << newPatchFaces.size()
207 << " out of " << selectedLocations_.size()
208 << " on local processor" << endl;
209 }
210
211 patchFaces.transfer(newPatchFaces);
212 }
213
214
215 // Shuffle and truncate if in random mode
216 label totalSize = returnReduce(patchFaces.size(), sumOp<label>());
217
218 if (maxPoints_ < totalSize)
219 {
220 // Check what fraction of maxPoints_ I need to generate locally.
221 label myMaxPoints =
222 label(scalar(patchFaces.size())/totalSize*maxPoints_);
223
224 labelList subset = identity(patchFaces.size());
225 for (label iter = 0; iter < 4; ++iter)
226 {
227 forAll(subset, i)
228 {
229 label j = rndGen.position<label>(0, subset.size()-1);
230 std::swap(subset[i], subset[j]);
231 }
232 }
233 // Truncate
234 subset.setSize(myMaxPoints);
235
236 // Subset patchFaces
237
238 if (debug)
239 {
240 Pout<< "In random mode : selected " << subset.size()
241 << " faces out of " << patchFaces.size() << endl;
242 }
243
244 patchFaces = labelUIndList(patchFaces, subset)();
245 }
246
247
248 // Get points on patchFaces.
249 globalIndex globalSampleNumbers(patchFaces.size());
250
251 samplingPts.setCapacity(patchFaces.size());
252 samplingCells.setCapacity(patchFaces.size());
253 samplingFaces.setCapacity(patchFaces.size());
254 samplingSegments.setCapacity(patchFaces.size());
255 samplingCurveDist.setCapacity(patchFaces.size());
256
257 // For calculation of min-decomp tet base points
258 (void)mesh().tetBasePtIs();
259
260 forAll(patchFaces, i)
261 {
262 label facei = patchFaces[i];
263
264 // Slightly shift point in since on warped face face-diagonal
265 // decomposition might be outside cell for face-centre decomposition!
267 (
268 mesh(),
269 facei,
271 );
272 label celli = mesh().faceOwner()[facei];
273
274 if (info.hit())
275 {
276 // Move the point into the cell
277 const point& cc = mesh().cellCentres()[celli];
278 samplingPts.append
279 (
280 info.hitPoint() + 1e-1*(cc-info.hitPoint())
281 );
282 }
283 else
284 {
285 samplingPts.append(info.rawPoint());
286 }
287 samplingCells.append(celli);
288 samplingFaces.append(facei);
289 samplingSegments.append(0);
290 samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
291 }
292}
293
294
295void Foam::patchSeedSet::genSamples()
296{
297 // Storage for sample points
298 DynamicList<point> samplingPts;
299 DynamicList<label> samplingCells;
300 DynamicList<label> samplingFaces;
301 DynamicList<label> samplingSegments;
302 DynamicList<scalar> samplingCurveDist;
303
304 calcSamples
305 (
306 samplingPts,
307 samplingCells,
308 samplingFaces,
309 samplingSegments,
310 samplingCurveDist
311 );
312
313 samplingPts.shrink();
314 samplingCells.shrink();
315 samplingFaces.shrink();
316 samplingSegments.shrink();
317 samplingCurveDist.shrink();
318
319 // Move into *this
320 setSamples
321 (
322 std::move(samplingPts),
323 std::move(samplingCells),
324 std::move(samplingFaces),
325 std::move(samplingSegments),
326 std::move(samplingCurveDist)
327 );
328
329 if (debug)
330 {
331 write(Info);
332 }
333}
334
335
336// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
337
339(
340 const word& name,
341 const polyMesh& mesh,
342 const meshSearch& searchEngine,
343 const dictionary& dict
344)
345:
346 sampledSet(name, mesh, searchEngine, dict),
347 patchSet_
348 (
349 mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
350 ),
351 maxPoints_(dict.get<label>("maxPoints")),
352 selectedLocations_
353 (
354 dict.getOrDefault<pointField>
355 (
356 "points",
357 pointField(0)
358 )
359 )
360{
361 genSamples();
362}
363
364
365// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Minimal example by using system/controlDict.functions:
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
const pointField & points() const noexcept
Return the points.
Definition: coordSet.H:146
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
Initialises points on or just off patch.
Definition: patchSeedSet.H:96
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const vectorField & cellCentres() const
int myProcNo() const noexcept
Return processor number.
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:86
const labelList & faces() const noexcept
Definition: sampledSet.H:339
const polyMesh & mesh() const noexcept
Definition: sampledSet.H:319
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:68
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
runTime write()
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23