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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faceOnlySet.H"
30 #include "meshSearch.H"
31 #include "DynamicList.H"
32 #include "polyMesh.H"
33 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(faceOnlySet, 0);
41  addToRunTimeSelectionTable(sampledSet, faceOnlySet, word);
42 }
43 
44 const Foam::scalar Foam::faceOnlySet::tol = ROOTSMALL;
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 bool Foam::faceOnlySet::trackToBoundary
50 (
51  passiveParticleCloud& particleCloud,
52  passiveParticle& singleParticle,
53  const scalar smallDist,
54  DynamicList<point>& samplingPts,
55  DynamicList<label>& samplingCells,
56  DynamicList<label>& samplingFaces,
57  DynamicList<scalar>& samplingCurveDist
58 ) const
59 {
60  const vector offset = (end_ - start_);
61 
62  particle::trackingData td(particleCloud);
63 
64  point trackPt = singleParticle.position();
65 
66  while (true)
67  {
68  point oldPoint = trackPt;
69 
70  singleParticle.trackToAndHitFace(end_ - start_, 0, particleCloud, td);
71 
72  trackPt = singleParticle.position();
73 
74  if (singleParticle.face() != -1 && mag(oldPoint - trackPt) > smallDist)
75  {
76  // Reached face. Sample.
77  samplingPts.append(trackPt);
78  samplingCells.append(singleParticle.cell());
79  samplingFaces.append(singleParticle.face());
80  samplingCurveDist.append(mag(trackPt - start_));
81  }
82 
83  if (-smallDist < ((trackPt - end_) & offset))
84  {
85  // Projected onto sampling vector
86  // - done when we are near or past the end of the sampling vector
87  return false;
88  }
89  else if (singleParticle.onBoundaryFace())
90  {
91  // Boundary reached
92  return true;
93  }
94  }
95 }
96 
97 
98 void Foam::faceOnlySet::calcSamples
99 (
100  DynamicList<point>& samplingPts,
101  DynamicList<label>& samplingCells,
102  DynamicList<label>& samplingFaces,
103  DynamicList<label>& samplingSegments,
104  DynamicList<scalar>& samplingCurveDist
105 ) const
106 {
107  // Distance vector between sampling points
108  if (mag(end_ - start_) < SMALL)
109  {
111  << "Incorrect sample specification :"
112  << " start equals end point." << endl
113  << " start:" << start_
114  << " end:" << end_
115  << exit(FatalError);
116  }
117 
118  const vector offset = (end_ - start_);
119  const vector normOffset = offset/mag(offset);
120  const vector smallVec = tol*offset;
121  const scalar smallDist = mag(smallVec);
122 
123  // Force calculation of cloud addressing on all processors
124  const bool oldMoving = const_cast<polyMesh&>(mesh()).moving(false);
125  passiveParticleCloud particleCloud(mesh());
126 
127  // Get all boundary intersections
128  List<pointIndexHit> bHits = searchEngine().intersections
129  (
130  start_ - smallVec,
131  end_ + smallVec
132  );
133 
134  point bPoint(GREAT, GREAT, GREAT);
135  label bFacei = -1;
136 
137  if (bHits.size())
138  {
139  bPoint = bHits[0].hitPoint();
140  bFacei = bHits[0].index();
141  }
142 
143  // Get first tracking point. Use bPoint, bFacei if provided.
144  point trackPt;
145  label trackCelli = -1;
146  label trackFacei = -1;
147 
148  // Pout<< "before getTrackingPoint : bPoint:" << bPoint
149  // << " bFacei:" << bFacei << endl;
150 
151  getTrackingPoint
152  (
153  start_,
154  bPoint,
155  bFacei,
156  smallDist,
157  trackPt,
158  trackCelli,
159  trackFacei
160  );
161 
162  // Pout<< "after getTrackingPoint : "
163  // << " trackPt:" << trackPt
164  // << " trackCelli:" << trackCelli
165  // << " trackFacei:" << trackFacei
166  // << endl;
167 
168  if (trackCelli == -1)
169  {
170  // Line start_ - end_ does not intersect domain at all.
171  // (or is along edge)
172  // Set points and cell/face labels to empty lists
173  //Info<< "calcSamples : Both start_ and end_ outside domain"
174  // << endl;
175 
176  const_cast<polyMesh&>(mesh()).moving(oldMoving);
177  return;
178  }
179 
180  if (trackFacei == -1)
181  {
182  // No boundary face. Check for nearish internal face
183  trackFacei = findNearFace(trackCelli, trackPt, smallDist);
184  }
185 
186  // Pout<< "calcSamples : got first point to track from :"
187  // << " trackPt:" << trackPt
188  // << " trackCell:" << trackCelli
189  // << " trackFace:" << trackFacei
190  // << endl;
191 
192  //
193  // Track until hit end of all boundary intersections
194  //
195 
196  // current segment number
197  label segmentI = 0;
198 
199  // starting index of current segment in samplePts
200  label startSegmentI = 0;
201 
202  // index in bHits; current boundary intersection
203  label bHitI = 1;
204 
205  while (true)
206  {
207  if (trackFacei != -1)
208  {
209  // Pout<< "trackPt:" << trackPt << " on face so use." << endl;
210  samplingPts.append(trackPt);
211  samplingCells.append(trackCelli);
212  samplingFaces.append(trackFacei);
213  samplingCurveDist.append(mag(trackPt - start_));
214  }
215 
216  // Initialize tracking starting from trackPt
217  passiveParticle singleParticle
218  (
219  mesh(),
220  trackPt,
221  trackCelli
222  );
223 
224  bool reachedBoundary = trackToBoundary
225  (
226  particleCloud,
227  singleParticle,
228  smallDist,
229  samplingPts,
230  samplingCells,
231  samplingFaces,
232  samplingCurveDist
233  );
234 
235  // Fill sampleSegments
236  for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
237  {
238  samplingSegments.append(segmentI);
239  }
240 
241  if (!reachedBoundary)
242  {
243  // Pout<< "calcSamples : Reached end of samples: "
244  // << " samplePt now:" << singleParticle.position()
245  // << endl;
246  break;
247  }
248 
249  bool foundValidB = false;
250 
251  while (bHitI < bHits.size())
252  {
253  scalar dist =
254  (bHits[bHitI].hitPoint() - singleParticle.position())
255  & normOffset;
256 
257  // Pout<< "Finding next boundary : "
258  // << "bPoint:" << bHits[bHitI].hitPoint()
259  // << " tracking:" << singleParticle.position()
260  // << " dist:" << dist
261  // << endl;
262 
263  if (dist > smallDist)
264  {
265  // Hit-point is past tracking position
266  foundValidB = true;
267  break;
268  }
269  else
270  {
271  ++bHitI;
272  }
273  }
274 
275  if (!foundValidB || bHitI == bHits.size() - 1)
276  {
277  // No valid boundary intersection found beyond tracking position
278  break;
279  }
280 
281  // Update starting point for tracking
282  trackFacei = bHits[bHitI].index();
283  trackPt = pushIn(bHits[bHitI].hitPoint(), trackFacei);
284  trackCelli = getBoundaryCell(trackFacei);
285 
286  ++segmentI;
287 
288  startSegmentI = samplingPts.size();
289  }
290 
291  const_cast<polyMesh&>(mesh()).moving(oldMoving);
292 }
293 
294 
295 void Foam::faceOnlySet::genSamples()
296 {
297  // Storage for sample points
298  DynamicList<point> samplingPts;
299  DynamicList<label> samplingCells;
300  DynamicList<label> samplingFaces;
301  DynamicList<label> samplingSegments;
302  DynamicList<scalar> samplingCurveDist;
303 
304  calcSamples
305  (
306  samplingPts,
307  samplingCells,
308  samplingFaces,
309  samplingSegments,
310  samplingCurveDist
311  );
312 
313  samplingPts.shrink();
314  samplingCells.shrink();
315  samplingFaces.shrink();
316  samplingSegments.shrink();
317  samplingCurveDist.shrink();
318 
319  // Move into *this
320  setSamples
321  (
322  std::move(samplingPts),
323  std::move(samplingCells),
324  std::move(samplingFaces),
325  std::move(samplingSegments),
326  std::move(samplingCurveDist)
327  );
328 
329  if (debug)
330  {
331  write(Info);
332  }
333 }
334 
335 
336 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
337 
339 (
340  const word& name,
341  const polyMesh& mesh,
342  const meshSearch& searchEngine,
343  const word& axis,
344  const point& start,
345  const point& end
346 )
347 :
348  sampledSet(name, mesh, searchEngine, axis),
349  start_(start),
350  end_(end)
351 {
352  genSamples();
353 }
354 
355 
357 (
358  const word& name,
359  const polyMesh& mesh,
360  const meshSearch& searchEngine,
361  const dictionary& dict
362 )
363 :
364  sampledSet(name, mesh, searchEngine, dict),
365  start_(dict.get<point>("start")),
366  end_(dict.get<point>("end"))
367 {
368  genSamples();
369 }
370 
371 
372 // ************************************************************************* //
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
faceOnlySet.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:339
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
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
Foam::faceOnlySet::tol
static const scalar tol
Tolerance when comparing points relative to difference between.
Definition: faceOnlySet.H:141
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
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
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::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)