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-2016 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  const scalar uniformSet::tol = 1e-3;
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 bool Foam::uniformSet::nextSample
47 (
48  const point& currentPt,
49  const vector& offset,
50  const scalar smallDist,
51  point& samplePt,
52  label& sampleI
53 ) const
54 {
55  bool pointFound = false;
56 
57  const vector normOffset = offset/mag(offset);
58 
59  samplePt += offset;
60  sampleI++;
61 
62  for (; sampleI < nPoints_; sampleI++)
63  {
64  scalar s = (samplePt - currentPt) & normOffset;
65 
66  if (s > -smallDist)
67  {
68  // samplePt is close to or beyond currentPt -> use it
69  pointFound = true;
70 
71  break;
72  }
73  samplePt += offset;
74  }
75 
76  return pointFound;
77 }
78 
79 
80 bool Foam::uniformSet::trackToBoundary
81 (
82  passiveParticleCloud& particleCloud,
83  passiveParticle& singleParticle,
84  point& samplePt,
85  label& sampleI,
86  DynamicList<point>& samplingPts,
87  DynamicList<label>& samplingCells,
88  DynamicList<label>& samplingFaces,
89  DynamicList<scalar>& samplingCurveDist
90 ) const
91 {
92  // distance vector between sampling points
93  const vector offset = (end_ - start_)/(nPoints_ - 1);
94  const vector smallVec = tol*offset;
95  const scalar smallDist = mag(smallVec);
96 
97  // Alias
98  const point& trackPt = singleParticle.position();
99 
100  particle::TrackingData<passiveParticleCloud> trackData(particleCloud);
101 
102  while(true)
103  {
104  // Find next samplePt on/after trackPt. Update samplePt, sampleI
105  if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
106  {
107  // no more samples.
108  if (debug)
109  {
110  Pout<< "trackToBoundary : Reached end : samplePt now:"
111  << samplePt << " sampleI now:" << sampleI << endl;
112  }
113  return false;
114  }
115 
116  if (mag(samplePt - trackPt) < smallDist)
117  {
118  // trackPt corresponds with samplePt. Store and use next samplePt
119  if (debug)
120  {
121  Pout<< "trackToBoundary : samplePt corresponds to trackPt : "
122  << " trackPt:" << trackPt << " samplePt:" << samplePt
123  << endl;
124  }
125 
126  samplingPts.append(trackPt);
127  samplingCells.append(singleParticle.cell());
128  samplingFaces.append(-1);
129  samplingCurveDist.append(mag(trackPt - start_));
130 
131  // go to next samplePt
132  if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
133  {
134  // no more samples.
135  if (debug)
136  {
137  Pout<< "trackToBoundary : Reached end : "
138  << " samplePt now:" << samplePt
139  << " sampleI now:" << sampleI
140  << endl;
141  }
142 
143  return false;
144  }
145  }
146 
147 
148  if (debug)
149  {
150  Pout<< "Searching along trajectory from "
151  << " trackPt:" << trackPt
152  << " trackCelli:" << singleParticle.cell()
153  << " to:" << samplePt << endl;
154  }
155 
156  point oldPos = trackPt;
157  label facei = -1;
158  do
159  {
160  singleParticle.stepFraction() = 0;
161  singleParticle.track(samplePt, trackData);
162 
163  if (debug)
164  {
165  Pout<< "Result of tracking "
166  << " trackPt:" << trackPt
167  << " trackCelli:" << singleParticle.cell()
168  << " trackFacei:" << singleParticle.face()
169  << " onBoundary:" << singleParticle.onBoundary()
170  << " samplePt:" << samplePt
171  << " smallDist:" << smallDist
172  << endl;
173  }
174  }
175  while
176  (
177  !singleParticle.onBoundary()
178  && (mag(trackPt - oldPos) < smallDist)
179  );
180 
181  if (singleParticle.onBoundary())
182  {
183  //Pout<< "trackToBoundary : reached boundary" << endl;
184  if (mag(trackPt - samplePt) < smallDist)
185  {
186  //Pout<< "trackToBoundary : boundary is also sampling point"
187  // << endl;
188  // Reached samplePt on boundary
189  samplingPts.append(trackPt);
190  samplingCells.append(singleParticle.cell());
191  samplingFaces.append(facei);
192  samplingCurveDist.append(mag(trackPt - start_));
193  }
194 
195  return true;
196  }
197 
198  //Pout<< "trackToBoundary : reached internal sampling point" << endl;
199  // Reached samplePt in cell or on internal face
200  samplingPts.append(trackPt);
201  samplingCells.append(singleParticle.cell());
202  samplingFaces.append(-1);
203  samplingCurveDist.append(mag(trackPt - start_));
204 
205  // go to next samplePt
206  }
207 }
208 
209 
210 void Foam::uniformSet::calcSamples
211 (
212  DynamicList<point>& samplingPts,
213  DynamicList<label>& samplingCells,
214  DynamicList<label>& samplingFaces,
215  DynamicList<label>& samplingSegments,
216  DynamicList<scalar>& samplingCurveDist
217 ) const
218 {
219  // distance vector between sampling points
220  if ((nPoints_ < 2) || (mag(end_ - start_) < SMALL))
221  {
223  << "Incorrect sample specification. Either too few points or"
224  << " start equals end point." << endl
225  << "nPoints:" << nPoints_
226  << " start:" << start_
227  << " end:" << end_
228  << exit(FatalError);
229  }
230 
231  const vector offset = (end_ - start_)/(nPoints_ - 1);
232  const vector normOffset = offset/mag(offset);
233  const vector smallVec = tol*offset;
234  const scalar smallDist = mag(smallVec);
235 
236  // Force calculation of cloud addressing on all processors
237  const bool oldMoving = const_cast<polyMesh&>(mesh()).moving(false);
238  passiveParticleCloud particleCloud(mesh());
239 
240  // Get all boundary intersections
241  List<pointIndexHit> bHits = searchEngine().intersections
242  (
243  start_ - smallVec,
244  end_ + smallVec
245  );
246 
247  point bPoint(GREAT, GREAT, GREAT);
248  label bFacei = -1;
249 
250  if (bHits.size())
251  {
252  bPoint = bHits[0].hitPoint();
253  bFacei = bHits[0].index();
254  }
255 
256  // Get first tracking point. Use bPoint, bFacei if provided.
257 
258  point trackPt;
259  label trackCelli = -1;
260  label trackFacei = -1;
261 
262  bool isSample =
263  getTrackingPoint
264  (
265  start_,
266  bPoint,
267  bFacei,
268  smallDist,
269 
270  trackPt,
271  trackCelli,
272  trackFacei
273  );
274 
275  if (trackCelli == -1)
276  {
277  // Line start_ - end_ does not intersect domain at all.
278  // (or is along edge)
279  // Set points and cell/face labels to empty lists
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 
394 void Foam::uniformSet::genSamples()
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 
403  calcSamples
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 // ************************************************************************* //
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:891
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual ~uniformSet()
Destructor.
Definition: uniformSet.C:480
Macros for easy insertion into run-time selection tables.
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))
dynamicFvMesh & mesh
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
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
A class for handling words, derived from string.
Definition: word.H:59
label readLabel(Istream &is)
Definition: label.H:64
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
vector point
Point is a vector.
Definition: point.H:41
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static const scalar tol
Tolerance when comparing points relative to difference between.
Definition: uniformSet.H:119
runTime write()
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451