coordSet.H
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 Class
25  Foam::coordSet
26 
27 Description
28  Holds list of sampling positions
29 
30 SourceFiles
31  coordSet.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef coordSet_H
36 #define coordSet_H
37 
38 #include "pointField.H"
39 #include "word.H"
40 #include "labelPair.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class coordSet Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 class coordSet
52 {
53 public:
54 
55  // Public data types
56 
57  //- Enumeration defining the output format for coordinates
58  enum class axisType
59  {
60  XYZ,
61  X,
62  Y,
63  Z,
64  DISTANCE,
65  DEFAULT
66  };
67 
68 
69  //- String representation of axis enums
71 
72 
73 protected:
74 
75  // Protected Data
76 
77  //- Connected segments
79 
80  //- Name of the positions
82 
83  //- Point positions
85 
86  //- Name of the distances
88 
89  //- Scalar distances
91 
92  //- Axis
94 
95 
96 public:
97 
98  // Constructors
99 
100  //- Construct null
101  coordSet();
102 
103  //- Construct from components
104  coordSet
105  (
106  const labelList& segments,
107  const word& positionName = word::null,
108  const pointField& positions = pointField::null(),
109  const word& distanceName = word::null,
110  const scalarField& distances = scalarField::null(),
111  const word& axis = axisTypeNames_[axisType::DEFAULT]
112  );
113 
114  //- Construct from positions
115  coordSet
116  (
117  const bool contiguous,
118  const word& positionName,
119  const pointField& positions,
120  const word& axis = axisTypeNames_[axisType::DEFAULT]
121  );
122 
123  //- Construct from distances
124  coordSet
125  (
126  const bool contiguous,
127  const word& distanceName,
128  const scalarField& distances,
129  const word& axis = axisTypeNames_[axisType::DEFAULT]
130  );
131 
132 
133  // Member Functions
134 
135  //- Return the size
136  inline label size() const
137  {
138  return segments_.size();
139  }
140 
141  //- Return the segments
142  inline const labelList& segments() const
143  {
144  return segments_;
145  }
146 
147  //- Return the axis name
148  inline word axis() const
149  {
150  return axisTypeNames_[axis_];
151  }
152 
153  //- Is the coordinate axis a scalar?
154  bool hasScalarAxis() const;
155 
156  //- Is the coordinate axis a point?
157  bool hasPointAxis() const;
158 
159  //- Get scalar coordinate (axis is x, y, z or distance)
160  scalar scalarCoord(const label index) const;
161 
162  //- Get scalar coordinates (axis is x, y, z or distance)
164 
165  //- Return the name of the scalar coordinates
166  word scalarName() const;
167 
168  //- Get vector coordinate (axis is xyz)
169  point pointCoord(const label index) const;
170 
171  //- Get vector coordinate (axis is xyz)
173 
174  //- Return the name of the point coordinates
175  word pointName() const;
176 
177  //- Return a list of isolated vertices. These are the points that are
178  // not adjacent to any points in the same segment.
179  labelList vertices() const;
180 
181  //- Return a list of edges. These are adjacent pairs of points which
182  // are in the same segment.
183  labelPairList edges() const;
184 
185  //- Return a list of lines. These are lists of points which are in the
186  // same segment.
187  labelListList lines() const;
188 
189  //- Combine coordinate sets onto the master. Return both the combined
190  // coordinate set, and an ordering to be used for gathering associated
191  // fields
193 
194  //- Combine a field using the ordering obtained from the coordinate set
195  // gather operation
196  template<class Type>
197  static tmp<Field<Type>> gather
198  (
199  const Field<Type>& values,
200  const labelList& order
201  );
202 };
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 #ifdef NoRepository
212  #include "coordSetTemplates.C"
213 #endif
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #endif
218 
219 // ************************************************************************* //
labelList segments_
Connected segments.
Definition: coordSet.H:77
word axis() const
Return the axis name.
Definition: coordSet.H:147
bool hasScalarAxis() const
Is the coordinate axis a scalar?
Definition: coordSet.C:138
autoPtr< scalarField > distances_
Scalar distances.
Definition: coordSet.H:89
label size() const
Return the size.
Definition: coordSet.H:135
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:65
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
word positionName_
Name of the positions.
Definition: coordSet.H:80
tmp< pointField > pointCoords() const
Get vector coordinate (axis is xyz)
Definition: coordSet.C:281
autoPtr< pointField > positions_
Point positions.
Definition: coordSet.H:83
static const NamedEnum< axisType, 6 > axisTypeNames_
String representation of axis enums.
Definition: coordSet.H:69
word scalarName() const
Return the name of the scalar coordinates.
Definition: coordSet.C:223
const labelList & segments() const
Return the segments.
Definition: coordSet.H:141
axisType axis_
Axis.
Definition: coordSet.H:92
bool hasPointAxis() const
Is the coordinate axis a point?
Definition: coordSet.C:152
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 class for handling words, derived from string.
Definition: word.H:59
point pointCoord(const label index) const
Get vector coordinate (axis is xyz)
Definition: coordSet.C:254
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
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
bool contiguous()
contiguous
Definition: contiguous.H:53
labelPairList edges() const
Return a list of edges. These are adjacent pairs of points which.
Definition: coordSet.C:359
word pointName() const
Return the name of the point coordinates.
Definition: coordSet.C:308
static const Field< Type > & null()
Return a null field.
Definition: Field.H:111
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
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