faceOnlySet.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 "faceOnlySet.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(faceOnlySet, 0);
38  addToRunTimeSelectionTable(sampledSet, faceOnlySet, word);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 (
46  passiveParticleCloud& particleCloud,
47  passiveParticle& singleParticle,
48  DynamicList<point>& samplingPts,
49  DynamicList<label>& samplingCells,
50  DynamicList<label>& samplingFaces,
51  DynamicList<scalar>& samplingCurveDist
52 ) const
53 {
54  // distance vector between sampling points
55  const vector offset = end_ - start_;
56  const vector smallVec = tol*offset;
57  const scalar smallDist = mag(smallVec);
58 
59  particle::TrackingData<passiveParticleCloud> trackData(particleCloud);
60 
61  // Alias
62  const point& trackPt = singleParticle.position();
63 
64  while(true)
65  {
66  point oldPoint = trackPt;
67 
68  singleParticle.trackToFace(end_, trackData);
69 
70  if (singleParticle.face() != -1 && mag(oldPoint - trackPt) > smallDist)
71  {
72  // Reached face. Sample.
73  samplingPts.append(trackPt);
74  samplingCells.append(singleParticle.cell());
75  samplingFaces.append(singleParticle.face());
76  samplingCurveDist.append(mag(trackPt - start_));
77  }
78 
79  if (mag(trackPt - end_) < smallDist)
80  {
81  // end reached
82  return false;
83  }
84  else if (singleParticle.onBoundary())
85  {
86  // Boundary reached.
87  return true;
88  }
89  }
90 }
91 
92 
94 (
95  DynamicList<point>& samplingPts,
96  DynamicList<label>& samplingCells,
97  DynamicList<label>& samplingFaces,
98  DynamicList<label>& samplingSegments,
99  DynamicList<scalar>& samplingCurveDist
100 ) const
101 {
102  // distance vector between sampling points
103  if (mag(end_ - start_) < SMALL)
104  {
106  << "Incorrect sample specification :"
107  << " start equals end point." << endl
108  << " start:" << start_
109  << " end:" << end_
110  << exit(FatalError);
111  }
112 
113  const vector offset = (end_ - start_);
114  const vector normOffset = offset/mag(offset);
115  const vector smallVec = tol*offset;
116  const scalar smallDist = mag(smallVec);
117 
118  // Force calculation of cloud addressing on all processors
119  const bool oldMoving = const_cast<polyMesh&>(mesh()).moving(false);
120  passiveParticleCloud particleCloud(mesh());
121 
122  // Get all boundary intersections
123  List<pointIndexHit> bHits = searchEngine().intersections
124  (
125  start_ - smallVec,
126  end_ + smallVec
127  );
128 
129  point bPoint(GREAT, GREAT, GREAT);
130  label bFaceI = -1;
131 
132  if (bHits.size())
133  {
134  bPoint = bHits[0].hitPoint();
135  bFaceI = bHits[0].index();
136  }
137 
138  // Get first tracking point. Use bPoint, bFaceI if provided.
139 
140  point trackPt;
141  label trackCellI = -1;
142  label trackFaceI = -1;
143 
144  //Info<< "before getTrackingPoint : bPoint:" << bPoint
145  // << " bFaceI:" << bFaceI << endl;
146 
147  getTrackingPoint
148  (
149  offset,
150  start_,
151  bPoint,
152  bFaceI,
153 
154  trackPt,
155  trackCellI,
156  trackFaceI
157  );
158 
159  //Info<< "after getTrackingPoint : "
160  // << " trackPt:" << trackPt
161  // << " trackCellI:" << trackCellI
162  // << " trackFaceI:" << trackFaceI
163  // << endl;
164 
165  if (trackCellI == -1)
166  {
167  // Line start_ - end_ does not intersect domain at all.
168  // (or is along edge)
169  // Set points and cell/face labels to empty lists
170  //Info<< "calcSamples : Both start_ and end_ outside domain"
171  // << endl;
172  const_cast<polyMesh&>(mesh()).moving(oldMoving);
173 
174  return;
175  }
176 
177  if (trackFaceI == -1)
178  {
179  // No boundary face. Check for nearish internal face
180  trackFaceI = findNearFace(trackCellI, trackPt, smallDist);
181  }
182 
183  //Info<< "calcSamples : got first point to track from :"
184  // << " trackPt:" << trackPt
185  // << " trackCell:" << trackCellI
186  // << " trackFace:" << trackFaceI
187  // << endl;
188 
189  //
190  // Track until hit end of all boundary intersections
191  //
192 
193  // current segment number
194  label segmentI = 0;
195 
196  // starting index of current segment in samplePts
197  label startSegmentI = 0;
198 
199  // index in bHits; current boundary intersection
200  label bHitI = 1;
201 
202  while(true)
203  {
204  if (trackFaceI != -1)
205  {
206  //Info<< "trackPt:" << trackPt << " on face so use." << endl;
207  samplingPts.append(trackPt);
208  samplingCells.append(trackCellI);
209  samplingFaces.append(trackFaceI);
210  samplingCurveDist.append(mag(trackPt - start_));
211  }
212 
213  // Initialize tracking starting from trackPt
214  passiveParticle singleParticle
215  (
216  mesh(),
217  trackPt,
218  trackCellI
219  );
220 
221  bool reachedBoundary = trackToBoundary
222  (
223  particleCloud,
224  singleParticle,
225  samplingPts,
226  samplingCells,
227  samplingFaces,
228  samplingCurveDist
229  );
230 
231  // fill sampleSegments
232  for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
233  {
234  samplingSegments.append(segmentI);
235  }
236 
237 
238  if (!reachedBoundary)
239  {
240  //Info<< "calcSamples : Reached end of samples: "
241  // << " samplePt now:" << singleParticle.position()
242  // << endl;
243  break;
244  }
245 
246 
247  // Go past boundary intersection where tracking stopped
248  // Use coordinate comparison instead of face comparison for
249  // accuracy reasons
250 
251  bool foundValidB = false;
252 
253  while (bHitI < bHits.size())
254  {
255  scalar dist =
256  (bHits[bHitI].hitPoint() - singleParticle.position())
257  & normOffset;
258 
259  //Info<< "Finding next boundary : "
260  // << "bPoint:" << bHits[bHitI].hitPoint()
261  // << " tracking:" << singleParticle.position()
262  // << " dist:" << dist
263  // << endl;
264 
265  if (dist > smallDist)
266  {
267  // hitpoint is past tracking position
268  foundValidB = true;
269  break;
270  }
271  else
272  {
273  bHitI++;
274  }
275  }
276 
277  if (!foundValidB)
278  {
279  // No valid boundary intersection found beyond tracking position
280  break;
281  }
282 
283  // Update starting point for tracking
284  trackFaceI = bHits[bHitI].index();
285  trackPt = pushIn(bHits[bHitI].hitPoint(), trackFaceI);
286  trackCellI = getBoundaryCell(trackFaceI);
287 
288  segmentI++;
289 
290  startSegmentI = samplingPts.size();
291  }
292 
293  const_cast<polyMesh&>(mesh()).moving(oldMoving);
294 }
295 
296 
298 {
299  // Storage for sample points
300  DynamicList<point> samplingPts;
301  DynamicList<label> samplingCells;
302  DynamicList<label> samplingFaces;
303  DynamicList<label> samplingSegments;
304  DynamicList<scalar> samplingCurveDist;
305 
307  (
308  samplingPts,
309  samplingCells,
310  samplingFaces,
311  samplingSegments,
312  samplingCurveDist
313  );
314 
315  samplingPts.shrink();
316  samplingCells.shrink();
317  samplingFaces.shrink();
318  samplingSegments.shrink();
319  samplingCurveDist.shrink();
320 
321  // Copy into *this
322  setSamples
323  (
324  samplingPts,
325  samplingCells,
326  samplingFaces,
327  samplingSegments,
328  samplingCurveDist
329  );
330 }
331 
332 
333 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
334 
336 (
337  const word& name,
338  const polyMesh& mesh,
339  const meshSearch& searchEngine,
340  const word& axis,
341  const point& start,
342  const point& end
343 )
344 :
345  sampledSet(name, mesh, searchEngine, axis),
346  start_(start),
347  end_(end)
348 {
349  genSamples();
350 
351  if (debug)
352  {
353  write(Info);
354  }
355 }
356 
357 
359 (
360  const word& name,
361  const polyMesh& mesh,
362  const meshSearch& searchEngine,
363  const dictionary& dict
364 )
365 :
366  sampledSet(name, mesh, searchEngine, dict),
367  start_(dict.lookup("start")),
368  end_(dict.lookup("end"))
369 {
370  genSamples();
371 
372  if (debug)
373  {
374  write(Info);
375  }
376 }
377 
378 
379 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
380 
382 {}
383 
384 
385 // ************************************************************************* //
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::particle::face
label & face()
Return current face particle is on otherwise -1.
Definition: particleI.H:658
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
faceOnlySet.H
Foam::particle::trackToFace
scalar trackToFace(const vector &endPosition, TrackData &td)
Track particle to a given position and returns 1.0 if the.
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::faceOnlySet::trackToBoundary
bool trackToBoundary(passiveParticleCloud &particleCloud, passiveParticle &singleParticle, DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< scalar > &samplingCurve) const
Samples from startTrackPt/CellI. Updates particle/samplePt/sampleI.
Definition: faceOnlySet.C:45
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::faceOnlySet::faceOnlySet
faceOnlySet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const point &start, const point &end)
Construct from components.
Definition: faceOnlySet.C:336
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::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
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::faceOnlySet::calcSamples
void calcSamples(DynamicList< point > &samplingPts, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces, DynamicList< label > &samplingSegments, DynamicList< scalar > &samplingCurveDist) const
Samples from start_ to end_. samplingSegments contains segmentNo.
Definition: faceOnlySet.C:94
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
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::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::position
const vector & position() const
Return current particle position.
Definition: particleI.H:586
Foam::faceOnlySet::genSamples
void genSamples()
Uses calcSamples to obtain samples. Copies them into *this.
Definition: faceOnlySet.C:297
Foam::faceOnlySet::~faceOnlySet
virtual ~faceOnlySet()
Destructor.
Definition: faceOnlySet.C:381
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