polyLineSet.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 "polyLineSet.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(polyLineSet, 0);
38  addToRunTimeSelectionTable(sampledSet, polyLineSet, word);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 (
46  passiveParticleCloud& particleCloud,
47  passiveParticle& singleParticle,
48  label& sampleI,
49  DynamicList<point>& samplingPts,
50  DynamicList<label>& samplingCells,
51  DynamicList<label>& samplingFaces,
52  DynamicList<scalar>& samplingCurveDist
53 ) const
54 {
55  particle::TrackingData<passiveParticleCloud> trackData(particleCloud);
56 
57  // Alias
58  const point& trackPt = singleParticle.position();
59 
60  while (true)
61  {
62  // Local geometry info
63  const vector offset = sampleCoords_[sampleI+1] - sampleCoords_[sampleI];
64  const scalar smallDist = mag(tol*offset);
65 
66  point oldPos = trackPt;
67  label facei = -1;
68  do
69  {
70  singleParticle.stepFraction() = 0;
71  singleParticle.track(sampleCoords_[sampleI+1], trackData);
72  }
73  while
74  (
75  !singleParticle.onBoundary()
76  && (mag(trackPt - oldPos) < smallDist)
77  );
78 
79  if (singleParticle.onBoundary())
80  {
81  //Info<< "trackToBoundary : reached boundary"
82  // << " trackPt:" << trackPt << endl;
83  if
84  (
85  mag(trackPt - sampleCoords_[sampleI+1])
86  < smallDist
87  )
88  {
89  // Reached samplePt on boundary
90  //Info<< "trackToBoundary : boundary. also sampling."
91  // << " trackPt:" << trackPt << " sampleI+1:" << sampleI+1
92  // << endl;
93  samplingPts.append(trackPt);
94  samplingCells.append(singleParticle.cell());
95  samplingFaces.append(facei);
96 
97  // trackPt is at sampleI+1
98  samplingCurveDist.append(1.0*(sampleI+1));
99  }
100  return true;
101  }
102 
103  // Reached samplePt in cell
104  samplingPts.append(trackPt);
105  samplingCells.append(singleParticle.cell());
106  samplingFaces.append(-1);
107 
108  // Convert trackPt to fraction inbetween sampleI and sampleI+1
109  scalar dist =
110  mag(trackPt - sampleCoords_[sampleI])
111  / mag(sampleCoords_[sampleI+1] - sampleCoords_[sampleI]);
112  samplingCurveDist.append(sampleI + dist);
113 
114  // go to next samplePt
115  sampleI++;
116 
117  if (sampleI == sampleCoords_.size() - 1)
118  {
119  // no more samples.
120  //Info<< "trackToBoundary : Reached end : sampleI now:" << sampleI
121  // << endl;
122  return false;
123  }
124  }
125 }
126 
127 
129 (
130  DynamicList<point>& samplingPts,
131  DynamicList<label>& samplingCells,
132  DynamicList<label>& samplingFaces,
133  DynamicList<label>& samplingSegments,
134  DynamicList<scalar>& samplingCurveDist
135 ) const
136 {
137  // Check sampling points
138  if (sampleCoords_.size() < 2)
139  {
141  << "Incorrect sample specification. Too few points:"
142  << sampleCoords_ << exit(FatalError);
143  }
144  point oldPoint = sampleCoords_[0];
145  for (label sampleI = 1; sampleI < sampleCoords_.size(); sampleI++)
146  {
147  if (mag(sampleCoords_[sampleI] - oldPoint) < SMALL)
148  {
150  << "Incorrect sample specification."
151  << " Point " << sampleCoords_[sampleI-1]
152  << " at position " << sampleI-1
153  << " and point " << sampleCoords_[sampleI]
154  << " at position " << sampleI
155  << " are too close" << exit(FatalError);
156  }
157  oldPoint = sampleCoords_[sampleI];
158  }
159 
160  // Force calculation of cloud addressing on all processors
161  const bool oldMoving = const_cast<polyMesh&>(mesh()).moving(false);
162  passiveParticleCloud particleCloud(mesh());
163 
164  // current segment number
165  label segmentI = 0;
166 
167  // starting index of current segment in samplePts
168  label startSegmentI = 0;
169 
170  label sampleI = 0;
171 
172  point lastSample(GREAT, GREAT, GREAT);
173  while (true)
174  {
175  // Get boundary intersection
176  point trackPt;
177  label trackCellI = -1;
178  label trackFaceI = -1;
179 
180  do
181  {
182  const vector offset =
183  sampleCoords_[sampleI+1] - sampleCoords_[sampleI];
184  const scalar smallDist = mag(tol*offset);
185 
186 
187  // Get all boundary intersections
188  List<pointIndexHit> bHits = searchEngine().intersections
189  (
190  sampleCoords_[sampleI],
191  sampleCoords_[sampleI+1]
192  );
193 
194  point bPoint(GREAT, GREAT, GREAT);
195  label bFaceI = -1;
196 
197  if (bHits.size())
198  {
199  bPoint = bHits[0].hitPoint();
200  bFaceI = bHits[0].index();
201  }
202 
203  // Get tracking point
204 
205  bool isSample =
206  getTrackingPoint
207  (
208  sampleCoords_[sampleI+1] - sampleCoords_[sampleI],
209  sampleCoords_[sampleI],
210  bPoint,
211  bFaceI,
212 
213  trackPt,
214  trackCellI,
215  trackFaceI
216  );
217 
218  if (isSample && (mag(lastSample - trackPt) > smallDist))
219  {
220  //Info<< "calcSamples : getTrackingPoint returned valid sample "
221  // << " trackPt:" << trackPt
222  // << " trackFaceI:" << trackFaceI
223  // << " trackCellI:" << trackCellI
224  // << " sampleI:" << sampleI
225  // << " dist:" << dist
226  // << endl;
227 
228  samplingPts.append(trackPt);
229  samplingCells.append(trackCellI);
230  samplingFaces.append(trackFaceI);
231 
232  // Convert sampling position to unique curve parameter. Get
233  // fraction of distance between sampleI and sampleI+1.
234  scalar dist =
235  mag(trackPt - sampleCoords_[sampleI])
236  / mag(sampleCoords_[sampleI+1] - sampleCoords_[sampleI]);
237  samplingCurveDist.append(sampleI + dist);
238 
239  lastSample = trackPt;
240  }
241 
242  if (trackCellI == -1)
243  {
244  // No intersection found. Go to next point
245  sampleI++;
246  }
247  } while ((trackCellI == -1) && (sampleI < sampleCoords_.size() - 1));
248 
249  if (sampleI == sampleCoords_.size() - 1)
250  {
251  //Info<< "calcSamples : Reached end of samples: "
252  // << " sampleI now:" << sampleI
253  // << endl;
254  break;
255  }
256 
257  //
258  // Segment sampleI .. sampleI+1 intersected by domain
259  //
260 
261  // Initialize tracking starting from sampleI
262  passiveParticle singleParticle
263  (
264  mesh(),
265  trackPt,
266  trackCellI
267  );
268 
269  bool bReached = trackToBoundary
270  (
271  particleCloud,
272  singleParticle,
273  sampleI,
274  samplingPts,
275  samplingCells,
276  samplingFaces,
277  samplingCurveDist
278  );
279 
280  // fill sampleSegments
281  for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
282  {
283  samplingSegments.append(segmentI);
284  }
285 
286  if (!bReached)
287  {
288  //Info<< "calcSamples : Reached end of samples: "
289  // << " sampleI now:" << sampleI
290  // << endl;
291  break;
292  }
293  lastSample = singleParticle.position();
294 
295 
296  // Find next boundary.
297  sampleI++;
298 
299  if (sampleI == sampleCoords_.size() - 1)
300  {
301  //Info<< "calcSamples : Reached end of samples: "
302  // << " sampleI now:" << sampleI
303  // << endl;
304  break;
305  }
306 
307  segmentI++;
308 
309  startSegmentI = samplingPts.size();
310  }
311 
312  const_cast<polyMesh&>(mesh()).moving(oldMoving);
313 }
314 
315 
317 {
318  // Storage for sample points
319  DynamicList<point> samplingPts;
320  DynamicList<label> samplingCells;
321  DynamicList<label> samplingFaces;
322  DynamicList<label> samplingSegments;
323  DynamicList<scalar> samplingCurveDist;
324 
326  (
327  samplingPts,
328  samplingCells,
329  samplingFaces,
330  samplingSegments,
331  samplingCurveDist
332  );
333 
334  samplingPts.shrink();
335  samplingCells.shrink();
336  samplingFaces.shrink();
337  samplingSegments.shrink();
338  samplingCurveDist.shrink();
339 
340  setSamples
341  (
342  samplingPts,
343  samplingCells,
344  samplingFaces,
345  samplingSegments,
346  samplingCurveDist
347  );
348 }
349 
350 
351 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
352 
354 (
355  const word& name,
356  const polyMesh& mesh,
357  const meshSearch& searchEngine,
358  const word& axis,
359  const List<point>& sampleCoords
360 )
361 :
362  sampledSet(name, mesh, searchEngine, axis),
363  sampleCoords_(sampleCoords)
364 {
365  genSamples();
366 
367  if (debug)
368  {
369  write(Info);
370  }
371 }
372 
373 
375 (
376  const word& name,
377  const polyMesh& mesh,
378  const meshSearch& searchEngine,
379  const dictionary& dict
380 )
381 :
382  sampledSet(name, mesh, searchEngine, dict),
383  sampleCoords_(dict.lookup("points"))
384 {
385  genSamples();
386 
387  if (debug)
388  {
389  write(Info);
390  }
391 }
392 
393 
394 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
395 
397 {}
398 
399 
400 // ************************************************************************* //
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::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::polyLineSet::calcSamples
void calcSamples(DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< label > &samplingSegments, DynamicList< scalar > &samplingCurveDist) const
Samples all point in sampleCoords_.
Definition: polyLineSet.C:129
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::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
Foam::particle::cell
label & cell()
Return current cell particle is in.
Definition: particleI.H:598
Foam::polyLineSet::trackToBoundary
bool trackToBoundary(passiveParticleCloud &particleCloud, passiveParticle &singleParticle, label &sampleI, DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< scalar > &samplingCurveDist) const
Sample till hits boundary. Called with singleParticle at position.
Definition: polyLineSet.C:45
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::Info
messageStream Info
polyLineSet.H
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
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
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
Foam::polyLineSet::polyLineSet
polyLineSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const List< point > &samplePoints)
Construct from components.
Definition: polyLineSet.C:354
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
Foam::polyLineSet::genSamples
void genSamples()
Uses calcSamples to obtain samples. Copies them into *this.
Definition: polyLineSet.C:316
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::polyLineSet::~polyLineSet
virtual ~polyLineSet()
Destructor.
Definition: polyLineSet.C:396
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::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::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