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-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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"
31 #include "treeBoundBox.H"
32 #include "treeDataFace.H"
33 #include "mappedPatchBase.H"
34 #include "indirectPrimitivePatch.H"
35 #include "triangulatedPatch.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(patchSeedSet, 0);
43  addToRunTimeSelectionTable(sampledSet, patchSeedSet, word);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::patchSeedSet::calcPatchSamples
50 (
51  const label nAvailable,
52  const label nPatchPoints,
53  DynamicList<point>& samplingPts,
54  DynamicList<label>& samplingCells,
55  DynamicList<label>& samplingFaces,
56  DynamicList<label>& samplingSegments,
57  DynamicList<scalar>& samplingCurveDist
58 )
59 {
60  if (nAvailable < 1)
61  {
62  return;
63  }
64 
65  Random& rndGen = *rndGenPtr_;
66 
67  globalIndex globalSampleNumbers(nAvailable);
68  label nGlobalPatchPoints = returnReduce(nPatchPoints, sumOp<label>());
69 
70  point pt;
71  label facei;
72  label celli;
73 
74  const bool perturb = true;
75 
76  for (const label patchi : patchSet_)
77  {
78  const polyPatch& pp = mesh().boundaryMesh()[patchi];
79  triangulatedPatch tp(pp, perturb);
80 
81  const label np = nAvailable*pp.size()/scalar(nGlobalPatchPoints);
82  for (label i = 0; i < np; ++i)
83  {
84  tp.randomLocalPoint(rndGen, pt, facei, celli);
85 
86  samplingPts.append(pt);
87  samplingCells.append(celli);
88  samplingFaces.append(pp.start() + facei);
89  samplingSegments.append(0);
90  samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
91  }
92  }
93 }
94 
95 
96 void Foam::patchSeedSet::calcSelectedLocations
97 (
98  const label nAvailable,
99  const label nPatchPoints,
100  DynamicList<point>& samplingPts,
101  DynamicList<label>& samplingCells,
102  DynamicList<label>& samplingFaces,
103  DynamicList<label>& samplingSegments,
104  DynamicList<scalar>& samplingCurveDist
105 )
106 {
107  if (nAvailable < 1)
108  {
109  return;
110  }
111 
112  Random& rndGen = *rndGenPtr_;
113 
114  labelList patchFaces(nPatchPoints);
115  label sz = 0;
116  for (const label patchi : patchSet_)
117  {
118  const polyPatch& pp = mesh().boundaryMesh()[patchi];
119  forAll(pp, localFacei)
120  {
121  patchFaces[sz++] = pp.start() + localFacei;
122  }
123  }
124 
125  {
126  DynamicList<label> newPatchFaces(patchFaces.size());
127 
128  // Find the nearest patch face
129  {
130  // 1. All processors find nearest local patch face for all
131  // selectedLocations
132 
133  // All the info for nearest. Construct to miss
134  List<mappedPatchBase::nearInfo> nearest(nAvailable);
135 
137  (
138  IndirectList<face>(mesh().faces(), patchFaces),
139  mesh().points()
140  );
141 
142  treeBoundBox patchBb
143  (
144  treeBoundBox(pp.points(), pp.meshPoints())
145  .extend(rndGen, 1e-4, ROOTVSMALL)
146  );
147 
148  indexedOctree<treeDataFace> boundaryTree
149  (
150  treeDataFace(mesh(), patchFaces), // boundary faces only
151 
152  patchBb, // overall search domain
153  8, // maxLevel
154  10, // leafsize
155  3.0 // duplicity
156  );
157 
158  // Get some global dimension so all points are equally likely
159  // to be found
160  const scalar globalDistSqr
161  (
162  //boundBox(pp.points(), pp.meshPoints(), true).magSqr()
163  GREAT
164  );
165 
166  for (label sampleI = 0; sampleI < nAvailable; ++sampleI)
167  {
168  const auto& treeData = boundaryTree.shapes();
169  const point& sample = selectedLocations_[sampleI];
170 
171  pointIndexHit& nearInfo = nearest[sampleI].first();
172  auto& distSqrProc = nearest[sampleI].second();
173 
174  nearInfo = boundaryTree.findNearest
175  (
176  sample,
177  globalDistSqr
178  );
179 
180  if (!nearInfo.hit())
181  {
182  distSqrProc.first() = Foam::sqr(GREAT);
183  distSqrProc.second() = Pstream::myProcNo();
184  }
185  else
186  {
187  nearInfo.setPoint(treeData.centre(nearInfo.index()));
188 
189  distSqrProc.first() = sample.distSqr(nearInfo.point());
190  distSqrProc.second() = Pstream::myProcNo();
191  }
192  }
193 
194 
195  // 2. Reduce on master. Select nearest processor.
196 
197  // Find nearest - globally consistent
198  Pstream::listCombineReduce(nearest, mappedPatchBase::nearestEqOp());
199 
200  // 3. Pick up my local faces that have won
201 
202  forAll(nearest, sampleI)
203  {
204  if (nearest[sampleI].first().hit())
205  {
206  label procI = nearest[sampleI].second().second();
207  label index = nearest[sampleI].first().index();
208 
209  if (procI == Pstream::myProcNo())
210  {
211  newPatchFaces.append(pp.addressing()[index]);
212  }
213  }
214  }
215  }
216 
217  if (debug)
218  {
219  Pout<< "Found " << newPatchFaces.size()
220  << " out of " << nAvailable
221  << " on local processor" << endl;
222  }
223 
224  patchFaces.transfer(newPatchFaces);
225  }
226 
227 
228  // Shuffle and truncate if in random mode
229  const label totalSize = returnReduce(patchFaces.size(), sumOp<label>());
230 
231  if (totalSize > nAvailable)
232  {
233  // Check what fraction of maxPoints_ I need to generate locally.
234  label myMaxPoints = scalar(patchFaces.size())/totalSize*nAvailable;
235 
236  labelList subset = identity(patchFaces.size());
237  for (label iter = 0; iter < 4; ++iter)
238  {
239  forAll(subset, i)
240  {
241  label j = rndGen.position<label>(0, subset.size()-1);
242  std::swap(subset[i], subset[j]);
243  }
244  }
245 
246  // Truncate
247  subset.setSize(myMaxPoints);
248 
249  // Subset patchFaces
250 
251  if (debug)
252  {
253  Pout<< "In random mode : selected " << subset.size()
254  << " faces out of " << patchFaces.size() << endl;
255  }
256 
257  patchFaces = labelUIndList(patchFaces, subset)();
258  }
259 
260 
261  // Get points on patchFaces.
262  globalIndex globalSampleNumbers(patchFaces.size());
263 
264  samplingPts.setCapacity(patchFaces.size());
265  samplingCells.setCapacity(patchFaces.size());
266  samplingFaces.setCapacity(patchFaces.size());
267  samplingSegments.setCapacity(patchFaces.size());
268  samplingCurveDist.setCapacity(patchFaces.size());
269 
270  // For calculation of min-decomp tet base points
271  (void)mesh().tetBasePtIs();
272 
273  forAll(patchFaces, i)
274  {
275  const label facei = patchFaces[i];
276 
277  // Slightly shift point in since on warped face face-diagonal
278  // decomposition might be outside cell for face-centre decomposition!
280  (
281  mesh(),
282  facei,
284  );
285  const label celli = mesh().faceOwner()[facei];
286 
287  if (info.hit())
288  {
289  // Move the point into the cell
290  const point& cc = mesh().cellCentres()[celli];
291  samplingPts.append
292  (
293  info.point() + 1e-1*(cc-info.point())
294  );
295  }
296  else
297  {
298  samplingPts.append(info.point());
299  }
300  samplingCells.append(celli);
301  samplingFaces.append(facei);
302  samplingSegments.append(0);
303  samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
304  }
305 }
306 
307 
308 void Foam::patchSeedSet::genSamples()
309 {
310  // Storage for sample points
311  DynamicList<point> samplingPts;
312  DynamicList<label> samplingCells;
313  DynamicList<label> samplingFaces;
314  DynamicList<label> samplingSegments;
315  DynamicList<scalar> samplingCurveDist;
316 
317  calcSamples
318  (
319  samplingPts,
320  samplingCells,
321  samplingFaces,
322  samplingSegments,
323  samplingCurveDist
324  );
325 
326  samplingPts.shrink();
327  samplingCells.shrink();
328  samplingFaces.shrink();
329  samplingSegments.shrink();
330  samplingCurveDist.shrink();
331 
332  // Move into *this
333  setSamples
334  (
335  std::move(samplingPts),
336  std::move(samplingCells),
337  std::move(samplingFaces),
338  std::move(samplingSegments),
339  std::move(samplingCurveDist)
340  );
341 
342  if (debug)
343  {
344  write(Info);
345  }
346 }
347 
348 
349 void Foam::patchSeedSet::calcSamples
350 (
351  DynamicList<point>& samplingPts,
352  DynamicList<label>& samplingCells,
353  DynamicList<label>& samplingFaces,
354  DynamicList<label>& samplingSegments,
355  DynamicList<scalar>& samplingCurveDist
356 )
357 {
358  DebugInfo << "patchSeedSet : sampling on patches :" << endl;
359 
360  if (!rndGenPtr_)
361  {
362  rndGenPtr_.reset(new Random(0));
363  }
364 
365  label nPatchPoints = 0;
366  for (const label patchi : patchSet_)
367  {
368  const polyPatch& pp = mesh().boundaryMesh()[patchi];
369  nPatchPoints += pp.size();
370 
371  DebugInfo << " " << pp.name() << " size " << pp.size() << endl;
372  }
373 
374  label nAvailable = min(maxPoints_, selectedLocations_.size());
375 
376  calcSelectedLocations
377  (
378  nAvailable,
379  nPatchPoints,
380  samplingPts,
381  samplingCells,
382  samplingFaces,
383  samplingSegments,
384  samplingCurveDist
385  );
386 
387  nAvailable = maxPoints_ - nAvailable;
388 
389  calcPatchSamples
390  (
391  nAvailable,
392  nPatchPoints,
393  samplingPts,
394  samplingCells,
395  samplingFaces,
396  samplingSegments,
397  samplingCurveDist
398  );
399 }
400 
401 
402 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
403 
405 (
406  const word& name,
407  const polyMesh& mesh,
408  const meshSearch& searchEngine,
409  const dictionary& dict
410 )
411 :
412  sampledSet(name, mesh, searchEngine, dict),
413  patchSet_
414  (
415  mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
416  ),
417  maxPoints_(dict.get<label>("maxPoints")),
418  selectedLocations_
419  (
420  dict.getOrDefault<pointField>
421  (
422  "points",
423  pointField(0)
424  )
425  )
426 {
427  genSamples();
428 }
429 
430 
431 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:56
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
dimensionedSymmTensor sqr(const dimensionedVector &dv)
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Random rndGen
Definition: createFields.H:23
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
patchSeedSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: patchSeedSet.C:398
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:79
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const Field< point_type > & points() const noexcept
Return reference to global points.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const vectorField & cellCentres() const
#define DebugInfo
Report an information message using Foam::Info.
const polyMesh & mesh() const noexcept
Definition: sampledSet.H:373
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
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:37
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.