lineFace.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 "lineFace.H"
27 #include "meshSearch.H"
28 #include "DynamicList.H"
29 #include "polyMesh.H"
30 #include "treeDataCell.H"
31 #include "sampledSetCloud.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace sampledSets
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
47 
48 void Foam::sampledSets::lineFace::calcSamples
49 (
50  const polyMesh& mesh,
51  const meshSearch& searchEngine,
52  const vector& start,
53  const vector& end,
54  const label storeFaces,
55  const bool storeCells,
56  DynamicList<point>& samplingPositions,
57  DynamicList<scalar>& samplingDistances,
58  DynamicList<label>& samplingSegments,
59  DynamicList<label>& samplingCells,
60  DynamicList<label>& samplingFaces
61 )
62 {
63  // Get all candidates for starting the tracks
64  List<DynamicList<label>> procCandidateCells(Pstream::nProcs());
65  List<DynamicList<scalar>> procCandidateTs(Pstream::nProcs());
66  {
67  const label startCelli = searchEngine.findCell(start);
68  if (startCelli != -1)
69  {
70  procCandidateCells[Pstream::myProcNo()].append(startCelli);
71  procCandidateTs[Pstream::myProcNo()].append(0);
72  }
73 
74  const label endCelli = searchEngine.findCell(end);
75  if (endCelli != -1)
76  {
77  procCandidateCells[Pstream::myProcNo()].append(endCelli);
78  procCandidateTs[Pstream::myProcNo()].append(1);
79  }
80 
81  const List<pointIndexHit> bHits =
82  searchEngine.intersections(start, end);
83  forAll(bHits, bHiti)
84  {
85  for (label bHitj = bHiti + 1; bHitj < bHits.size(); ++ bHitj)
86  {
87  const point midP =
88  (bHits[bHiti].hitPoint() + bHits[bHitj].hitPoint())/2;
89 
90  const label midCelli = searchEngine.findCell(midP);
91  const scalar midT = mag(midP - start)/mag(end - start);
92 
93  if (midCelli != -1)
94  {
95  procCandidateCells[Pstream::myProcNo()].append(midCelli);
96  procCandidateTs[Pstream::myProcNo()].append(midT);
97  }
98  }
99  }
100  }
101  Pstream::gatherList(procCandidateCells);
102  Pstream::scatterList(procCandidateCells);
103  Pstream::gatherList(procCandidateTs);
104  Pstream::scatterList(procCandidateTs);
105 
106  // Tracking data
107  const List<point> startEnd({start, end});
108  const List<point> endStart({end, start});
109  DynamicList<point> halfSegmentPositions;
110  DynamicList<scalar> halfSegmentDistances;
111  DynamicList<label> halfSegmentCells;
112  DynamicList<label> halfSegmentFaces;
113 
114  // Create a cloud with which to track segments
115  sampledSetCloud particles
116  (
117  mesh,
118  lineFace::typeName,
120  );
121 
122  // Create each segment in turn
123  label segmenti = 0;
124  label nLocateBoundaryHits = 0;
125  forAll(procCandidateCells, proci)
126  {
127  forAll(procCandidateCells[proci], candidatei)
128  {
129  const label celli = procCandidateCells[proci][candidatei];
130 
131  if (celli == -1) continue;
132 
133  const scalar t = procCandidateTs[proci][candidatei];
134  const point p = (1 - t)*start + t*end;
135 
136  // Track in both directions to form parts of the segment either
137  // side of the candidate
138  Pair<scalar> segmentT;
139  forAll(segmentT, i)
140  {
142  (
143  particles,
144  i == 0 ? endStart : startEnd,
145  false,
146  storeFaces,
147  storeCells,
148  halfSegmentPositions,
149  halfSegmentDistances,
150  halfSegmentCells,
151  halfSegmentFaces
152  );
153 
154  particles.clear();
155  if (proci == Pstream::myProcNo())
156  {
157  particles.addParticle
158  (
160  (
161  mesh,
162  p,
163  celli,
164  nLocateBoundaryHits,
165  0,
166  i == 0 ? t : 1 - t,
167  0
168  )
169  );
170  }
171 
172  particles.move(particles, tdBwd);
173 
174  segmentT[i] =
176  (
177  particles.size()
178  ? mag(particles.first()->position(mesh) - start)
179  /mag(end - start)
180  : i == 0 ? vGreat : -vGreat,
181  [i](const scalar a, const scalar b)
182  {
183  return (i == 0) == (a < b) ? a : b;
184  }
185  );
186 
187  if (i == 0)
188  {
189  reverse(halfSegmentPositions);
190  reverse(halfSegmentDistances);
191  reverse(halfSegmentCells);
192  reverse(halfSegmentFaces);
193  }
194 
195  const label n = halfSegmentPositions.size();
196 
197  samplingPositions.append(halfSegmentPositions);
198  samplingSegments.append(labelList(n, segmenti));
199  samplingCells.append(halfSegmentCells);
200  samplingFaces.append(halfSegmentFaces);
201 
202  halfSegmentPositions.clear();
203  halfSegmentDistances.clear();
204  halfSegmentCells.clear();
205  halfSegmentFaces.clear();
206 
207  // If storing cells we need to store the starting cells between
208  // the tracks
209  if (proci == Pstream::myProcNo() && i == 0 && storeCells)
210  {
211  particle trackBwd(mesh, p, celli, nLocateBoundaryHits);
212  particle trackFwd(trackBwd);
213  trackBwd.trackToFace(mesh, start - p, 0);
214  trackFwd.trackToFace(mesh, end - p, 0);
215  if (trackBwd.onFace() && trackFwd.onFace())
216  {
217  samplingPositions.append
218  (
219  (
220  trackBwd.position(mesh)
221  + trackFwd.position(mesh)
222  )/2
223  );
224  samplingSegments.append(segmenti);
225  samplingCells.append(celli);
226  samplingFaces.append(-1);
227  }
228  }
229  }
230 
231  // Disable all candidates that fall within the bounds of the
232  // computed segment
233  forAll(procCandidateCells, procj)
234  {
235  forAll(procCandidateCells[procj], candidatej)
236  {
237  const label cellj = procCandidateCells[procj][candidatej];
238 
239  if (cellj == -1) continue;
240 
241  const scalar t = procCandidateTs[procj][candidatej];
242 
243  if
244  (
245  t > segmentT.first() - rootSmall
246  && t < segmentT.second() + rootSmall
247  )
248  {
249  procCandidateCells[procj][candidatej] = -1;
250  procCandidateTs[procj][candidatej] = NaN;
251  }
252  }
253  }
254 
255  // Move onto the next segment
256  segmenti ++;
257  }
258  }
259 
260  // Set the distances
261  samplingDistances.resize(samplingPositions.size());
262  forAll(samplingPositions, i)
263  {
264  samplingDistances[i] = mag(samplingPositions[i] - start);
265  }
266 }
267 
268 
269 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
270 
271 void Foam::sampledSets::lineFace::calcSamples
272 (
273  DynamicList<point>& samplingPositions,
274  DynamicList<scalar>& samplingDistances,
275  DynamicList<label>& samplingSegments,
276  DynamicList<label>& samplingCells,
277  DynamicList<label>& samplingFaces
278 ) const
279 {
280  calcSamples
281  (
282  mesh(),
283  searchEngine(),
284  start_,
285  end_,
286  1,
287  false,
288  samplingPositions,
289  samplingDistances,
290  samplingSegments,
291  samplingCells,
292  samplingFaces
293  );
294 }
295 
296 
297 void Foam::sampledSets::lineFace::genSamples()
298 {
299  DynamicList<point> samplingPositions;
300  DynamicList<scalar> samplingDistances;
301  DynamicList<label> samplingSegments;
302  DynamicList<label> samplingCells;
303  DynamicList<label> samplingFaces;
304 
305  calcSamples
306  (
307  samplingPositions,
308  samplingDistances,
309  samplingSegments,
310  samplingCells,
311  samplingFaces
312  );
313 
314  samplingPositions.shrink();
315  samplingDistances.shrink();
316  samplingSegments.shrink();
317  samplingCells.shrink();
318  samplingFaces.shrink();
319 
320  setSamples
321  (
322  samplingPositions,
323  samplingDistances,
324  samplingSegments,
325  samplingCells,
326  samplingFaces
327  );
328 }
329 
330 
331 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
332 
334 (
335  const word& name,
336  const polyMesh& mesh,
337  const meshSearch& searchEngine,
338  const dictionary& dict
339 )
340 :
341  sampledSet(name, mesh, searchEngine, dict),
342  start_(dict.lookup("start")),
343  end_(dict.lookup("end"))
344 {
345  genSamples();
346 }
347 
348 
349 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
350 
352 {}
353 
354 
355 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:191
void clear()
Definition: Cloud.H:239
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:219
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:281
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
void resize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:216
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
Template class for intrusive linked lists.
Definition: ILList.H:67
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
const Type & second() const
Return second.
Definition: Pair.H:110
const Type & first() const
Return first.
Definition: Pair.H:98
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
T * first()
Return the first entry.
Definition: UILList.H:109
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:58
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:861
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:750
Base particle class.
Definition: particle.H:83
scalar trackToFace(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:565
bool onFace() const
Is the particle on a face?
Definition: particleI.H:231
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:255
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A Cloud of sampledSet particles.
Particle for generating line-type sampled sets.
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:64
const meshSearch & searchEngine() const
Access the search engine.
Definition: sampledSet.H:216
const polyMesh & mesh() const
Access the mesh.
Definition: sampledSet.H:210
Face-intersections along a line.
Definition: lineFace.H:92
virtual ~lineFace()
Destructor.
Definition: lineFace.C:351
lineFace(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: lineFace.C:334
Set of sets to sample. Call sampledSets.write() to sample&write files.
A class for handling words, derived from string.
Definition: word.H:62
volScalarField & b
Definition: createFields.H:25
defineTypeNameAndDebug(arcUniform, 0)
addToRunTimeSelectionTable(sampledSet, arcUniform, word)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
dictionary dict
volScalarField & p