coordSet.H
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-2012 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 Class
28  Foam::coordSet
29 
30 Description
31  Holds list of sampling positions
32 
33 SourceFiles
34  coordSet.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef Foam_coordSet_H
39 #define Foam_coordSet_H
40 
41 #include "pointField.H"
42 #include "scalarList.H"
43 #include "Enum.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class coordSet Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class coordSet
55 :
56  public pointField
57 {
58 public:
59 
60  // Public Data Types
61 
62  //- Enumeration defining the output format for coordinates
63  enum class coordFormat
64  {
65  X,
66  Y,
67  Z,
68  RADIUS,
69  DISTANCE,
70  XYZ,
71  };
73 
74  //- String representation of coordFormat enum
76 
77 
78 protected:
79 
80  // Protected Data
81 
82  //- Name
83  word name_;
84 
85  //- Cumulative distance for "distance" write specifier.
87 
88  //- Axis type
90 
91 
92  // Protected Member Functions
93 
94  //- Check for consistent dimensions of points and curve distance
95  void checkDimensions() const;
96 
97 
98 public:
99 
100  // Constructors
101 
104 
105  //- Default construct with name and axis type
106  coordSet(const word& name, const coordFormat axisType);
107 
108  //- Default construct with name and axis type
109  coordSet(const word& name, const word& axis);
110 
111  //- Copy construct from components
112  coordSet
113  (
114  const word& name,
115  const word& axis,
116  const List<point>& points,
117  const scalarList& dist
118  );
119 
120  //- Move construct from components
121  coordSet
122  (
123  const word& name,
124  const word& axis,
126  scalarList&& dist
127  );
128 
129 
130  // Member Functions
131 
132  // Access
133 
134  //- The coord-set name
135  const word& name() const noexcept
136  {
137  return name_;
138  }
139 
140  //- The sort axis name
141  const word& axis() const
142  {
143  return coordFormatNames[axis_];
144  }
145 
146  //- Return the points
147  const pointField& points() const noexcept
148  {
149  return static_cast<const pointField&>(*this);
150  }
151 
152  //- Return the cumulative distance
153  const scalarList& distance() const noexcept
154  {
155  return distance_;
156  }
157 
158  //- Return the number of points
159  label nPoints() const noexcept
160  {
161  return pointField::size();
162  }
163 
164 
165  // Edit
166 
167  //- Rename the coordinate set
168  void rename(const word& newName)
169  {
170  name_ = newName;
171  }
172 
173  //- Copy assign new points
174  void setPoints(const List<point>& newPoints)
175  {
176  static_cast<pointField&>(*this) = newPoints;
177  }
178 
179  //- Move assign new points
180  void setPoints(List<point>&& newPoints)
181  {
182  static_cast<pointField&>(*this) = std::move(newPoints);
183  }
185  //- Copy assign the cumulative distance
186  void setDistance(const scalarList& dist, const bool check=true)
187  {
188  distance_ = dist;
189  if (check) checkDimensions();
190  }
191 
192  //- Move assign the cumulative distance
193  void setDistance(scalarList&& dist, const bool check=true)
194  {
195  distance_ = std::move(dist);
196  if (check) checkDimensions();
197  }
198 
199 
200  // Output-related
201 
202  //- True if axis specification is a vector
203  bool hasVectorAxis() const noexcept;
204 
205  //- Get coordinate of point according to axis specification.
206  // If axis="distance" is the distance[index]
207  scalar scalarCoord(const label index) const;
208 
209  //- Get point according to axis="xyz" specification
210  const vector& vectorCoord(const label index) const;
212  //- Write to stream
213  Ostream& write(Ostream& os) const;
214 
215 
216  // Other
217 
218  //- Gather and sort.
219  // \return (on master) gathered set and overall sort order
220  autoPtr<coordSet> gatherSort(labelList& sortOrder) const;
221 
222 
223  // Housekeeping
224 
225  //- Return the cumulative distance
226  const scalarList& curveDist() const noexcept
227  {
228  return distance_;
229  }
230 
231  //- Copy assign the cumulative distance
232  void setCurveDist(const scalarList& dist)
233  {
234  setDistance(dist);
235  }
236 };
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 } // End namespace Foam
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 #endif
246 
247 // ************************************************************************* //
Use &#39;x&#39; component of points for (scalar) axis.
const word & axis() const
The sort axis name.
Definition: coordSet.H:160
Use &#39;y&#39; component of points for (scalar) axis.
Use &#39;z&#39; component of points for (scalar) axis.
const scalarList & curveDist() const noexcept
Return the cumulative distance.
Definition: coordSet.H:275
word name_
Name.
Definition: coordSet.H:84
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
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
void setDistance(const scalarList &dist, const bool check=true)
Copy assign the cumulative distance.
Definition: coordSet.H:219
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Holds list of sampling positions.
Definition: coordSet.H:49
label nPoints() const noexcept
Return the number of points.
Definition: coordSet.H:184
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
scalarList distance_
Cumulative distance for "distance" write specifier.
Definition: coordSet.H:89
const pointField & points() const noexcept
Return the points.
Definition: coordSet.H:168
void setCurveDist(const scalarList &dist)
Copy assign the cumulative distance.
Definition: coordSet.H:283
coordFormat
Enumeration defining the output format for coordinates.
Definition: coordSet.H:60
const word & name() const noexcept
The coord-set name.
Definition: coordSet.H:152
void rename(const word &newName)
Rename the coordinate set.
Definition: coordSet.H:195
coordFormat axis_
Axis type.
Definition: coordSet.H:94
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)
static void check(const int retVal, const char *what)
Use mag of points for (scalar) axis.
coordSet(const word &name, const coordFormat axisType)
Default construct with name and axis type.
Definition: coordSet.C:70
Use additional distance field for (scalar) axis.
void setPoints(const List< point > &newPoints)
Copy assign new points.
Definition: coordSet.H:203
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Namespace for OpenFOAM.
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