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