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 -------------------------------------------------------------------------------
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"
32 #include "treeBoundBox.H"
33 #include "treeDataFace.H"
34 #include "Time.H"
35 #include "meshTools.H"
36 #include "mappedPatchBase.H"
37 #include "indirectPrimitivePatch.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(patchSeedSet, 0);
44  addToRunTimeSelectionTable(sampledSet, patchSeedSet, word);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void 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 
103  const indirectPrimitivePatch pp
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. Combine on master.
181  Pstream::listCombineGather(nearest, mappedPatchBase::nearestEqOp());
183 
184 
185  // 3. Pick up my local faces that have won
186 
187  forAll(nearest, sampleI)
188  {
189  if (nearest[sampleI].first().hit())
190  {
191  label procI = nearest[sampleI].second().second();
192  label index = nearest[sampleI].first().index();
193 
194  if (procI == Pstream::myProcNo())
195  {
196  newPatchFaces.append(pp.addressing()[index]);
197  }
198  }
199  }
200  }
201 
202  if (debug)
203  {
204  Pout<< "Found " << newPatchFaces.size()
205  << " out of " << selectedLocations_.size()
206  << " on local processor" << endl;
207  }
208 
209  patchFaces.transfer(newPatchFaces);
210  }
211 
212 
213  // Shuffle and truncate if in random mode
214  label totalSize = returnReduce(patchFaces.size(), sumOp<label>());
215 
216  if (maxPoints_ < totalSize)
217  {
218  // Check what fraction of maxPoints_ I need to generate locally.
219  label myMaxPoints =
220  label(scalar(patchFaces.size())/totalSize*maxPoints_);
221 
222  labelList subset = identity(patchFaces.size());
223  for (label iter = 0; iter < 4; ++iter)
224  {
225  forAll(subset, i)
226  {
227  label j = rndGen.position<label>(0, subset.size()-1);
228  std::swap(subset[i], subset[j]);
229  }
230  }
231  // Truncate
232  subset.setSize(myMaxPoints);
233 
234  // Subset patchFaces
235  patchFaces = labelUIndList(patchFaces, subset)();
236 
237  if (debug)
238  {
239  Pout<< "In random mode : selected " << patchFaces.size()
240  << " faces out of " << patchFaces.size() << endl;
241  }
242  }
243 
244 
245  // Get points on patchFaces.
246  globalIndex globalSampleNumbers(patchFaces.size());
247 
248  samplingPts.setCapacity(patchFaces.size());
249  samplingCells.setCapacity(patchFaces.size());
250  samplingFaces.setCapacity(patchFaces.size());
251  samplingSegments.setCapacity(patchFaces.size());
252  samplingCurveDist.setCapacity(patchFaces.size());
253 
254  // For calculation of min-decomp tet base points
255  (void)mesh().tetBasePtIs();
256 
257  forAll(patchFaces, i)
258  {
259  label facei = patchFaces[i];
260 
261  // Slightly shift point in since on warped face face-diagonal
262  // decomposition might be outside cell for face-centre decomposition!
264  (
265  mesh(),
266  facei,
268  );
269  label celli = mesh().faceOwner()[facei];
270 
271  if (info.hit())
272  {
273  // Move the point into the cell
274  const point& cc = mesh().cellCentres()[celli];
275  samplingPts.append
276  (
277  info.hitPoint() + 1e-1*(cc-info.hitPoint())
278  );
279  }
280  else
281  {
282  samplingPts.append(info.rawPoint());
283  }
284  samplingCells.append(celli);
285  samplingFaces.append(facei);
286  samplingSegments.append(0);
287  samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
288  }
289 }
290 
291 
292 void Foam::patchSeedSet::genSamples()
293 {
294  // Storage for sample points
295  DynamicList<point> samplingPts;
296  DynamicList<label> samplingCells;
297  DynamicList<label> samplingFaces;
298  DynamicList<label> samplingSegments;
299  DynamicList<scalar> samplingCurveDist;
300 
301  calcSamples
302  (
303  samplingPts,
304  samplingCells,
305  samplingFaces,
306  samplingSegments,
307  samplingCurveDist
308  );
309 
310  samplingPts.shrink();
311  samplingCells.shrink();
312  samplingFaces.shrink();
313  samplingSegments.shrink();
314  samplingCurveDist.shrink();
315 
316  // Move into *this
317  setSamples
318  (
319  std::move(samplingPts),
320  std::move(samplingCells),
321  std::move(samplingFaces),
322  std::move(samplingSegments),
323  std::move(samplingCurveDist)
324  );
325 
326  if (debug)
327  {
328  write(Info);
329  }
330 }
331 
332 
333 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
334 
336 (
337  const word& name,
338  const polyMesh& mesh,
339  const meshSearch& searchEngine,
340  const dictionary& dict
341 )
342 :
343  sampledSet(name, mesh, searchEngine, dict),
344  patchSet_
345  (
346  mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
347  ),
348  maxPoints_(dict.get<label>("maxPoints")),
349  selectedLocations_
350  (
352  (
353  "points",
354  pointField(0)
355  )
356  )
357 {
358  genSamples();
359 }
360 
361 
362 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
meshTools.H
Foam::polyMesh::FACE_DIAG_TRIS
Definition: polyMesh.H:107
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
patchSeedSet.H
Foam::subset
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
Foam::patchSeedSet::patchSeedSet
patchSeedSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: patchSeedSet.C:336
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
treeDataFace.H
Foam::Field< vector >
treeBoundBox.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
indirectPrimitivePatch.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::sampledSet::faces
const labelList & faces() const
Definition: sampledSet.H:301
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::indirectPrimitivePatch
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Definition: indirectPrimitivePatch.H:49
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Time.H
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::sampledSet::mesh
const polyMesh & mesh() const
Definition: sampledSet.H:281
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:864
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
rndGen
Random rndGen
Definition: createFields.H:23
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::mappedPatchBase::facePoint
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
Definition: mappedPatchBase.C:1723
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
sample
Minimal example by using system/controlDict.functions:
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58
mappedPatchBase.H
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:892