coordSet.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "coordSet.H"
30 #include "globalIndex.H"
31 
32 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
33 
34 const Foam::Enum
35 <
37 >
39 ({
40  { coordFormat::X, "x" },
41  { coordFormat::Y, "y" },
42  { coordFormat::Z, "z" },
43  { coordFormat::RADIUS, "radius" },
44  { coordFormat::DISTANCE, "distance" },
45  { coordFormat::XYZ, "xyz" },
47 });
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 {
54  if (points().size() != distance().size())
55  {
57  << "Size not equal :" << nl
58  << " points:" << points().size()
59  << " distance:" << distance().size()
60  << abort(FatalError);
61  }
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
67 // Foam::coordSet::coordSet()
68 // :
69 // pointField(),
70 // name_(),
71 // distance_(),
72 // axis_(coordFormat::DEFAULT)
73 // {}
74 
75 
77 (
78  const word& name,
79  const coordFormat axisType
80 )
81 :
82  pointField(),
83  name_(name),
84  distance_(),
85  axis_(axisType)
86 {}
87 
88 
90 (
91  const word& name,
92  const word& axis
93 )
94 :
95  coordSet(name, coordFormatNames.get(axis))
96 {}
97 
98 
100 (
101  const word& name,
102  const word& axis,
103  const List<point>& points,
104  const scalarList& dist
105 )
106 :
108  name_(name),
109  distance_(dist),
110  axis_(coordFormatNames[axis])
111 {
112  checkDimensions();
113 }
114 
115 
117 (
118  const word& name,
119  const word& axis,
121  scalarList&& dist
122 )
123 :
124  pointField(std::move(points)),
125  name_(name),
126  distance_(std::move(dist)),
127  axis_(coordFormatNames.get(axis))
128 {
129  checkDimensions();
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 {
137  return axis_ == coordFormat::XYZ;
138 }
139 
140 
141 Foam::scalar Foam::coordSet::scalarCoord(const label index) const
142 {
143  switch (axis_)
144  {
145  case coordFormat::X:
146  {
147  return points()[index].x();
148  }
149  case coordFormat::Y:
150  {
151  return points()[index].y();
152  }
153  case coordFormat::Z:
154  {
155  return points()[index].z();
156  }
157  case coordFormat::RADIUS:
158  {
159  return mag(points()[index]);
160  }
161  case coordFormat::DISTANCE:
162  {
163  // Note: the distance will unset it constructed from
164  // 'name' and 'axis' only
165 
166  if (distance().empty())
167  {
169  << "Axis type '" << coordFormatNames[axis_]
170  << "' requested but curve distance has not been set"
171  << abort(FatalError);
172  }
173  return distance()[index];
174  }
175  default:
176  {
178  << "Illegal axis specification '" << coordFormatNames[axis_]
179  << "' for sampling " << name_
180  << exit(FatalError);
181  }
182  }
183 
184  return 0;
185 }
186 
188 const Foam::point& Foam::coordSet::vectorCoord(const label index) const
189 {
190  return points()[index];
191 }
192 
193 
195 {
196  os << "name:" << name_ << " axis:" << coordFormatNames[axis_] << nl
197  << nl
198  << "\t(coord)" << nl;
199 
200  for (const point& p : *this)
201  {
202  os << '\t' << p << nl;
203  }
204 
205  return os;
206 }
207 
208 
210 (
211  labelList& sortOrder
212 ) const
213 {
214  // Combine sampleSet from processors. Sort by distance.
215  // Return ordering in indexSet.
216  // Note: only master results are valid
217 
218  List<point> allPoints(globalIndex::gatherOp(points()));
219  List<scalar> allDistance(globalIndex::gatherOp(distance()));
220 
221  if (Pstream::master() && allDistance.empty())
222  {
224  << "Gathered empty coordSet: " << name() << endl;
225  }
226 
227  // Sort according to distance
228  Foam::sortedOrder(allDistance, sortOrder); // uses stable sort
229 
230  // Repopulate gathered points/distances in the correct order
231  allPoints = List<point>(allPoints, sortOrder);
232  allDistance = List<scalar>(allDistance, sortOrder);
233 
235  (
236  name(),
237  axis(),
238  std::move(allPoints),
239  std::move(allDistance)
240  );
241 }
242 
243 
244 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const scalarList & distance() const noexcept
Return the cumulative distance.
Definition: coordSet.H:176
const vector & vectorCoord(const label index) const
Get point according to axis="xyz" specification.
Definition: coordSet.C:181
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
autoPtr< coordSet > gatherSort(labelList &sortOrder) const
Gather and sort.
Definition: coordSet.C:203
bool hasVectorAxis() const noexcept
True if axis specification is a vector.
Definition: coordSet.C:128
Holds list of sampling positions.
Definition: coordSet.H:49
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
void checkDimensions() const
Check for consistent dimensions of points and curve distance.
Definition: coordSet.C:45
A class for handling words, derived from Foam::string.
Definition: word.H:63
Ostream & write(Ostream &os) const
Write to stream.
Definition: coordSet.C:187
const pointField & points() const noexcept
Return the points.
Definition: coordSet.H:168
coordFormat
Enumeration defining the output format for coordinates.
Definition: coordSet.H:60
errorManip< error > abort(error &err)
Definition: errorManip.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
OBJstream os(runTime.globalPath()/outputName)
coordSet(const word &name, const coordFormat axisType)
Default construct with name and axis type.
Definition: coordSet.C:70
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
PtrList< volScalarField > & Y
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
static void gatherOp(const UList< Type > &sendData, List< Type > &allData, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking, const label comm=UPstream::worldComm)
Collect data in processor order on master (in serial: performs a simple copy).
scalar scalarCoord(const label index) const
Get coordinate of point according to axis specification.
Definition: coordSet.C:134
static const Enum< coordFormat > coordFormatNames
String representation of coordFormat enum.
Definition: coordSet.H:74