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-2022 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  scalar distance = 0;
70  while (pointi < points.size())
71  {
72  // Sum the distance to the start of the track
73  for (label pointj = pointi0; pointj < pointi; ++ pointj)
74  {
75  distance += mag(points[pointj + 1] - points[pointj]);
76  }
77 
78  // Update the old point index
79  pointi0 = pointi;
80 
81  // Get unique processor and cell that this sample point is in
82  const labelPair procAndCelli = returnReduce
83  (
84  labelPair
85  (
88  ),
89  [](const labelPair& a, const labelPair& b)
90  {
91  return
92  a.second() != -1 && b.second() != -1
93  ? a.first() < b.first() ? a : b
94  : a.second() != -1 ? a : b;
95  }
96  );
97 
98  // Skip this point if it is not in the global mesh
99  if (procAndCelli.second() == -1)
100  {
101  ++ pointi;
102  }
103 
104  // If the point is in the global mesh then track to create a segment
105  else
106  {
108  (
109  particles,
110  points,
111  true,
112  false,
113  false,
114  samplingPositions,
115  samplingDistances,
116  samplingCells,
117  samplingFaces
118  );
119 
120  // Clear the cloud, then, if the point is in this local mesh,
121  // initialise a particle at the point
122  particles.clear();
123  if (procAndCelli.first() == Pstream::myProcNo())
124  {
125  particles.addParticle
126  (
128  (
129  mesh,
130  points[pointi],
131  procAndCelli.second(),
132  pointi,
133  1,
134  distance
135  )
136  );
137 
138  particles.first()->store(particles, td);
139  }
140 
141  // Track to create this segment
142  particles.move(particles, td);
143 
144  // Set the segment indices
145  samplingSegments.append
146  (
147  labelList
148  (
149  samplingPositions.size() - samplingSegments.size(),
150  segmenti
151  )
152  );
153 
154  // Move on to the next segment
155  ++ segmenti;
156 
157  // Determine the global number of samples completed
158  const label samplei0 = samplei;
159  samplei = returnReduce(samplingPositions.size(), sumOp<label>());
160 
161  // Move to the next unsampled point
162  pointi += samplei - samplei0;
163  }
164  }
165 }
166 
167 
168 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
169 
170 void Foam::sampledSets::points::calcSamplesUnordered
171 (
172  DynamicList<point>& samplingPositions,
173  DynamicList<label>& samplingSegments,
174  DynamicList<label>& samplingCells,
175  DynamicList<label>& samplingFaces
176 ) const
177 {
178  forAll(points_, i)
179  {
180  const point& pt = points_[i];
181  const label celli = searchEngine().findCell(pt);
182 
183  if (celli != -1)
184  {
185  samplingPositions.append(pt);
186  samplingSegments.append(i);
187  samplingCells.append(celli);
188  samplingFaces.append(-1);
189  }
190  }
191 }
192 
193 
194 void Foam::sampledSets::points::calcSamplesOrdered
195 (
196  DynamicList<point>& samplingPositions,
197  DynamicList<scalar>& samplingDistances,
198  DynamicList<label>& samplingSegments,
199  DynamicList<label>& samplingCells,
200  DynamicList<label>& samplingFaces
201 ) const
202 {
203  // Calculate the sampling topology
204  calcSamples
205  (
206  mesh(),
207  searchEngine(),
208  pointField(points_),
209  samplingPositions,
210  samplingDistances,
211  samplingSegments,
212  samplingCells,
213  samplingFaces
214  );
215 }
216 
217 
218 void Foam::sampledSets::points::genSamples()
219 {
220  DynamicList<point> samplingPositions;
221  DynamicList<scalar> samplingDistances;
222  DynamicList<label> samplingSegments;
223  DynamicList<label> samplingCells;
224  DynamicList<label> samplingFaces;
225 
226  if (!ordered_)
227  {
228  calcSamplesUnordered
229  (
230  samplingPositions,
231  samplingSegments,
232  samplingCells,
233  samplingFaces
234  );
235  }
236  else
237  {
238  calcSamplesOrdered
239  (
240  samplingPositions,
241  samplingDistances,
242  samplingSegments,
243  samplingCells,
244  samplingFaces
245  );
246  }
247 
248  samplingPositions.shrink();
249  samplingDistances.shrink();
250  samplingSegments.shrink();
251  samplingCells.shrink();
252  samplingFaces.shrink();
253 
254  if (!ordered_)
255  {
256  setSamples
257  (
258  samplingPositions,
259  samplingSegments,
260  samplingCells,
261  samplingFaces
262  );
263  }
264  else
265  {
266  setSamples
267  (
268  samplingPositions,
269  samplingDistances,
270  samplingSegments,
271  samplingCells,
272  samplingFaces
273  );
274  }
275 }
276 
277 
278 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
279 
281 (
282  const word& name,
283  const polyMesh& mesh,
284  const meshSearch& searchEngine,
285  const dictionary& dict
286 )
287 :
288  sampledSet(name, mesh, searchEngine, dict),
289  points_(dict.lookup("points")),
290  ordered_(dict.lookup("ordered"))
291 {
292  genSamples();
293 }
294 
295 
296 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
297 
299 {}
300 
301 
302 // ************************************************************************* //
#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:231
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:218
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:275
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:160
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:107
virtual ~points()
Destructor.
Definition: points.C:298
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:281
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:27
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
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dictionary dict