sampledSet.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) 2011-2017 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 "sampledSet.H"
30 #include "polyMesh.H"
31 #include "primitiveMesh.H"
32 #include "meshSearch.H"
33 #include "writer.H"
34 #include "particle.H"
35 #include "globalIndex.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(sampledSet, 0);
42  defineRunTimeSelectionTable(sampledSet, word);
43 }
44 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
50  if
51  (
52  (cells_.size() != size())
53  || (faces_.size() != size())
54  || (segments_.size() != size())
55  || (curveDist_.size() != size())
56  )
57  {
59  << "sizes not equal : "
60  << " points:" << size()
61  << " cells:" << cells_.size()
62  << " faces:" << faces_.size()
63  << " segments:" << segments_.size()
64  << " curveDist:" << curveDist_.size()
65  << abort(FatalError);
66  }
67 }
68 
69 
70 Foam::label Foam::sampledSet::getBoundaryCell(const label facei) const
71 {
72  return mesh().faceOwner()[facei];
73 }
74 
75 
76 Foam::label Foam::sampledSet::getNeighbourCell(const label facei) const
77 {
78  if (facei >= mesh().nInternalFaces())
79  {
80  return mesh().faceOwner()[facei];
81  }
82  else
83  {
84  return mesh().faceNeighbour()[facei];
85  }
86 }
87 
88 
90 (
91  const point& p,
92  const label samplei
93 ) const
94 {
95  // Collect the face owner and neighbour cells of the sample into an array
96  // for convenience
97  const label cells[4] =
98  {
99  mesh().faceOwner()[faces_[samplei]],
100  getNeighbourCell(faces_[samplei]),
101  mesh().faceOwner()[faces_[samplei+1]],
102  getNeighbourCell(faces_[samplei+1])
103  };
104 
105  // Find the sampled cell by checking the owners and neighbours of the
106  // sampled faces
107  label cellm =
108  (cells[0] == cells[2] || cells[0] == cells[3]) ? cells[0]
109  : (cells[1] == cells[2] || cells[1] == cells[3]) ? cells[1]
110  : -1;
111 
112  if (cellm != -1)
113  {
114  // If found the sampled cell check the point is in the cell
115  // otherwise ignore
116  if (!mesh().pointInCell(p, cellm, searchEngine_.decompMode()))
117  {
118  cellm = -1;
119 
120  if (debug)
121  {
123  << "Could not find mid-point " << p
124  << " cell " << cellm << endl;
125  }
126  }
127  }
128  else
129  {
130  // If the sample does not pass through a single cell check if the point
131  // is in any of the owners or neighbours otherwise ignore
132  for (label i=0; i<4; ++i)
133  {
134  if (mesh().pointInCell(p, cells[i], searchEngine_.decompMode()))
135  {
136  return cells[i];
137  }
138  }
139 
140  if (debug)
141  {
143  << "Could not find cell for mid-point" << nl
144  << " samplei: " << samplei
145  << " pts[samplei]: " << operator[](samplei)
146  << " face[samplei]: " << faces_[samplei]
147  << " pts[samplei+1]: " << operator[](samplei+1)
148  << " face[samplei+1]: " << faces_[samplei+1]
149  << " cellio: " << cells[0]
150  << " cellin: " << cells[1]
151  << " celljo: " << cells[2]
152  << " celljn: " << cells[3]
153  << endl;
154  }
155  }
156 
157  return cellm;
158 }
159 
160 
161 Foam::scalar Foam::sampledSet::calcSign
162 (
163  const label facei,
164  const point& sample
165 ) const
166 {
167  vector vec = sample - mesh().faceCentres()[facei];
168 
169  scalar magVec = mag(vec);
170 
171  if (magVec < VSMALL)
172  {
173  // Sample on face centre. Regard as inside
174  return -1;
175  }
176 
177  vec /= magVec;
178 
179  const vector n = normalised(mesh().faceAreas()[facei]);
180 
181  return n & vec;
182 }
183 
184 
186 (
187  const label celli,
188  const point& sample,
189  const scalar smallDist
190 ) const
191 {
192  const cell& myFaces = mesh().cells()[celli];
193 
194  forAll(myFaces, myFacei)
195  {
196  const face& f = mesh().faces()[myFaces[myFacei]];
197 
198  pointHit inter = f.nearestPoint(sample, mesh().points());
199 
200  scalar dist;
201 
202  if (inter.hit())
203  {
204  dist = mag(inter.hitPoint() - sample);
205  }
206  else
207  {
208  dist = mag(inter.missPoint() - sample);
209  }
210 
211  if (dist < smallDist)
212  {
213  return myFaces[myFacei];
214  }
215  }
216  return -1;
217 }
218 
219 
221 (
222  const point& facePt,
223  const label facei
224 ) const
225 {
226  label celli = mesh().faceOwner()[facei];
227  const point& cC = mesh().cellCentres()[celli];
228 
229  point newPosition = facePt;
230 
231  // Taken from particle::initCellFacePt()
232  label tetFacei;
233  label tetPtI;
234  mesh().findTetFacePt(celli, facePt, tetFacei, tetPtI);
235 
236  // This is the tolerance that was defined as a static constant of the
237  // particle class. It is no longer used by particle, following the switch to
238  // barycentric tracking. The only place in which the tolerance is now used
239  // is here. I'm not sure what the purpose of this code is, but it probably
240  // wants removing. It is doing tet-searches for a particle position. This
241  // should almost certainly be left to the particle class.
242  const scalar trackingCorrectionTol = 1e-5;
243 
244  if (tetFacei == -1 || tetPtI == -1)
245  {
246  newPosition = facePt;
247 
248  label trap(1.0/trackingCorrectionTol + 1);
249 
250  label iterNo = 0;
251 
252  do
253  {
254  newPosition += trackingCorrectionTol*(cC - facePt);
255 
257  (
258  celli,
259  newPosition,
260  tetFacei,
261  tetPtI
262  );
263 
264  ++iterNo;
265 
266  } while (tetFacei < 0 && iterNo <= trap);
267  }
268 
269  if (tetFacei == -1)
270  {
272  << "After pushing " << facePt << " to " << newPosition
273  << " it is still outside face " << facei
274  << " at " << mesh().faceCentres()[facei]
275  << " of cell " << celli
276  << " at " << cC << endl
277  << "Please change your starting point"
278  << abort(FatalError);
279  }
280 
281  return newPosition;
282 }
283 
284 
286 (
287  const point& samplePt,
288  const point& bPoint,
289  const label bFacei,
290  const scalar smallDist,
291 
292  point& trackPt,
293  label& trackCelli,
294  label& trackFacei
295 ) const
296 {
297  bool isGoodSample = false;
298 
299  if (bFacei == -1)
300  {
301  // No boundary intersection. Try and find cell samplePt is in
302  trackCelli = mesh().findCell(samplePt, searchEngine_.decompMode());
303 
304  if
305  (
306  (trackCelli == -1)
307  || !mesh().pointInCell
308  (
309  samplePt,
310  trackCelli,
311  searchEngine_.decompMode()
312  )
313  )
314  {
315  // Line samplePt - end_ does not intersect domain at all.
316  // (or is along edge)
317 
318  trackCelli = -1;
319  trackFacei = -1;
320 
321  isGoodSample = false;
322  }
323  else
324  {
325  // Start is inside. Use it as tracking point
326 
327  trackPt = samplePt;
328  trackFacei = -1;
329 
330  isGoodSample = true;
331  }
332  }
333  else if (mag(samplePt - bPoint) < smallDist)
334  {
335  // samplePt close to bPoint. Snap to it
336  trackPt = pushIn(bPoint, bFacei);
337  trackFacei = bFacei;
338  trackCelli = getBoundaryCell(trackFacei);
339 
340  isGoodSample = true;
341  }
342  else
343  {
344  scalar sign = calcSign(bFacei, samplePt);
345 
346  if (sign < 0)
347  {
348  // samplePt inside or marginally outside.
349  trackPt = samplePt;
350  trackFacei = -1;
351  trackCelli = mesh().findCell(trackPt, searchEngine_.decompMode());
352 
353  isGoodSample = true;
354  }
355  else
356  {
357  // samplePt outside. use bPoint
358  trackPt = pushIn(bPoint, bFacei);
359  trackFacei = bFacei;
360  trackCelli = getBoundaryCell(trackFacei);
361 
362  isGoodSample = false;
363  }
364  }
365 
367  << " samplePt:" << samplePt
368  << " bPoint:" << bPoint
369  << " bFacei:" << bFacei << nl
370  << " Calculated first tracking point :"
371  << " trackPt:" << trackPt
372  << " trackCelli:" << trackCelli
373  << " trackFacei:" << trackFacei
374  << " isGoodSample:" << isGoodSample
375  << endl;
376 
377  return isGoodSample;
378 }
379 
380 
382 (
383  const List<point>& samplingPts,
384  const labelList& samplingCells,
385  const labelList& samplingFaces,
386  const labelList& samplingSegments,
387  const scalarList& samplingCurveDist
388 )
389 {
390  setPoints(samplingPts);
391  curveDist_ = samplingCurveDist;
392 
393  segments_ = samplingSegments;
394  cells_ = samplingCells;
395  faces_ = samplingFaces;
396 
397  checkDimensions();
398 }
399 
400 
402 (
403  List<point>&& samplingPts,
404  labelList&& samplingCells,
405  labelList&& samplingFaces,
406  labelList&& samplingSegments,
407  scalarList&& samplingCurveDist
408 )
409 {
410  setPoints(std::move(samplingPts));
411  curveDist_ = std::move(samplingCurveDist);
412 
413  segments_ = std::move(samplingSegments);
414  cells_ = std::move(samplingCells);
415  faces_ = std::move(samplingFaces);
416 
417  checkDimensions();
418 }
419 
420 
422 (
423  labelList& indexSet,
424  labelList& allSegments
425 ) const
426 {
427  // Combine sampleSet from processors. Sort by curveDist. Return
428  // ordering in indexSet.
429  // Note: only master results are valid
430 
431  List<point> allPts;
432  globalIndex::gatherOp(*this, allPts);
433 
434  globalIndex::gatherOp(segments(), allSegments);
435 
436  scalarList allCurveDist;
437  globalIndex::gatherOp(curveDist(), allCurveDist);
438 
439 
440  if (Pstream::master() && allCurveDist.empty())
441  {
443  << "Sample set " << name()
444  << " has zero points." << endl;
445  }
446 
447  // Sort curveDist and use to fill masterSamplePts
448  Foam::sortedOrder(allCurveDist, indexSet); // uses stable sort
449  scalarList sortedDist(allCurveDist, indexSet); // with indices for mapping
450 
451  allSegments = UIndirectList<label>(allSegments, indexSet)();
452 
454  (
455  name(),
456  axis(),
457  List<point>(UIndirectList<point>(allPts, indexSet)),
458  sortedDist
459  );
460 }
461 
462 
463 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
464 
466 (
467  const word& name,
468  const polyMesh& mesh,
469  const meshSearch& searchEngine,
470  const coordSet::coordFormat axisType
471 )
472 :
473  coordSet(name, axisType),
474  mesh_(mesh),
475  searchEngine_(searchEngine),
476  segments_(),
477  cells_(),
478  faces_()
479 {}
480 
481 
483 (
484  const word& name,
485  const polyMesh& mesh,
486  const meshSearch& searchEngine,
487  const word& axis
488 )
489 :
490  coordSet(name, axis),
491  mesh_(mesh),
492  searchEngine_(searchEngine),
493  segments_(),
494  cells_(),
495  faces_()
496 {}
497 
498 
500 (
501  const word& name,
502  const polyMesh& mesh,
503  const meshSearch& searchEngine,
504  const dictionary& dict
505 )
506 :
507  coordSet(name, dict.get<word>("axis")),
508  mesh_(mesh),
509  searchEngine_(searchEngine),
510  segments_(),
511  cells_(),
512  faces_()
513 {}
514 
515 
516 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
517 
519 (
520  const word& name,
521  const polyMesh& mesh,
522  const meshSearch& searchEngine,
523  const dictionary& dict
524 )
525 {
526  const word sampleType(dict.get<word>("type"));
527 
528  auto* ctorPtr = wordConstructorTable(sampleType);
529 
530  if (!ctorPtr)
531  {
533  (
534  dict,
535  "sample",
536  sampleType,
537  *wordConstructorTablePtr_
538  ) << exit(FatalIOError);
539  }
540 
541  return autoPtr<sampledSet>
542  (
543  ctorPtr
544  (
545  name,
546  mesh,
547  searchEngine,
548  dict.optionalSubDict(sampleType + "Coeffs")
549  )
550  );
551 }
552 
553 
555 {
557 
558  os << nl << "\t(celli)\t(facei)" << nl;
559 
560  forAll(*this, samplei)
561  {
562  os << '\t' << cells_[samplei]
563  << '\t' << faces_[samplei]
564  << nl;
565  }
566 
567  return os;
568 }
569 
570 
571 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::sampledSet::New
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition: sampledSet.C:519
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::sampledSet::checkDimensions
void checkDimensions() const
Check for consistent sizing.
Definition: sampledSet.C:48
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::sampledSet::write
Ostream & write(Ostream &) const
Output for debugging.
Definition: sampledSet.C:554
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::sampledSet::gather
autoPtr< coordSet > gather(labelList &indexSet, labelList &allSegments) const
Helper: gather onto master and sort.
Definition: sampledSet.C:422
Foam::sampledSet::getNeighbourCell
label getNeighbourCell(const label) const
Returns the neighbour cell or the owner if face in on the boundary.
Definition: sampledSet.C:76
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::coordSet::coordFormat
coordFormat
Enumeration defining the output format for coordinates.
Definition: coordSet.H:62
globalIndex.H
Foam::sampledSet::calcSign
scalar calcSign(const label facei, const point &sample) const
Calculates inproduct of face normal and vector sample-face centre.
Definition: sampledSet.C:162
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::sampledSet::findNearFace
label findNearFace(const label celli, const point &sample, const scalar smallDist) const
Returns face label (or -1) of face which is close to sample.
Definition: sampledSet.C:186
Foam::PointHit::hitPoint
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
Definition: PointHit.H:145
Foam::FatalIOError
IOerror FatalIOError
primitiveMesh.H
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::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
polyMesh.H
Foam::sampledSet::sampledSet
sampledSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const coordSet::coordFormat axisType)
Construct from components.
Definition: sampledSet.C:466
Foam::polyMesh::pointInCell
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1397
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
n
label n
Definition: TABSMDCalcMethod2.H:31
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
sampledSet.H
Foam::PointHit::missPoint
const point_type & missPoint() const
Return the miss point. Fatal if hit.
Definition: PointHit.H:158
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::sampledSet::segments_
labelList segments_
Segment numbers.
Definition: sampledSet.H:99
Foam::globalIndex::gatherOp
static void gatherOp(const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const Pstream::commsTypes=Pstream::commsTypes::nonBlocking)
Collect data in processor order on master.
Definition: globalIndexTemplates.C:472
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::sampledSet::getBoundaryCell
label getBoundaryCell(const label) const
Returns cell next to boundary face.
Definition: sampledSet.C:70
Foam::sampledSet::pushIn
point pushIn(const point &sample, const label facei) const
Moves sample in direction of -n to it is 'inside' of facei.
Definition: sampledSet.C:221
Foam::sampledSet::cells_
labelList cells_
Cell numbers.
Definition: sampledSet.H:102
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::coordSet::curveDist_
scalarList curveDist_
Cumulative distance "distance" write specifier.
Definition: coordSet.H:89
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1507
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
meshSearch.H
Foam::sampledSet::setSamples
void setSamples(const List< point > &samplingPts, const labelList &samplingCells, const labelList &samplingFaces, const labelList &samplingSegments, const scalarList &samplingCurveDist)
Set sample data. Copy list contents.
Definition: sampledSet.C:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::coordSet::write
Ostream & write(Ostream &os) const
Write to stream.
Definition: coordSet.C:186
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::sampledSet::getTrackingPoint
bool getTrackingPoint(const point &samplePt, const point &bPoint, const label bFacei, const scalar smallDist, point &trackPt, label &trackCelli, label &trackFacei) const
Calculates start of tracking given samplePt and first boundary.
Definition: sampledSet.C:286
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::dictionary::optionalSubDict
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:577
particle.H
Foam::sampledSet::faces_
labelList faces_
Face numbers (-1 if not known)
Definition: sampledSet.H:105
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
writer.H
Foam::polyMesh::findTetFacePt
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Find the tetFacei and tetPti for point p in celli.
Definition: polyMesh.C:1381
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::PointHit::hit
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::sampledSet::pointInCell
label pointInCell(const point &p, const label samplei) const
Return the cell in which the point on the sample line.
Definition: sampledSet.C:90
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
sample
Minimal example by using system/controlDict.functions: