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 | Copyright (C) 2011-2012 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 "sampledSet.H"
27 #include "polyMesh.H"
28 #include "primitiveMesh.H"
29 #include "meshSearch.H"
30 #include "writer.H"
31 #include "particle.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  const scalar sampledSet::tol = 1e-6;
38 
39  defineTypeNameAndDebug(sampledSet, 0);
40  defineRunTimeSelectionTable(sampledSet, word);
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 {
48  return mesh().faceOwner()[faceI];
49 }
50 
51 
53 (
54  const label faceI,
55  const point& sample
56 ) const
57 {
58  if (faceI == -1)
59  {
61  (
62  "sampledSet::getCell(const label, const point&)"
63  ) << "Illegal face label " << faceI
64  << abort(FatalError);
65  }
66 
67  if (faceI >= mesh().nInternalFaces())
68  {
69  label cellI = getBoundaryCell(faceI);
70 
71  if (!mesh().pointInCell(sample, cellI, searchEngine_.decompMode()))
72  {
74  (
75  "sampledSet::getCell(const label, const point&)"
76  ) << "Found cell " << cellI << " using face " << faceI
77  << ". But cell does not contain point " << sample
78  << abort(FatalError);
79  }
80  return cellI;
81  }
82  else
83  {
84  // Try owner and neighbour to see which one contains sample
85 
86  label cellI = mesh().faceOwner()[faceI];
87 
88  if (mesh().pointInCell(sample, cellI, searchEngine_.decompMode()))
89  {
90  return cellI;
91  }
92  else
93  {
94  cellI = mesh().faceNeighbour()[faceI];
95 
96  if (mesh().pointInCell(sample, cellI, searchEngine_.decompMode()))
97  {
98  return cellI;
99  }
100  else
101  {
103  (
104  "sampledSet::getCell(const label, const point&)"
105  ) << "None of the neighbours of face "
106  << faceI << " contains point " << sample
107  << abort(FatalError);
108 
109  return -1;
110  }
111  }
112  }
113 }
114 
115 
116 Foam::scalar Foam::sampledSet::calcSign
117 (
118  const label faceI,
119  const point& sample
120 ) const
121 {
122  vector vec = sample - mesh().faceCentres()[faceI];
123 
124  scalar magVec = mag(vec);
125 
126  if (magVec < VSMALL)
127  {
128  // sample on face centre. Regard as inside
129  return -1;
130  }
131 
132  vec /= magVec;
133 
134  vector n = mesh().faceAreas()[faceI];
135 
136  n /= mag(n) + VSMALL;
137 
138  return n & vec;
139 }
140 
141 
142 // Return face (or -1) of face which is within smallDist of sample
144 (
145  const label cellI,
146  const point& sample,
147  const scalar smallDist
148 ) const
149 {
150  const cell& myFaces = mesh().cells()[cellI];
151 
152  forAll(myFaces, myFaceI)
153  {
154  const face& f = mesh().faces()[myFaces[myFaceI]];
155 
156  pointHit inter = f.nearestPoint(sample, mesh().points());
157 
158  scalar dist;
159 
160  if (inter.hit())
161  {
162  dist = mag(inter.hitPoint() - sample);
163  }
164  else
165  {
166  dist = mag(inter.missPoint() - sample);
167  }
168 
169  if (dist < smallDist)
170  {
171  return myFaces[myFaceI];
172  }
173  }
174  return -1;
175 }
176 
177 
178 // 'Pushes' point facePt (which is almost on face) in direction of cell centre
179 // so it is clearly inside.
181 (
182  const point& facePt,
183  const label faceI
184 ) const
185 {
186  label cellI = mesh().faceOwner()[faceI];
187  const point& cC = mesh().cellCentres()[cellI];
188 
189  point newPosition = facePt;
190 
191  // Taken from particle::initCellFacePt()
192  label tetFaceI;
193  label tetPtI;
194  mesh().findTetFacePt(cellI, facePt, tetFaceI, tetPtI);
195 
196  if (tetFaceI == -1 || tetPtI == -1)
197  {
198  newPosition = facePt;
199 
201 
202  label iterNo = 0;
203 
204  do
205  {
206  newPosition += particle::trackingCorrectionTol*(cC - facePt);
207 
209  (
210  cellI,
211  newPosition,
212  tetFaceI,
213  tetPtI
214  );
215 
216  iterNo++;
217 
218  } while (tetFaceI < 0 && iterNo <= trap);
219  }
220 
221  if (tetFaceI == -1)
222  {
224  (
225  "sampledSet::pushIn(const point&, const label)"
226  ) << "After pushing " << facePt << " to " << newPosition
227  << " it is still outside face " << faceI
228  << " at " << mesh().faceCentres()[faceI]
229  << " of cell " << cellI
230  << " at " << cC << endl
231  << "Please change your starting point"
232  << abort(FatalError);
233  }
234 
235  //Info<< "pushIn : moved " << facePt << " to " << newPosition
236  // << endl;
237 
238  return newPosition;
239 }
240 
241 
242 // Calculates start of tracking given samplePt and first boundary intersection
243 // (bPoint, bFaceI). bFaceI == -1 if no boundary intersection.
244 // Returns true if trackPt is sampling point
246 (
247  const vector& offset,
248  const point& samplePt,
249  const point& bPoint,
250  const label bFaceI,
251 
252  point& trackPt,
253  label& trackCellI,
254  label& trackFaceI
255 ) const
256 {
257  const scalar smallDist = mag(tol*offset);
258 
259  bool isGoodSample = false;
260 
261  if (bFaceI == -1)
262  {
263  // No boundary intersection. Try and find cell samplePt is in
264  trackCellI = mesh().findCell(samplePt, searchEngine_.decompMode());
265 
266  if
267  (
268  (trackCellI == -1)
269  || !mesh().pointInCell
270  (
271  samplePt,
272  trackCellI,
273  searchEngine_.decompMode()
274  )
275  )
276  {
277  // Line samplePt - end_ does not intersect domain at all.
278  // (or is along edge)
279  //Info<< "getTrackingPoint : samplePt outside domain : "
280  // << " samplePt:" << samplePt
281  // << endl;
282 
283  trackCellI = -1;
284  trackFaceI = -1;
285 
286  isGoodSample = false;
287  }
288  else
289  {
290  // start is inside. Use it as tracking point
291  //Info<< "getTrackingPoint : samplePt inside :"
292  // << " samplePt:" << samplePt
293  // << " trackCellI:" << trackCellI
294  // << endl;
295 
296  trackPt = samplePt;
297  trackFaceI = -1;
298 
299  isGoodSample = true;
300  }
301  }
302  else if (mag(samplePt - bPoint) < smallDist)
303  {
304  //Info<< "getTrackingPoint : samplePt:" << samplePt
305  // << " close to bPoint:"
306  // << bPoint << endl;
307 
308  // samplePt close to bPoint. Snap to it
309  trackPt = pushIn(bPoint, bFaceI);
310  trackFaceI = bFaceI;
311  trackCellI = getBoundaryCell(trackFaceI);
312 
313  isGoodSample = true;
314  }
315  else
316  {
317  scalar sign = calcSign(bFaceI, samplePt);
318 
319  if (sign < 0)
320  {
321  // samplePt inside or marginally outside.
322  trackPt = samplePt;
323  trackFaceI = -1;
324  trackCellI = mesh().findCell(trackPt, searchEngine_.decompMode());
325 
326  isGoodSample = true;
327  }
328  else
329  {
330  // samplePt outside. use bPoint
331  trackPt = pushIn(bPoint, bFaceI);
332  trackFaceI = bFaceI;
333  trackCellI = getBoundaryCell(trackFaceI);
334 
335  isGoodSample = false;
336  }
337  }
338 
339  if (debug)
340  {
341  Info<< "sampledSet::getTrackingPoint :"
342  << " offset:" << offset
343  << " samplePt:" << samplePt
344  << " bPoint:" << bPoint
345  << " bFaceI:" << bFaceI
346  << endl << " Calculated first tracking point :"
347  << " trackPt:" << trackPt
348  << " trackCellI:" << trackCellI
349  << " trackFaceI:" << trackFaceI
350  << " isGoodSample:" << isGoodSample
351  << endl;
352  }
353 
354  return isGoodSample;
355 }
356 
357 
359 (
360  const List<point>& samplingPts,
361  const labelList& samplingCells,
362  const labelList& samplingFaces,
363  const labelList& samplingSegments,
364  const scalarList& samplingCurveDist
365 )
366 {
367  setSize(samplingPts.size());
368  cells_.setSize(samplingCells.size());
369  faces_.setSize(samplingFaces.size());
370  segments_.setSize(samplingSegments.size());
371  curveDist_.setSize(samplingCurveDist.size());
372 
373  if
374  (
375  (cells_.size() != size())
376  || (faces_.size() != size())
377  || (segments_.size() != size())
378  || (curveDist_.size() != size())
379  )
380  {
381  FatalErrorIn("sampledSet::setSamples()")
382  << "sizes not equal : "
383  << " points:" << size()
384  << " cells:" << cells_.size()
385  << " faces:" << faces_.size()
386  << " segments:" << segments_.size()
387  << " curveDist:" << curveDist_.size()
388  << abort(FatalError);
389  }
390 
391  forAll(samplingPts, sampleI)
392  {
393  operator[](sampleI) = samplingPts[sampleI];
394  }
395  curveDist_ = samplingCurveDist;
396 
397  cells_ = samplingCells;
398  faces_ = samplingFaces;
399  segments_ = samplingSegments;
400 }
401 
402 
403 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
404 
406 (
407  const word& name,
408  const polyMesh& mesh,
409  const meshSearch& searchEngine,
410  const word& axis
411 )
412 :
413  coordSet(name, axis),
414  mesh_(mesh),
415  searchEngine_(searchEngine),
416  segments_(0),
417  cells_(0),
418  faces_(0)
419 {}
420 
421 
423 (
424  const word& name,
425  const polyMesh& mesh,
426  const meshSearch& searchEngine,
427  const dictionary& dict
428 )
429 :
430  coordSet(name, dict.lookup("axis")),
431  mesh_(mesh),
432  searchEngine_(searchEngine),
433  segments_(0),
434  cells_(0),
435  faces_(0)
436 {}
437 
438 
439 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
440 
442 {}
443 
444 
445 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
446 
448 (
449  const word& name,
450  const polyMesh& mesh,
451  const meshSearch& searchEngine,
452  const dictionary& dict
453 )
454 {
455  const word sampleType(dict.lookup("type"));
456 
457  wordConstructorTable::iterator cstrIter =
458  wordConstructorTablePtr_->find(sampleType);
459 
460  if (cstrIter == wordConstructorTablePtr_->end())
461  {
463  (
464  "sampledSet::New"
465  "(const word&, const polyMesh&, const meshSearch&"
466  ", const dictionary&)"
467  ) << "Unknown sample type "
468  << sampleType << nl << nl
469  << "Valid sample types : " << endl
470  << wordConstructorTablePtr_->sortedToc()
471  << exit(FatalError);
472  }
473 
474  return autoPtr<sampledSet>
475  (
476  cstrIter()
477  (
478  name,
479  mesh,
480  searchEngine,
481  dict
482  )
483  );
484 }
485 
486 
488 {
489  coordSet::write(os);
490 
491  os << endl << "\t(cellI)\t(faceI)" << endl;
492 
493  forAll(*this, sampleI)
494  {
495  os << '\t' << cells_[sampleI]
496  << '\t' << faces_[sampleI]
497  << endl;
498  }
499 
500  return os;
501 }
502 
503 
504 // ************************************************************************* //
const vectorField & faceAreas() const
virtual ~sampledSet()
Destructor.
Definition: sampledSet.C:441
const pointField & points
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:1341
dimensioned< scalar > mag(const dimensioned< Type > &)
label getCell(const label faceI, const point &sample) const
Returns cell using face and containing sample.
Definition: sampledSet.C:53
Holds list of sampling positions.
Definition: coordSet.H:49
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
labelList f(nPoints)
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:448
const Point & missPoint() const
Return miss point.
Definition: PointHit.H:145
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:144
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
A class for handling words, derived from string.
Definition: word.H:59
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
static const scalar tol
Tolerance when comparing points. Usually relative to difference.
Definition: sampledSet.H:199
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const cellList & cells() const
dimensionedScalar sign(const dimensionedScalar &ds)
const vectorField & cellCentres() const
messageStream Info
point pushIn(const point &sample, const label faceI) const
Moves sample in direction of -n to it is &#39;inside&#39; of faceI.
Definition: sampledSet.C:181
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
label n
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
dictionary dict
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static const scalar trackingCorrectionTol
Fraction of distance to tet centre to move a particle to.
Definition: particle.H:322
label getBoundaryCell(const label) const
Returns cell next to boundary face.
Definition: sampledSet.C:46
#define forAll(list, i)
Definition: UList.H:421
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1476
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1073
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
points setSize(newPointi)
Ostream & write(Ostream &os) const
Definition: coordSet.C:138
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
bool hit() const
Is there a hit.
Definition: PointHit.H:120
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
Ostream & write(Ostream &) const
Output for debugging.
Definition: sampledSet.C:487
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1079
sampledSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis)
Construct from components.
Definition: sampledSet.C:406
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1357
scalar calcSign(const label faceI, const point &sample) const
Calculates inproduct of face normal and vector sample-face centre.
Definition: sampledSet.C:117
bool getTrackingPoint(const vector &offset, const point &samplePt, const point &bPoint, const label bFaceI, point &trackPt, label &trackCellI, label &trackFaceI) const
Calculates start of tracking given samplePt and first boundary.
Definition: sampledSet.C:246
void setSamples(const List< point > &samplingPts, const labelList &samplingCells, const labelList &samplingFaces, const labelList &samplingSegments, const scalarList &samplingCurveDist)
Sets sample data.
Definition: sampledSet.C:359
defineTypeNameAndDebug(combustionModel, 0)
const vectorField & faceCentres() const