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-2022 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())
112  .extend(rndGen, 1e-4, ROOTVSMALL)
113  );
114 
115  indexedOctree<treeDataFace> boundaryTree
116  (
117  treeDataFace(mesh(), patchFaces), // boundary faces only
118 
119  patchBb, // overall search domain
120  8, // maxLevel
121  10, // leafsize
122  3.0 // duplicity
123  );
124 
125  // Get some global dimension so all points are equally likely
126  // to be found
127  const scalar globalDistSqr
128  (
129  //boundBox(pp.points(), pp.meshPoints(), true).magSqr()
130  GREAT
131  );
132 
133  forAll(selectedLocations_, sampleI)
134  {
135  const auto& treeData = boundaryTree.shapes();
136  const point& sample = selectedLocations_[sampleI];
137 
138  pointIndexHit& nearInfo = nearest[sampleI].first();
139  auto& distSqrProc = nearest[sampleI].second();
140 
141  nearInfo = boundaryTree.findNearest
142  (
143  sample,
144  globalDistSqr
145  );
146 
147  if (!nearInfo.hit())
148  {
149  distSqrProc.first() = Foam::sqr(GREAT);
150  distSqrProc.second() = Pstream::myProcNo();
151  }
152  else
153  {
154  nearInfo.setPoint(treeData.centre(nearInfo.index()));
155 
156  distSqrProc.first() = sample.distSqr(nearInfo.point());
157  distSqrProc.second() = Pstream::myProcNo();
158  }
159  }
160 
161 
162  // 2. Reduce on master. Select nearest processor.
163 
164  // Find nearest - globally consistent
165  Pstream::listCombineReduce(nearest, mappedPatchBase::nearestEqOp());
166 
167  // 3. Pick up my local faces that have won
168 
169  forAll(nearest, sampleI)
170  {
171  if (nearest[sampleI].first().hit())
172  {
173  label procI = nearest[sampleI].second().second();
174  label index = nearest[sampleI].first().index();
175 
176  if (procI == Pstream::myProcNo())
177  {
178  newPatchFaces.append(pp.addressing()[index]);
179  }
180  }
181  }
182  }
183 
184  if (debug)
185  {
186  Pout<< "Found " << newPatchFaces.size()
187  << " out of " << selectedLocations_.size()
188  << " on local processor" << endl;
189  }
190 
191  patchFaces.transfer(newPatchFaces);
192  }
193 
194 
195  // Shuffle and truncate if in random mode
196  label totalSize = returnReduce(patchFaces.size(), sumOp<label>());
197 
198  if (maxPoints_ < totalSize)
199  {
200  // Check what fraction of maxPoints_ I need to generate locally.
201  label myMaxPoints =
202  label(scalar(patchFaces.size())/totalSize*maxPoints_);
203 
204  labelList subset = identity(patchFaces.size());
205  for (label iter = 0; iter < 4; ++iter)
206  {
207  forAll(subset, i)
208  {
209  label j = rndGen.position<label>(0, subset.size()-1);
210  std::swap(subset[i], subset[j]);
211  }
212  }
213  // Truncate
214  subset.setSize(myMaxPoints);
215 
216  // Subset patchFaces
217 
218  if (debug)
219  {
220  Pout<< "In random mode : selected " << subset.size()
221  << " faces out of " << patchFaces.size() << endl;
222  }
223 
224  patchFaces = labelUIndList(patchFaces, subset)();
225  }
226 
227 
228  // Get points on patchFaces.
229  globalIndex globalSampleNumbers(patchFaces.size());
230 
231  samplingPts.setCapacity(patchFaces.size());
232  samplingCells.setCapacity(patchFaces.size());
233  samplingFaces.setCapacity(patchFaces.size());
234  samplingSegments.setCapacity(patchFaces.size());
235  samplingCurveDist.setCapacity(patchFaces.size());
236 
237  // For calculation of min-decomp tet base points
238  (void)mesh().tetBasePtIs();
239 
240  forAll(patchFaces, i)
241  {
242  label facei = patchFaces[i];
243 
244  // Slightly shift point in since on warped face face-diagonal
245  // decomposition might be outside cell for face-centre decomposition!
247  (
248  mesh(),
249  facei,
251  );
252  label celli = mesh().faceOwner()[facei];
253 
254  if (info.hit())
255  {
256  // Move the point into the cell
257  const point& cc = mesh().cellCentres()[celli];
258  samplingPts.append
259  (
260  info.point() + 1e-1*(cc-info.point())
261  );
262  }
263  else
264  {
265  samplingPts.append(info.point());
266  }
267  samplingCells.append(celli);
268  samplingFaces.append(facei);
269  samplingSegments.append(0);
270  samplingCurveDist.append(globalSampleNumbers.toGlobal(i));
271  }
272 }
273 
274 
275 void Foam::patchSeedSet::genSamples()
276 {
277  // Storage for sample points
278  DynamicList<point> samplingPts;
279  DynamicList<label> samplingCells;
280  DynamicList<label> samplingFaces;
281  DynamicList<label> samplingSegments;
282  DynamicList<scalar> samplingCurveDist;
283 
284  calcSamples
285  (
286  samplingPts,
287  samplingCells,
288  samplingFaces,
289  samplingSegments,
290  samplingCurveDist
291  );
292 
293  samplingPts.shrink();
294  samplingCells.shrink();
295  samplingFaces.shrink();
296  samplingSegments.shrink();
297  samplingCurveDist.shrink();
298 
299  // Move into *this
300  setSamples
301  (
302  std::move(samplingPts),
303  std::move(samplingCells),
304  std::move(samplingFaces),
305  std::move(samplingSegments),
306  std::move(samplingCurveDist)
307  );
308 
309  if (debug)
310  {
311  write(Info);
312  }
313 }
314 
315 
316 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
317 
319 (
320  const word& name,
321  const polyMesh& mesh,
322  const meshSearch& searchEngine,
323  const dictionary& dict
324 )
325 :
326  sampledSet(name, mesh, searchEngine, dict),
327  patchSet_
328  (
329  mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
330  ),
331  maxPoints_(dict.get<label>("maxPoints")),
332  selectedLocations_
333  (
334  dict.getOrDefault<pointField>
335  (
336  "points",
337  pointField(0)
338  )
339  )
340 {
341  genSamples();
342 }
343 
344 
345 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:491
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:893
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Random rndGen
Definition: createFields.H:23
T & first()
Access first element of the list, position [0].
Definition: UList.H:798
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
const labelList & faces() const noexcept
Definition: sampledSet.H:393
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
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.
patchSeedSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: patchSeedSet.C:312
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
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 INVALID.
Definition: exprTraits.C:52
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1104
const pointField & points() const noexcept
Return the points.
Definition: coordSet.H:168
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:
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
vector point
Point is a vector.
Definition: point.H:37
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
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.