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