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 | Copyright (C) 2012-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "patchSeedSet.H"
27 #include "polyMesh.H"
29 #include "treeBoundBox.H"
30 #include "treeDataFace.H"
31 #include "Time.H"
32 #include "meshTools.H"
33 #include "mappedPatchBase.H"
34 #include "indirectPrimitivePatch.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(patchSeedSet, 0);
41  addToRunTimeSelectionTable(sampledSet, patchSeedSet, word);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
48 (
49  DynamicList<point>& samplingPts,
50  DynamicList<label>& samplingCells,
51  DynamicList<label>& samplingFaces,
52  DynamicList<label>& samplingSegments,
53  DynamicList<scalar>& samplingCurveDist
54 )
55 {
56  if (debug)
57  {
58  Info<< "patchSeedSet : sampling on patches :" << endl;
59  }
60 
61  // Construct search tree for all patch faces.
62  label sz = 0;
63  forAllConstIter(labelHashSet, patchSet_, iter)
64  {
65  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
66 
67  sz += pp.size();
68 
69  if (debug)
70  {
71  Info<< " " << pp.name() << " size " << pp.size() << endl;
72  }
73  }
74 
76  sz = 0;
77  forAllConstIter(labelHashSet, patchSet_, iter)
78  {
79  const polyPatch& pp = mesh().boundaryMesh()[iter.key()];
80  forAll(pp, i)
81  {
82  patchFaces[sz++] = pp.start()+i;
83  }
84  }
85 
86 
87  if (!rndGenPtr_.valid())
88  {
89  rndGenPtr_.reset(new Random(0));
90  }
91  Random& rndGen = rndGenPtr_();
92 
93 
94  if (selectedLocations_.size())
95  {
96  DynamicList<label> newPatchFaces(patchFaces.size());
97 
98  // Find the nearest patch face
99  {
100  // 1. All processors find nearest local patch face for all
101  // selectedLocations
102 
103  // All the info for nearest. Construct to miss
104  List<mappedPatchBase::nearInfo> nearest(selectedLocations_.size());
105 
106  const indirectPrimitivePatch pp
107  (
108  IndirectList<face>(mesh().faces(), patchFaces),
109  mesh().points()
110  );
111 
112  treeBoundBox patchBb
113  (
115  (
116  rndGen,
117  1e-4
118  )
119  );
120  patchBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
121  patchBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
122 
123  indexedOctree<treeDataFace> boundaryTree
124  (
125  treeDataFace // all information needed to search faces
126  (
127  false, // do not cache bb
128  mesh(),
129  patchFaces // boundary faces only
130  ),
131  patchBb, // overall search domain
132  8, // maxLevel
133  10, // leafsize
134  3.0 // duplicity
135  );
136 
137  // Get some global dimension so all points are equally likely
138  // to be found
139  const scalar globalDistSqr
140  (
141  //magSqr
142  //(
143  // boundBox
144  // (
145  // pp.points(),
146  // pp.meshPoints(),
147  // true
148  // ).span()
149  //)
150  GREAT
151  );
152 
153  forAll(selectedLocations_, sampleI)
154  {
155  const point& sample = selectedLocations_[sampleI];
156 
157  pointIndexHit& nearInfo = nearest[sampleI].first();
158  nearInfo = boundaryTree.findNearest
159  (
160  sample,
161  globalDistSqr
162  );
163 
164  if (!nearInfo.hit())
165  {
166  nearest[sampleI].second().first() = Foam::sqr(GREAT);
167  nearest[sampleI].second().second() =
168  Pstream::myProcNo();
169  }
170  else
171  {
172  point fc(pp[nearInfo.index()].centre(pp.points()));
173  nearInfo.setPoint(fc);
174  nearest[sampleI].second().first() = magSqr(fc-sample);
175  nearest[sampleI].second().second() =
176  Pstream::myProcNo();
177  }
178  }
179 
180 
181  // 2. Reduce on master. Select nearest processor.
182 
183  // Find nearest. Combine on master.
184  Pstream::listCombineGather(nearest, mappedPatchBase::nearestEqOp());
185  Pstream::listCombineScatter(nearest);
186 
187 
188  // 3. Pick up my local faces that have won
189 
190  forAll(nearest, sampleI)
191  {
192  if (nearest[sampleI].first().hit())
193  {
194  label procI = nearest[sampleI].second().second();
195  label index = nearest[sampleI].first().index();
196 
197  if (procI == Pstream::myProcNo())
198  {
199  newPatchFaces.append(pp.addressing()[index]);
200  }
201  }
202  }
203  }
204 
205  if (debug)
206  {
207  Pout<< "Found " << newPatchFaces.size()
208  << " out of " << selectedLocations_.size()
209  << " on local processor" << endl;
210  }
211 
212  patchFaces.transfer(newPatchFaces);
213  }
214 
215 
216  // Shuffle and truncate if in random mode
217  label totalSize = returnReduce(patchFaces.size(), sumOp<label>());
218 
219  if (maxPoints_ < totalSize)
220  {
221  // Check what fraction of maxPoints_ I need to generate locally.
222  label myMaxPoints =
223  label(scalar(patchFaces.size())/totalSize*maxPoints_);
224 
226  for (label iter = 0; iter < 4; iter++)
227  {
228  forAll(subset, i)
229  {
230  label j = rndGen.integer(0, subset.size()-1);
231  Swap(subset[i], subset[j]);
232  }
233  }
234  // Truncate
235  subset.setSize(myMaxPoints);
236 
237  // Subset patchFaces
239 
240  if (debug)
241  {
242  Pout<< "In random mode : selected " << patchFaces.size()
243  << " faces out of " << patchFaces.size() << endl;
244  }
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,
270  polyMesh::FACE_DIAG_TRIS
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 
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 
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  setSamples
320  (
321  samplingPts,
322  samplingCells,
323  samplingFaces,
324  samplingSegments,
325  samplingCurveDist
326  );
327 }
328 
329 
330 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
331 
333 (
334  const word& name,
335  const polyMesh& mesh,
336  const meshSearch& searchEngine,
337  const dictionary& dict
338 )
339 :
340  sampledSet(name, mesh, searchEngine, dict),
341  patchSet_
342  (
344  (
345  wordReList(dict.lookup("patches"))
346  )
347  ),
348  maxPoints_(readLabel(dict.lookup("maxPoints"))),
349  selectedLocations_
350  (
352  (
353  "points",
354  pointField(0)
355  )
356  )
357 {
358  genSamples();
359 
360  if (debug)
361  {
362  write(Info);
363  }
364 }
365 
366 
367 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
368 
370 {}
371 
372 
373 // ************************************************************************* //
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
Foam::Random
Simple random number generator.
Definition: Random.H:49
meshTools.H
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:57
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
patchSeedSet.H
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::patchSeedSet::patchSeedSet
patchSeedSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: patchSeedSet.C:333
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mappedPatchBase::nearestEqOp
Definition: mappedPatchBase.H:139
polyMesh.H
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
treeDataFace.H
Foam::facePoint
label facePoint(const int facei, const block &block, const label i, const label j)
Definition: blockMeshMergeFast.C:202
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::PointIndexHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointIndexHit.H:143
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::treeBoundBox::extend
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:317
treeBoundBox.H
Foam::Info
messageStream Info
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::DynamicList::setCapacity
void setCapacity(const label)
Alter the size of the underlying storage.
Definition: DynamicListI.H:109
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:102
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::patchSeedSet::genSamples
void genSamples()
Uses calcSamples to obtain samples. Copies them into *this.
Definition: patchSeedSet.C:295
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::subset
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
Definition: ListOpsTemplates.C:290
Foam::patchSeedSet::~patchSeedSet
virtual ~patchSeedSet()
Destructor.
Definition: patchSeedSet.C:369
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sampledSet::setSamples
void setSamples(const List< point > &samplingPts, const labelList &samplingCells, const labelList &samplingFaces, const labelList &samplingSegments, const scalarList &samplingCurveDist)
Sets sample data.
Definition: sampledSet.C:351
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::indexedOctree::findNearest
void findNearest(const label nodeI, const linePointRef &ln, treeBoundBox &tightest, label &nearestShapeI, point &linePoint, point &nearestPoint, const FindNearestOp &fnOp) const
Find nearest point to line.
Definition: indexedOctree.C:590
Foam::patchSeedSet::calcSamples
void calcSamples(DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< label > &samplingSegments, DynamicList< scalar > &samplingCurveDist)
Samples all points in sampleCoords.
Definition: patchSeedSet.C:48
Foam::Tuple2
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:750
Foam::Swap
void Swap(T &a, T &b)
Definition: Swap.H:43
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
write
Tcoeff write()
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
mappedPatchBase.H