coordSet.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-2025 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 "coordSet.H"
27 #include "ListListOps.H"
28 
29 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
30 
33 {
34  "xyz",
35  "x",
36  "y",
37  "z",
38  "distance",
39  "default"
40 };
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 :
47  segments_(0),
48  positionName_(word::null),
49  positions_(nullptr),
50  distanceName_(word::null),
51  distances_(nullptr),
52  axis_(axisType::DEFAULT)
53 {}
54 
55 
57 (
58  const labelList& segments,
59  const word& positionName,
60  const pointField& positions,
61  const word& distanceName,
62  const scalarField& distances,
63  const word& axis
64 )
65 :
66  segments_(segments),
67  positionName_(positionName),
68  positions_
69  (
70  isNull<pointField>(positions)
71  ? nullptr
72  : new pointField(positions)
73  ),
74  distanceName_(distanceName),
75  distances_
76  (
77  isNull<scalarField>(distances)
78  ? nullptr
79  : new scalarField(distances)
80  ),
81  axis_(axisTypeNames_[axis])
82 {}
83 
84 
86 (
87  const bool contiguous,
88  const word& positionName,
89  const pointField& positions,
90  const word& axis
91 )
92 :
93  coordSet
94  (
96  ? labelList(positions.size(), 0)
97  : identityMap(positions.size()),
98  positionName,
99  positions,
100  word::null,
101  scalarField::null(),
102  axis
103  )
104 {}
105 
106 
108 (
109  const bool contiguous,
110  const word& distanceName,
111  const scalarField& distances,
112  const word& axis
113 )
114 :
115  coordSet
116  (
117  contiguous
118  ? labelList(distances.size(), 0)
119  : identityMap(distances.size()),
120  word::null,
121  pointField::null(),
122  distanceName,
123  distances,
124  axis
125  )
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
132 {
133  // If axis is default and both positions and distances are valid, take
134  // scalar distances in preference
135 
136  return
137  (axis_ == axisType::X && positions_.valid())
138  || (axis_ == axisType::Y && positions_.valid())
139  || (axis_ == axisType::Z && positions_.valid())
140  || (axis_ == axisType::DISTANCE && distances_.valid())
141  || (axis_ == axisType::DEFAULT && distances_.valid());
142 }
143 
144 
146 {
147  // If axis is default and both positions and distances are valid, take
148  // scalar distances in preference
149 
150  return
151  (axis_ == axisType::XYZ && positions_.valid())
152  || (axis_ == axisType::DEFAULT && positions_.valid());
153 }
154 
155 
156 Foam::scalar Foam::coordSet::scalarCoord(const label index) const
157 {
158  switch (axis_)
159  {
160  case axisType::XYZ:
162  << "Scalar coordinate requested from coordinate set with "
163  << axisTypeNames_[axis_] << " axis" << exit(FatalError);
164  break;
165  case axisType::X:
166  return positions_()[index].x();
167  case axisType::Y:
168  return positions_()[index].y();
169  case axisType::Z:
170  return positions_()[index].z();
171  case axisType::DISTANCE:
172  return distances_()[index];
173  case axisType::DEFAULT:
174  if (distances_.valid()) return distances_()[index];
176  << "Scalar coordinate requested from coordinate set with "
177  << axisTypeNames_[axis_] << " axis, but with no valid distances"
178  << exit(FatalError);
179  break;
180  }
181 
182  return NaN;
183 }
184 
185 
187 {
188  switch (axis_)
189  {
190  case axisType::XYZ:
192  << "Scalar coordinate requested from coordinate set with "
193  << axisTypeNames_[axis_] << " axis" << exit(FatalError);
194  break;
195  case axisType::X:
196  return positions_().component(point::X);
197  case axisType::Y:
198  return positions_().component(point::Y);
199  case axisType::Z:
200  return positions_().component(point::Z);
201  case axisType::DISTANCE:
202  return distances_();
203  case axisType::DEFAULT:
204  if (distances_.valid()) return distances_();
206  << "Scalar coordinate requested from coordinate set with "
207  << axisTypeNames_[axis_] << " axis, but with no valid distances"
208  << exit(FatalError);
209  break;
210  }
211 
212  return scalarField::null();
213 }
214 
215 
217 {
218  switch (axis_)
219  {
220  case axisType::XYZ:
222  << "Scalar name requested from coordinate set with "
223  << axisTypeNames_[axis_] << " axis" << exit(FatalError);
224  break;
225  case axisType::X:
226  case axisType::Y:
227  case axisType::Z:
228  return
229  positionName_
230  + (positionName_ != word::null ? "_" : "")
231  + axis();
232  case axisType::DISTANCE:
233  return distanceName_;
234  case axisType::DEFAULT:
235  if (distances_.valid()) return distanceName_;
237  << "Scalar coordinate requested from coordinate set with "
238  << axisTypeNames_[axis_] << " axis, but with no valid distances"
239  << exit(FatalError);
240  break;
241  }
242 
243  return word::null;
244 }
245 
246 
248 {
249  switch (axis_)
250  {
251  case axisType::XYZ:
252  return positions_()[index];
253  case axisType::X:
254  case axisType::Y:
255  case axisType::Z:
256  case axisType::DISTANCE:
258  << "Point coordinate requested from coordinate set with "
259  << axisTypeNames_[axis_] << " axis" << exit(FatalError);
260  break;
261  case axisType::DEFAULT:
262  if (positions_.valid()) return positions_()[index];
264  << "Point coordinate requested from coordinate set with "
265  << axisTypeNames_[axis_] << " axis, but with no valid point"
266  << exit(FatalError);
267  break;
268  }
269 
270  return point::uniform(NaN);
271 }
272 
273 
275 {
276  switch (axis_)
277  {
278  case axisType::XYZ:
279  return positions_();
280  case axisType::X:
281  case axisType::Y:
282  case axisType::Z:
283  case axisType::DISTANCE:
285  << "Point coordinate requested from coordinate set with "
286  << axisTypeNames_[axis_] << " axis" << exit(FatalError);
287  break;
288  case axisType::DEFAULT:
289  if (positions_.valid()) return positions_();
291  << "Point coordinate requested from coordinate set with "
292  << axisTypeNames_[axis_] << " axis, but with no valid points"
293  << exit(FatalError);
294  break;
295  }
296 
297  return pointField::null();
298 }
299 
300 
302 {
303  if (!positions_.valid())
304  {
306  << "Positions requested from coordinate set with no valid points"
307  << exit(FatalError);
308  }
309 
310  return positions_();
311 }
312 
313 
315 {
316  switch (axis_)
317  {
318  case axisType::XYZ:
319  return positionName_;
320  case axisType::X:
321  case axisType::Y:
322  case axisType::Z:
323  case axisType::DISTANCE:
325  << "Point name requested from coordinate set with "
326  << axisTypeNames_[axis_] << " axis" << exit(FatalError);
327  break;
328  case axisType::DEFAULT:
329  if (positions_.valid()) return positionName_;
331  << "Point name requested from coordinate set with "
332  << axisTypeNames_[axis_] << " axis, but with no valid point"
333  << exit(FatalError);
334  break;
335  }
336 
337  return word::null;
338 }
339 
340 
342 {
343  label nVertices = 0;
344  labelList vertices(size());
345 
346  forAll(*this, pointi)
347  {
348  const label s = segments_[pointi];
349  if
350  (
351  (pointi == 0 || s != segments_[pointi - 1])
352  && (pointi == size() - 1 || s != segments_[pointi + 1])
353  )
354  {
355  vertices[nVertices ++] = pointi;
356  }
357  }
358 
359  vertices.resize(nVertices);
360 
361  return vertices;
362 }
363 
364 
366 {
367  label nEdges = 0;
368  labelPairList edges(size() - 1);
369 
370  for (label pointi = 0; pointi < size() - 1; ++ pointi)
371  {
372  if (segments_[pointi] == segments_[pointi + 1])
373  {
374  edges[nEdges ++] = labelPair(pointi, pointi + 1);
375  }
376  }
377 
378  edges.resize(nEdges);
379 
380  return edges;
381 }
382 
383 
385 {
388 
389  forAll(*this, pointi)
390  {
391  line.append(pointi);
392 
393  if
394  (
395  pointi == size() - 1
396  || segments_[pointi] != segments_[pointi + 1]
397  )
398  {
399  if (line.size() > 1)
400  {
401  lines.append(line);
402  }
403 
404  line.clear();
405  }
406  }
407 
408  labelListList linesNonDynamic;
409  linesNonDynamic.transfer(lines);
410 
411  return linesNonDynamic;
412 }
413 
414 
416 {
417  // Collect data from all processors
418  List<List<point>> gatheredPositions;
419  if (positions_.valid())
420  {
421  gatheredPositions.resize(Pstream::nProcs());
422  gatheredPositions[Pstream::myProcNo()] = positions_();
423  Pstream::gatherList(gatheredPositions);
424  }
425 
426  List<scalarList> gatheredDistances;
427  if (distances_.valid())
428  {
429  gatheredDistances.resize(Pstream::nProcs());
430  gatheredDistances[Pstream::myProcNo()] = distances_();
431  Pstream::gatherList(gatheredDistances);
432  }
433 
434  List<labelField> gatheredSegments(Pstream::nProcs());
435  gatheredSegments[Pstream::myProcNo()] = segments_;
436  Pstream::gatherList(gatheredSegments);
437 
438  // Combine processor lists into one big list.
439  List<point> allPositions;
440  if (positions_.valid())
441  {
442  allPositions =
443  ListListOps::combine<List<point>>
444  (
445  gatheredPositions,
447  );
448  };
449 
450  scalarList allDistances;
451  if (distances_.valid())
452  {
453  allDistances =
454  ListListOps::combine<scalarList>
455  (
456  gatheredDistances,
458  );
459  }
460 
461  labelList allSegments
462  (
463  ListListOps::combine<labelList>
464  (
465  gatheredSegments,
467  )
468  );
469 
470  // Construct a result tuple
472  (
473  coordSet(),
474  identityMap(allSegments.size())
475  );
476  coordSet& set = result.first();
477  labelList& order = result.second();
478 
479  // Sort by segment then by distance
480  if (distances_.valid())
481  {
482  stableSort
483  (
484  order,
485  [&](const label a, const label b)
486  {
487  return
488  allSegments[a] < allSegments[b] ? true
489  : allSegments[a] > allSegments[b] ? false
490  : allDistances[a] < allDistances[b];
491  }
492  );
493  }
494  else
495  {
496  stableSort
497  (
498  order,
499  [&](const label a, const label b)
500  {
501  return allSegments[a] < allSegments[b];
502  }
503  );
504  }
505 
506  // Set the data in the coordinate set
507  set.segments_ = labelField(allSegments, order);
508  set.positionName_ = positionName_;
509  if (positions_.valid())
510  {
511  set.positions_.set(new pointField(allPositions, order));
512  }
513  set.distanceName_ = distanceName_;
514  if (distances_.valid())
515  {
516  set.distances_.set(new scalarField(allDistances, order));
517  }
518  set.axis_ = axis_;
519 
520  return result;
521 }
522 
523 
524 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
static const Field< scalar > & null()
Return a null field.
Definition: Field.H:112
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
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
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
Holds list of sampling positions.
Definition: coordSet.H:51
word distanceName_
Name of the distances.
Definition: coordSet.H:86
bool hasPointAxis() const
Is the coordinate axis a point?
Definition: coordSet.C:145
labelListList lines() const
Return a list of lines. These are lists of points which are in the.
Definition: coordSet.C:384
labelList vertices() const
Return a list of isolated vertices. These are the points that are.
Definition: coordSet.C:341
Tuple2< coordSet, labelList > gather() const
Combine coordinate sets onto the master. Return both the combined.
Definition: coordSet.C:415
labelPairList edges() const
Return a list of edges. These are adjacent pairs of points which.
Definition: coordSet.C:365
autoPtr< pointField > positions_
Point positions.
Definition: coordSet.H:83
static const NamedEnum< axisType, 6 > axisTypeNames_
String representation of axis enums.
Definition: coordSet.H:69
axisType
Enumeration defining the output format for coordinates.
Definition: coordSet.H:58
labelList segments_
Connected segments.
Definition: coordSet.H:77
word positionName_
Name of the positions.
Definition: coordSet.H:80
axisType axis_
Axis.
Definition: coordSet.H:92
coordSet()
Construct null.
Definition: coordSet.C:45
tmp< scalarField > scalarCoords() const
Get scalar coordinates (axis is x, y, z or distance)
Definition: coordSet.C:186
point pointCoord(const label index) const
Get vector coordinate (axis is xyz)
Definition: coordSet.C:247
word scalarName() const
Return the name of the scalar coordinates.
Definition: coordSet.C:216
autoPtr< scalarField > distances_
Scalar distances.
Definition: coordSet.H:89
word pointName() const
Return the name of the point coordinates.
Definition: coordSet.C:314
tmp< pointField > pointCoords() const
Get vector coordinate (axis is xyz)
Definition: coordSet.C:274
bool hasScalarAxis() const
Is the coordinate axis a scalar?
Definition: coordSet.C:131
scalar scalarCoord(const label index) const
Get scalar coordinate (axis is x, y, z or distance)
Definition: coordSet.C:156
const pointField & positions() const
Get the positions.
Definition: coordSet.C:301
A line primitive.
Definition: line.H:71
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
pointField vertices(const blockVertexList &bvl)
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
int order(const scalar s)
bool contiguous()
contiguous
Definition: contiguous.H:53
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
void stableSort(UList< T > &)
Definition: UList.C:129
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
PtrList< volScalarField > & Y