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-2017 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  point trackPt = singleParticle.position();
98 
99  particle::TrackingData<passiveParticleCloud> trackData(particleCloud);
100 
101  while(true)
102  {
103  // Find next samplePt on/after trackPt. Update samplePt, sampleI
104  if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
105  {
106  // no more samples.
107  if (debug)
108  {
109  Pout<< "trackToBoundary : Reached end : samplePt now:"
110  << samplePt << " sampleI now:" << sampleI << endl;
111  }
112  return false;
113  }
114 
115  if (mag(samplePt - trackPt) < smallDist)
116  {
117  // trackPt corresponds with samplePt. Store and use next samplePt
118  if (debug)
119  {
120  Pout<< "trackToBoundary : samplePt corresponds to trackPt : "
121  << " trackPt:" << trackPt << " samplePt:" << samplePt
122  << endl;
123  }
124 
125  samplingPts.append(trackPt);
126  samplingCells.append(singleParticle.cell());
127  samplingFaces.append(-1);
128  samplingCurveDist.append(mag(trackPt - start_));
129 
130  // go to next samplePt
131  if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
132  {
133  // no more samples.
134  if (debug)
135  {
136  Pout<< "trackToBoundary : Reached end : "
137  << " samplePt now:" << samplePt
138  << " sampleI now:" << sampleI
139  << endl;
140  }
141 
142  return false;
143  }
144  }
145 
146 
147  if (debug)
148  {
149  Pout<< "Searching along trajectory from "
150  << " trackPt:" << trackPt
151  << " trackCelli:" << singleParticle.cell()
152  << " to:" << samplePt << endl;
153  }
154 
155  singleParticle.track(samplePt - trackPt, 0);
156  trackPt = singleParticle.position();
157 
158  if (singleParticle.onBoundaryFace())
159  {
160  //Pout<< "trackToBoundary : reached boundary" << endl;
161  if (mag(trackPt - samplePt) < smallDist)
162  {
163  //Pout<< "trackToBoundary : boundary is also sampling point"
164  // << endl;
165  // Reached samplePt on boundary
166  samplingPts.append(trackPt);
167  samplingCells.append(singleParticle.cell());
168  samplingFaces.append(singleParticle.face());
169  samplingCurveDist.append(mag(trackPt - start_));
170  }
171 
172  return true;
173  }
174 
175  //Pout<< "trackToBoundary : reached internal sampling point" << endl;
176  // Reached samplePt in cell or on internal face
177  samplingPts.append(trackPt);
178  samplingCells.append(singleParticle.cell());
179  samplingFaces.append(-1);
180  samplingCurveDist.append(mag(trackPt - start_));
181 
182  // go to next samplePt
183  }
184 }
185 
186 
187 void Foam::uniformSet::calcSamples
188 (
189  DynamicList<point>& samplingPts,
190  DynamicList<label>& samplingCells,
191  DynamicList<label>& samplingFaces,
192  DynamicList<label>& samplingSegments,
193  DynamicList<scalar>& samplingCurveDist
194 ) const
195 {
196  // distance vector between sampling points
197  if ((nPoints_ < 2) || (mag(end_ - start_) < SMALL))
198  {
200  << "Incorrect sample specification. Either too few points or"
201  << " start equals end point." << endl
202  << "nPoints:" << nPoints_
203  << " start:" << start_
204  << " end:" << end_
205  << exit(FatalError);
206  }
207 
208  const vector offset = (end_ - start_)/(nPoints_ - 1);
209  const vector normOffset = offset/mag(offset);
210  const vector smallVec = tol*offset;
211  const scalar smallDist = mag(smallVec);
212 
213  // Force calculation of cloud addressing on all processors
214  const bool oldMoving = const_cast<polyMesh&>(mesh()).moving(false);
215  passiveParticleCloud particleCloud(mesh());
216 
217  // Get all boundary intersections
218  List<pointIndexHit> bHits = searchEngine().intersections
219  (
220  start_ - smallVec,
221  end_ + smallVec
222  );
223 
224  point bPoint(GREAT, GREAT, GREAT);
225  label bFacei = -1;
226 
227  if (bHits.size())
228  {
229  bPoint = bHits[0].hitPoint();
230  bFacei = bHits[0].index();
231  }
232 
233  // Get first tracking point. Use bPoint, bFacei if provided.
234 
235  point trackPt;
236  label trackCelli = -1;
237  label trackFacei = -1;
238 
239  bool isSample =
240  getTrackingPoint
241  (
242  start_,
243  bPoint,
244  bFacei,
245  smallDist,
246 
247  trackPt,
248  trackCelli,
249  trackFacei
250  );
251 
252  if (trackCelli == -1)
253  {
254  // Line start_ - end_ does not intersect domain at all.
255  // (or is along edge)
256  // Set points and cell/face labels to empty lists
257 
258  const_cast<polyMesh&>(mesh()).moving(oldMoving);
259  return;
260  }
261 
262  if (isSample)
263  {
264  samplingPts.append(start_);
265  samplingCells.append(trackCelli);
266  samplingFaces.append(trackFacei);
267  samplingCurveDist.append(0.0);
268  }
269 
270  //
271  // Track until hit end of all boundary intersections
272  //
273 
274  // current segment number
275  label segmentI = 0;
276 
277  // starting index of current segment in samplePts
278  label startSegmentI = 0;
279 
280  label sampleI = 0;
281  point samplePt = start_;
282 
283  // index in bHits; current boundary intersection
284  label bHitI = 1;
285 
286  while(true)
287  {
288  // Initialize tracking starting from trackPt
289  passiveParticle singleParticle(mesh(), trackPt, trackCelli);
290 
291  bool reachedBoundary = trackToBoundary
292  (
293  particleCloud,
294  singleParticle,
295  samplePt,
296  sampleI,
297  samplingPts,
298  samplingCells,
299  samplingFaces,
300  samplingCurveDist
301  );
302 
303  // fill sampleSegments
304  for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
305  {
306  samplingSegments.append(segmentI);
307  }
308 
309 
310  if (!reachedBoundary)
311  {
312  if (debug)
313  {
314  Pout<< "calcSamples : Reached end of samples: "
315  << " samplePt now:" << samplePt
316  << " sampleI now:" << sampleI
317  << endl;
318  }
319  break;
320  }
321 
322 
323  bool foundValidB = false;
324 
325  while (bHitI < bHits.size())
326  {
327  scalar dist =
328  (bHits[bHitI].hitPoint() - singleParticle.position())
329  & normOffset;
330 
331  if (debug)
332  {
333  Pout<< "Finding next boundary : "
334  << "bPoint:" << bHits[bHitI].hitPoint()
335  << " tracking:" << singleParticle.position()
336  << " dist:" << dist
337  << endl;
338  }
339 
340  if (dist > smallDist)
341  {
342  // hitpoint is past tracking position
343  foundValidB = true;
344  break;
345  }
346  else
347  {
348  bHitI++;
349  }
350  }
351 
352  if (!foundValidB)
353  {
354  // No valid boundary intersection found beyond tracking position
355  break;
356  }
357 
358  // Update starting point for tracking
359  trackFacei = bFacei;
360  trackPt = pushIn(bPoint, trackFacei);
361  trackCelli = getBoundaryCell(trackFacei);
362 
363  segmentI++;
364 
365  startSegmentI = samplingPts.size();
366  }
367 
368  const_cast<polyMesh&>(mesh()).moving(oldMoving);
369 }
370 
371 
372 void Foam::uniformSet::genSamples()
373 {
374  // Storage for sample points
375  DynamicList<point> samplingPts;
376  DynamicList<label> samplingCells;
377  DynamicList<label> samplingFaces;
378  DynamicList<label> samplingSegments;
379  DynamicList<scalar> samplingCurveDist;
380 
381  calcSamples
382  (
383  samplingPts,
384  samplingCells,
385  samplingFaces,
386  samplingSegments,
387  samplingCurveDist
388  );
389 
390  samplingPts.shrink();
391  samplingCells.shrink();
392  samplingFaces.shrink();
393  samplingSegments.shrink();
394  samplingCurveDist.shrink();
395 
396  setSamples
397  (
398  samplingPts,
399  samplingCells,
400  samplingFaces,
401  samplingSegments,
402  samplingCurveDist
403  );
404 }
405 
406 
407 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
408 
410 (
411  const word& name,
412  const polyMesh& mesh,
413  const meshSearch& searchEngine,
414  const word& axis,
415  const point& start,
416  const point& end,
417  const label nPoints
418 )
419 :
420  sampledSet(name, mesh, searchEngine, axis),
421  start_(start),
422  end_(end),
423  nPoints_(nPoints)
424 {
425  genSamples();
426 
427  if (debug)
428  {
429  write(Pout);
430  }
431 }
432 
433 
435 (
436  const word& name,
437  const polyMesh& mesh,
438  const meshSearch& searchEngine,
439  const dictionary& dict
440 )
441 :
442  sampledSet(name, mesh, searchEngine, dict),
443  start_(dict.lookup("start")),
444  end_(dict.lookup("end")),
445  nPoints_(readLabel(dict.lookup("nPoints")))
446 {
447  genSamples();
448 
449  if (debug)
450  {
451  write(Pout);
452  }
453 }
454 
455 
456 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
457 
459 {}
460 
461 
462 // ************************************************************************* //
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
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:458
Macros for easy insertion into run-time selection tables.
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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:410
A class for handling words, derived from string.
Definition: word.H:59
label readLabel(Istream &is)
Definition: label.H:64
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
vector point
Point is a vector.
Definition: point.H:41
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:805
static const scalar tol
Tolerance when comparing points relative to difference between.
Definition: uniformSet.H:119
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576