points.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-2024 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 "points.H"
27 #include "meshSearch.H"
28 #include "DynamicList.H"
29 #include "polyMesh.H"
30 #include "sampledSetCloud.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace sampledSets
38 {
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
46 
48 (
49  const polyMesh& mesh,
50  const meshSearch& searchEngine,
51  const pointField& points,
52  DynamicList<point>& samplingPositions,
53  DynamicList<scalar>& samplingDistances,
54  DynamicList<label>& samplingSegments,
55  DynamicList<label>& samplingCells,
56  DynamicList<label>& samplingFaces
57 )
58 {
59  // Create a cloud with which to track segments
60  sampledSetCloud particles
61  (
62  mesh,
63  points::typeName,
65  );
66 
67  // Consider each point
68  label segmenti = 0, samplei = 0, pointi0 = labelMax, pointi = 0;
69  label nLocateBoundaryHits = 0;
70  scalar distance = 0;
71  while (pointi < points.size())
72  {
73  // Sum the distance to the start of the track
74  for (label pointj = pointi0; pointj < pointi; ++ pointj)
75  {
76  distance += mag(points[pointj + 1] - points[pointj]);
77  }
78 
79  // Update the old point index
80  pointi0 = pointi;
81 
82  // Get unique processor and cell that this sample point is in
83  const labelPair procAndCelli = returnReduce
84  (
85  labelPair
86  (
89  ),
90  [](const labelPair& a, const labelPair& b)
91  {
92  return
93  a.second() != -1 && b.second() != -1
94  ? a.first() < b.first() ? a : b
95  : a.second() != -1 ? a : b;
96  }
97  );
98 
99  // Skip this point if it is not in the global mesh
100  if (procAndCelli.second() == -1)
101  {
102  ++ pointi;
103  }
104 
105  // If the point is in the global mesh then track to create a segment
106  else
107  {
109  (
110  particles,
111  points,
112  true,
113  false,
114  false,
115  samplingPositions,
116  samplingDistances,
117  samplingCells,
118  samplingFaces
119  );
120 
121  // Clear the cloud, then, if the point is in this local mesh,
122  // initialise a particle at the point
123  particles.clear();
124  if (procAndCelli.first() == Pstream::myProcNo())
125  {
126  particles.addParticle
127  (
129  (
130  mesh,
131  points[pointi],
132  procAndCelli.second(),
133  nLocateBoundaryHits,
134  pointi,
135  1,
136  distance
137  )
138  );
139 
140  particles.first()->store(particles, td);
141  }
142 
143  // Track to create this segment
144  particles.move(particles, td);
145 
146  // Set the segment indices
147  samplingSegments.append
148  (
149  labelList
150  (
151  samplingPositions.size() - samplingSegments.size(),
152  segmenti
153  )
154  );
155 
156  // Move on to the next segment
157  ++ segmenti;
158 
159  // Determine the global number of samples completed
160  const label samplei0 = samplei;
161  samplei = returnReduce(samplingPositions.size(), sumOp<label>());
162 
163  // Move to the next unsampled point
164  pointi += samplei - samplei0;
165  }
166  }
167 }
168 
169 
170 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
171 
172 void Foam::sampledSets::points::calcSamplesUnordered
173 (
174  DynamicList<point>& samplingPositions,
175  DynamicList<label>& samplingSegments,
176  DynamicList<label>& samplingCells,
177  DynamicList<label>& samplingFaces
178 ) const
179 {
180  forAll(points_, i)
181  {
182  const point& pt = points_[i];
183  const label celli = searchEngine().findCell(pt);
184 
185  if (celli != -1)
186  {
187  samplingPositions.append(pt);
188  samplingSegments.append(i);
189  samplingCells.append(celli);
190  samplingFaces.append(-1);
191  }
192  }
193 }
194 
195 
196 void Foam::sampledSets::points::calcSamplesOrdered
197 (
198  DynamicList<point>& samplingPositions,
199  DynamicList<scalar>& samplingDistances,
200  DynamicList<label>& samplingSegments,
201  DynamicList<label>& samplingCells,
202  DynamicList<label>& samplingFaces
203 ) const
204 {
205  // Calculate the sampling topology
206  calcSamples
207  (
208  mesh(),
209  searchEngine(),
210  pointField(points_),
211  samplingPositions,
212  samplingDistances,
213  samplingSegments,
214  samplingCells,
215  samplingFaces
216  );
217 }
218 
219 
220 void Foam::sampledSets::points::genSamples()
221 {
222  DynamicList<point> samplingPositions;
223  DynamicList<scalar> samplingDistances;
224  DynamicList<label> samplingSegments;
225  DynamicList<label> samplingCells;
226  DynamicList<label> samplingFaces;
227 
228  if (!ordered_)
229  {
230  calcSamplesUnordered
231  (
232  samplingPositions,
233  samplingSegments,
234  samplingCells,
235  samplingFaces
236  );
237  }
238  else
239  {
240  calcSamplesOrdered
241  (
242  samplingPositions,
243  samplingDistances,
244  samplingSegments,
245  samplingCells,
246  samplingFaces
247  );
248  }
249 
250  samplingPositions.shrink();
251  samplingDistances.shrink();
252  samplingSegments.shrink();
253  samplingCells.shrink();
254  samplingFaces.shrink();
255 
256  if (!ordered_)
257  {
258  setSamples
259  (
260  samplingPositions,
261  samplingSegments,
262  samplingCells,
263  samplingFaces
264  );
265  }
266  else
267  {
268  setSamples
269  (
270  samplingPositions,
271  samplingDistances,
272  samplingSegments,
273  samplingCells,
274  samplingFaces
275  );
276  }
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
281 
283 (
284  const word& name,
285  const polyMesh& mesh,
286  const meshSearch& searchEngine,
287  const dictionary& dict
288 )
289 :
290  sampledSet(name, mesh, searchEngine, dict),
291  points_(dict.lookup("points")),
292  ordered_(dict.lookup<bool>("ordered"))
293 {
294  genSamples();
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
299 
301 {}
302 
303 
304 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
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
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Template class for intrusive linked lists.
Definition: ILList.H:67
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const Type & second() const
Return second.
Definition: Pair.H:110
const Type & first() const
Return first.
Definition: Pair.H:98
T * first()
Return the first entry.
Definition: UILList.H:109
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the size.
Definition: coordSet.H:135
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
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:750
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
Specified point samples. Optionally ordered into a continuous path. Ordering is an optimisation; it e...
Definition: points.H:106
virtual ~points()
Destructor.
Definition: points.C:300
static void calcSamples(const polyMesh &mesh, const meshSearch &searchEngine, const pointField &points, DynamicList< point > &samplingPositons, DynamicList< scalar > &samplingDistances, DynamicList< label > &samplingSegments, DynamicList< label > &samplingCells, DynamicList< label > &samplingFaces)
Calculate all the sampling points.
Definition: points.C:48
points(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: points.C:283
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.
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
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
static const label labelMax
Definition: label.H:62
dictionary dict