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