uniformSet.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-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 "uniformSet.H"
27 #include "meshSearch.H"
28 #include "DynamicList.H"
29 #include "polyMesh.H"
30 
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(uniformSet, 0);
38  addToRunTimeSelectionTable(sampledSet, uniformSet, word);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 (
46  const point& currentPt,
47  const vector& offset,
48  const scalar smallDist,
49  point& samplePt,
50  label& sampleI
51 ) const
52 {
53  bool pointFound = false;
54 
55  const vector normOffset = offset/mag(offset);
56 
57  samplePt += offset;
58  sampleI++;
59 
60  for (; sampleI < nPoints_; sampleI++)
61  {
62  scalar s = (samplePt - currentPt) & normOffset;
63 
64  if (s > -smallDist)
65  {
66  // samplePt is close to or beyond currentPt -> use it
67  pointFound = true;
68 
69  break;
70  }
71  samplePt += offset;
72  }
73 
74  return pointFound;
75 }
76 
77 
79 (
80  passiveParticleCloud& particleCloud,
81  passiveParticle& singleParticle,
82  point& samplePt,
83  label& sampleI,
84  DynamicList<point>& samplingPts,
85  DynamicList<label>& samplingCells,
86  DynamicList<label>& samplingFaces,
87  DynamicList<scalar>& samplingCurveDist
88 ) const
89 {
90  // distance vector between sampling points
91  const vector offset = (end_ - start_)/(nPoints_ - 1);
92  const vector smallVec = tol*offset;
93  const scalar smallDist = mag(smallVec);
94 
95  // Alias
96  const point& trackPt = singleParticle.position();
97 
98  particle::TrackingData<passiveParticleCloud> trackData(particleCloud);
99 
100  while(true)
101  {
102  // Find next samplePt on/after trackPt. Update samplePt, sampleI
103  if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
104  {
105  // no more samples.
106  if (debug)
107  {
108  Pout<< "trackToBoundary : Reached end : samplePt now:"
109  << samplePt << " sampleI now:" << sampleI << endl;
110  }
111  return false;
112  }
113 
114  if (mag(samplePt - trackPt) < smallDist)
115  {
116  // trackPt corresponds with samplePt. Store and use next samplePt
117  if (debug)
118  {
119  Pout<< "trackToBoundary : samplePt corresponds to trackPt : "
120  << " trackPt:" << trackPt << " samplePt:" << samplePt
121  << endl;
122  }
123 
124  samplingPts.append(trackPt);
125  samplingCells.append(singleParticle.cell());
126  samplingFaces.append(-1);
127  samplingCurveDist.append(mag(trackPt - start_));
128 
129  // go to next samplePt
130  if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
131  {
132  // no more samples.
133  if (debug)
134  {
135  Pout<< "trackToBoundary : Reached end : "
136  << " samplePt now:" << samplePt
137  << " sampleI now:" << sampleI
138  << endl;
139  }
140 
141  return false;
142  }
143  }
144 
145 
146  if (debug)
147  {
148  Pout<< "Searching along trajectory from "
149  << " trackPt:" << trackPt
150  << " trackCellI:" << singleParticle.cell()
151  << " to:" << samplePt << endl;
152  }
153 
154  point oldPos = trackPt;
155  label facei = -1;
156  do
157  {
158  singleParticle.stepFraction() = 0;
159  singleParticle.track(samplePt, trackData);
160 
161  if (debug)
162  {
163  Pout<< "Result of tracking "
164  << " trackPt:" << trackPt
165  << " trackCellI:" << singleParticle.cell()
166  << " trackFaceI:" << singleParticle.face()
167  << " onBoundary:" << singleParticle.onBoundary()
168  << " samplePt:" << samplePt
169  << " smallDist:" << smallDist
170  << endl;
171  }
172  }
173  while
174  (
175  !singleParticle.onBoundary()
176  && (mag(trackPt - oldPos) < smallDist)
177  );
178 
179  if (singleParticle.onBoundary())
180  {
181  //Pout<< "trackToBoundary : reached boundary" << endl;
182  if (mag(trackPt - samplePt) < smallDist)
183  {
184  //Pout<< "trackToBoundary : boundary is also sampling point"
185  // << endl;
186  // Reached samplePt on boundary
187  samplingPts.append(trackPt);
188  samplingCells.append(singleParticle.cell());
189  samplingFaces.append(facei);
190  samplingCurveDist.append(mag(trackPt - start_));
191  }
192 
193  return true;
194  }
195 
196  //Pout<< "trackToBoundary : reached internal sampling point" << endl;
197  // Reached samplePt in cell or on internal face
198  samplingPts.append(trackPt);
199  samplingCells.append(singleParticle.cell());
200  samplingFaces.append(-1);
201  samplingCurveDist.append(mag(trackPt - start_));
202 
203  // go to next samplePt
204  }
205 }
206 
207 
209 (
210  DynamicList<point>& samplingPts,
211  DynamicList<label>& samplingCells,
212  DynamicList<label>& samplingFaces,
213  DynamicList<label>& samplingSegments,
214  DynamicList<scalar>& samplingCurveDist
215 ) const
216 {
217  // distance vector between sampling points
218  if ((nPoints_ < 2) || (mag(end_ - start_) < SMALL))
219  {
221  << "Incorrect sample specification. Either too few points or"
222  << " start equals end point." << endl
223  << "nPoints:" << nPoints_
224  << " start:" << start_
225  << " end:" << end_
226  << exit(FatalError);
227  }
228 
229  const vector offset = (end_ - start_)/(nPoints_ - 1);
230  const vector normOffset = offset/mag(offset);
231  const vector smallVec = tol*offset;
232  const scalar smallDist = mag(smallVec);
233 
234  // Force calculation of cloud addressing on all processors
235  const bool oldMoving = const_cast<polyMesh&>(mesh()).moving(false);
236  passiveParticleCloud particleCloud(mesh());
237 
238  // Get all boundary intersections
239  List<pointIndexHit> bHits = searchEngine().intersections
240  (
241  start_ - smallVec,
242  end_ + smallVec
243  );
244 
245  point bPoint(GREAT, GREAT, GREAT);
246  label bFaceI = -1;
247 
248  if (bHits.size())
249  {
250  bPoint = bHits[0].hitPoint();
251  bFaceI = bHits[0].index();
252  }
253 
254  // Get first tracking point. Use bPoint, bFaceI if provided.
255 
256  point trackPt;
257  label trackCellI = -1;
258  label trackFaceI = -1;
259 
260  bool isSample =
261  getTrackingPoint
262  (
263  offset,
264  start_,
265  bPoint,
266  bFaceI,
267 
268  trackPt,
269  trackCellI,
270  trackFaceI
271  );
272 
273  if (trackCellI == -1)
274  {
275  // Line start_ - end_ does not intersect domain at all.
276  // (or is along edge)
277  // Set points and cell/face labels to empty lists
278 
279  const_cast<polyMesh&>(mesh()).moving(oldMoving);
280 
281  return;
282  }
283 
284  if (isSample)
285  {
286  samplingPts.append(start_);
287  samplingCells.append(trackCellI);
288  samplingFaces.append(trackFaceI);
289  samplingCurveDist.append(0.0);
290  }
291 
292  //
293  // Track until hit end of all boundary intersections
294  //
295 
296  // current segment number
297  label segmentI = 0;
298 
299  // starting index of current segment in samplePts
300  label startSegmentI = 0;
301 
302  label sampleI = 0;
303  point samplePt = start_;
304 
305  // index in bHits; current boundary intersection
306  label bHitI = 1;
307 
308  while(true)
309  {
310  // Initialize tracking starting from trackPt
311  passiveParticle singleParticle(mesh(), trackPt, trackCellI);
312 
313  bool reachedBoundary = trackToBoundary
314  (
315  particleCloud,
316  singleParticle,
317  samplePt,
318  sampleI,
319  samplingPts,
320  samplingCells,
321  samplingFaces,
322  samplingCurveDist
323  );
324 
325  // fill sampleSegments
326  for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
327  {
328  samplingSegments.append(segmentI);
329  }
330 
331 
332  if (!reachedBoundary)
333  {
334  if (debug)
335  {
336  Pout<< "calcSamples : Reached end of samples: "
337  << " samplePt now:" << samplePt
338  << " sampleI now:" << sampleI
339  << endl;
340  }
341  break;
342  }
343 
344 
345  bool foundValidB = false;
346 
347  while (bHitI < bHits.size())
348  {
349  scalar dist =
350  (bHits[bHitI].hitPoint() - singleParticle.position())
351  & normOffset;
352 
353  if (debug)
354  {
355  Pout<< "Finding next boundary : "
356  << "bPoint:" << bHits[bHitI].hitPoint()
357  << " tracking:" << singleParticle.position()
358  << " dist:" << dist
359  << endl;
360  }
361 
362  if (dist > smallDist)
363  {
364  // hitpoint is past tracking position
365  foundValidB = true;
366  break;
367  }
368  else
369  {
370  bHitI++;
371  }
372  }
373 
374  if (!foundValidB)
375  {
376  // No valid boundary intersection found beyond tracking position
377  break;
378  }
379 
380  // Update starting point for tracking
381  trackFaceI = bFaceI;
382  trackPt = pushIn(bPoint, trackFaceI);
383  trackCellI = getBoundaryCell(trackFaceI);
384 
385  segmentI++;
386 
387  startSegmentI = samplingPts.size();
388  }
389 
390  const_cast<polyMesh&>(mesh()).moving(oldMoving);
391 }
392 
393 
395 {
396  // Storage for sample points
397  DynamicList<point> samplingPts;
398  DynamicList<label> samplingCells;
399  DynamicList<label> samplingFaces;
400  DynamicList<label> samplingSegments;
401  DynamicList<scalar> samplingCurveDist;
402 
404  (
405  samplingPts,
406  samplingCells,
407  samplingFaces,
408  samplingSegments,
409  samplingCurveDist
410  );
411 
412  samplingPts.shrink();
413  samplingCells.shrink();
414  samplingFaces.shrink();
415  samplingSegments.shrink();
416  samplingCurveDist.shrink();
417 
418  setSamples
419  (
420  samplingPts,
421  samplingCells,
422  samplingFaces,
423  samplingSegments,
424  samplingCurveDist
425  );
426 }
427 
428 
429 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
430 
432 (
433  const word& name,
434  const polyMesh& mesh,
435  const meshSearch& searchEngine,
436  const word& axis,
437  const point& start,
438  const point& end,
439  const label nPoints
440 )
441 :
442  sampledSet(name, mesh, searchEngine, axis),
443  start_(start),
444  end_(end),
445  nPoints_(nPoints)
446 {
447  genSamples();
448 
449  if (debug)
450  {
451  write(Pout);
452  }
453 }
454 
455 
457 (
458  const word& name,
459  const polyMesh& mesh,
460  const meshSearch& searchEngine,
461  const dictionary& dict
462 )
463 :
464  sampledSet(name, mesh, searchEngine, dict),
465  start_(dict.lookup("start")),
466  end_(dict.lookup("end")),
467  nPoints_(readLabel(dict.lookup("nPoints")))
468 {
469  genSamples();
470 
471  if (debug)
472  {
473  write(Pout);
474  }
475 }
476 
477 
478 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
479 
481 {}
482 
483 
484 // ************************************************************************* //
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
Foam::particle::onBoundary
bool onBoundary() const
Is the particle on the boundary/(or outside the domain)?
Definition: particleI.H:810
Foam::particle::face
label & face()
Return current face particle is on otherwise -1.
Definition: particleI.H:658
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::particle::stepFraction
scalar & stepFraction()
Return the fraction of time-step completed.
Definition: particleI.H:816
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::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::passiveParticle
Copy of base particle.
Definition: passiveParticle.H:50
Foam::uniformSet::calcSamples
void calcSamples(DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< label > &samplingSegments, DynamicList< scalar > &samplingCurveDist) const
Samples from start_ to end_. samplingSegments contains segmentNo.
Definition: uniformSet.C:209
Foam::uniformSet::~uniformSet
virtual ~uniformSet()
Destructor.
Definition: uniformSet.C:480
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
polyMesh.H
Foam::meshSearch::intersections
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:899
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::particle::cell
label & cell()
Return current cell particle is in.
Definition: particleI.H:598
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::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Foam::uniformSet::trackToBoundary
bool trackToBoundary(passiveParticleCloud &particleCloud, passiveParticle &singleParticle, point &samplePt, label &sampleI, DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< scalar > &samplingCurveDist) const
Samples from startTrackPt/CellI. Updates particle/samplePt/sampleI.
Definition: uniformSet.C:79
Foam::passiveParticleCloud
A Cloud of passive particles.
Definition: passiveParticleCloud.H:49
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::uniformSet::genSamples
void genSamples()
Uses calcSamples to obtain samples. Copies them into *this.
Definition: uniformSet.C:394
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
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
meshSearch.H
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::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::particle::track
label track(const vector &endPosition, TrackData &td)
Track particle to end of trajectory.
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::uniformSet::uniformSet
uniformSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const point &start, const point &end, const label nPoints)
Construct from components.
Definition: uniformSet.C:432
Foam::particle::position
const vector & position() const
Return current particle position.
Definition: particleI.H:586
DynamicList.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::uniformSet::nextSample
bool nextSample(const point &currentPt, const vector &offset, const scalar smallDist, point &samplePt, label &sampleI) const
Calculates - starting at samplePt - the first sampling point.
Definition: uniformSet.C:45
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::particle::TrackingData
Definition: particle.H:94
write
Tcoeff write()
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
uniformSet.H