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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "uniformSet.H"
29 #include "meshSearch.H"
30 #include "DynamicList.H"
31 #include "polyMesh.H"
32 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(uniformSet, 0);
40  addToRunTimeSelectionTable(sampledSet, uniformSet, word);
41 }
42 
43 const Foam::scalar Foam::uniformSet::tol = 1e-3;
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 bool Foam::uniformSet::nextSample
49 (
50  const point& currentPt,
51  const vector& offset,
52  const scalar smallDist,
53  point& samplePt,
54  label& sampleI
55 ) const
56 {
57  bool pointFound = false;
58 
59  const vector normOffset = offset/mag(offset);
60 
61  samplePt += offset;
62  ++sampleI;
63 
64  for (; sampleI < nPoints_; ++sampleI)
65  {
66  scalar s = (samplePt - currentPt) & normOffset;
67 
68  if (s > -smallDist)
69  {
70  // samplePt is close to or beyond currentPt -> use it
71  pointFound = true;
72 
73  break;
74  }
75  samplePt += offset;
76  }
77 
78  return pointFound;
79 }
80 
81 
82 bool Foam::uniformSet::trackToBoundary
83 (
84  passiveParticleCloud& particleCloud,
85  passiveParticle& singleParticle,
86  point& samplePt,
87  label& sampleI,
88  DynamicList<point>& samplingPts,
89  DynamicList<label>& samplingCells,
90  DynamicList<label>& samplingFaces,
91  DynamicList<scalar>& samplingCurveDist
92 ) const
93 {
94  // distance vector between sampling points
95  const vector offset = (end_ - start_)/(nPoints_ - 1);
96  const vector smallVec = tol*offset;
97  const scalar smallDist = mag(smallVec);
98 
99  point trackPt = singleParticle.position();
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 
260  return;
261  }
262 
263  if (isSample)
264  {
265  samplingPts.append(start_);
266  samplingCells.append(trackCelli);
267  samplingFaces.append(trackFacei);
268  samplingCurveDist.append(0.0);
269  }
270 
271  //
272  // Track until hit end of all boundary intersections
273  //
274 
275  // current segment number
276  label segmentI = 0;
277 
278  // starting index of current segment in samplePts
279  label startSegmentI = 0;
280 
281  label sampleI = 0;
282  point samplePt = start_;
283 
284  // index in bHits; current boundary intersection
285  label bHitI = 1;
286 
287  while (true)
288  {
289  // Initialize tracking starting from trackPt
290  passiveParticle singleParticle(mesh(), trackPt, trackCelli);
291 
292  bool reachedBoundary = trackToBoundary
293  (
294  particleCloud,
295  singleParticle,
296  samplePt,
297  sampleI,
298  samplingPts,
299  samplingCells,
300  samplingFaces,
301  samplingCurveDist
302  );
303 
304  // fill sampleSegments
305  for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
306  {
307  samplingSegments.append(segmentI);
308  }
309 
310 
311  if (!reachedBoundary)
312  {
313  if (debug)
314  {
315  Pout<< "calcSamples : Reached end of samples: "
316  << " samplePt now:" << samplePt
317  << " sampleI now:" << sampleI
318  << endl;
319  }
320  break;
321  }
322 
323 
324  bool foundValidB = false;
325 
326  while (bHitI < bHits.size())
327  {
328  scalar dist =
329  (bHits[bHitI].hitPoint() - singleParticle.position())
330  & normOffset;
331 
332  if (debug)
333  {
334  Pout<< "Finding next boundary : "
335  << "bPoint:" << bHits[bHitI].hitPoint()
336  << " tracking:" << singleParticle.position()
337  << " dist:" << dist
338  << endl;
339  }
340 
341  if (dist > smallDist)
342  {
343  // hitpoint is past tracking position
344  foundValidB = true;
345  break;
346  }
347  else
348  {
349  ++bHitI;
350  }
351  }
352 
353  if (!foundValidB)
354  {
355  // No valid boundary intersection found beyond tracking position
356  break;
357  }
358 
359  // Update starting point for tracking
360  trackFacei = bFacei;
361  trackPt = pushIn(bPoint, trackFacei);
362  trackCelli = getBoundaryCell(trackFacei);
363 
364  ++segmentI;
365 
366  startSegmentI = samplingPts.size();
367  }
368 
369  const_cast<polyMesh&>(mesh()).moving(oldMoving);
370 }
371 
372 
373 void Foam::uniformSet::genSamples()
374 {
375  // Storage for sample points
376  DynamicList<point> samplingPts;
377  DynamicList<label> samplingCells;
378  DynamicList<label> samplingFaces;
379  DynamicList<label> samplingSegments;
380  DynamicList<scalar> samplingCurveDist;
381 
382  calcSamples
383  (
384  samplingPts,
385  samplingCells,
386  samplingFaces,
387  samplingSegments,
388  samplingCurveDist
389  );
390 
391  samplingPts.shrink();
392  samplingCells.shrink();
393  samplingFaces.shrink();
394  samplingSegments.shrink();
395  samplingCurveDist.shrink();
396 
397  // Move into *this
398  setSamples
399  (
400  std::move(samplingPts),
401  std::move(samplingCells),
402  std::move(samplingFaces),
403  std::move(samplingSegments),
404  std::move(samplingCurveDist)
405  );
406 
407  if (debug)
408  {
409  write(Pout);
410  }
411 }
412 
413 
414 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
415 
417 (
418  const word& name,
419  const polyMesh& mesh,
420  const meshSearch& searchEngine,
421  const word& axis,
422  const point& start,
423  const point& end,
424  const label nPoints
425 )
426 :
427  sampledSet(name, mesh, searchEngine, axis),
428  start_(start),
429  end_(end),
430  nPoints_(nPoints)
431 {
432  genSamples();
433 }
434 
435 
437 (
438  const word& name,
439  const polyMesh& mesh,
440  const meshSearch& searchEngine,
441  const dictionary& dict
442 )
443 :
444  sampledSet(name, mesh, searchEngine, dict),
445  start_(dict.get<point>("start")),
446  end_(dict.get<point>("end")),
447  nPoints_(dict.get<label>("nPoints"))
448 {
449  genSamples();
450 }
451 
452 
453 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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:888
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:121
Foam::uniformSet::tol
static const scalar tol
Tolerance when comparing points relative to difference between.
Definition: uniformSet.H:160
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
meshSearch.H
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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:417
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
DynamicList.H
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
uniformSet.H