lumpedPointState.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) 2016-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::lumpedPointState
28 
29 Description
30  The \a state of lumped points corresponds to positions and rotations.
31 
32  This class encapsulates the response from the external application and
33  serves as the entry point for applying relaxation, sub-stepping etc.
34 
35  \heading Dictionary input format
36  \table
37  Property | Description | Required | Default
38  points | List of points | yes |
39  angles | List of Euler rotation angles | yes |
40  rotationOrder | The Euler-angle rotation order | no | zxz
41  degrees | Rotation angles in degrees | no | false
42  \endtable
43 
44  \heading Plain input format.
45  Blank and comment lines starting with a '#' character are ignored.
46  The angles are always in radians.
47  \verbatim
48  NumPoints
49  x0 y0 z0 eulerz0 eulerx'0 eulerz''0
50  x1 y1 z1 eulerz1 eulerx'1 eulerz''1
51  ...
52  \endverbatim
53 
54 SeeAlso
55  Foam::coordinateRotations::euler, Foam::quaternion
56 
57 SourceFiles
58  lumpedPointState.C
59  lumpedPointStateI.H
60 
61 \*---------------------------------------------------------------------------*/
62 
63 #ifndef lumpedPointState_H
64 #define lumpedPointState_H
65 
66 #include "dictionary.H"
67 #include "scalarList.H"
68 #include "pointField.H"
69 #include "scalarField.H"
70 #include "vectorField.H"
71 #include "tensorField.H"
72 #include "quaternion.H"
73 #include "Enum.H"
74 
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 
77 namespace Foam
78 {
79 
80 // Forward Declarations
81 class Istream;
82 class Ostream;
83 
84 /*---------------------------------------------------------------------------*\
85  Class lumpedPointState Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 class lumpedPointState
89 {
90 public:
91 
92  // Data Types
93 
94  //- Input format types
95  enum class inputFormatType
96  {
97  PLAIN,
98  DICTIONARY
99  };
100 
101 
102  // Static Data
103 
104  //- Names for the input format types
105  static const Enum<inputFormatType> formatNames;
106 
107 
108 private:
109 
110  // Private Data
111 
112  //- Positions of lumped points
113  pointField points_;
114 
115  //- Orientation of lumped points (as Euler angles)
116  vectorField angles_;
117 
118  //- The Euler-angle rotation order (default: zxz)
120 
121  //- Euler angles measured in degrees
122  bool degrees_;
123 
124  //- Tensor rotation of lumped points
125  mutable unique_ptr<tensorField> rotationPtr_;
126 
127 
128  // Private Member Functions
129 
130  void calcRotations() const;
131 
132  void readDict
133  (
134  const dictionary& dict,
135  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
136  const bool degrees = false
137  );
138 
139 public:
140 
141  // Public Data
142 
143  //- Enable/disable visualization of unused points
144  static bool visUnused;
145 
146  //- The length for visualization triangles
147  static scalar visLength;
148 
149 
150  // Constructors
151 
152  //- Default construct
154 
155  //- Copy construct
157 
158  //- Copy construct from points and angles
160  (
161  const pointField& pts,
162  const vectorField& ang,
163  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
164  const bool degrees = false
165  );
166 
167  //- Copy construct from points with zero-rotation
169  (
170  const pointField& pts,
171  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
172  const bool degrees = false
173  );
174 
175  //- Construct from points with zero-rotation
176  explicit lumpedPointState
177  (
178  tmp<pointField>& pts,
179  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
180  const bool degrees = false
181  );
182 
183  //- Construct from dictionary
184  explicit lumpedPointState
185  (
186  const dictionary& dict,
187  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
188  const bool degrees = false
189  );
190 
191 
192  //- Destructor
193  virtual ~lumpedPointState() = default;
194 
195 
196  // Member Functions
197 
198  //- Has positions and consistent number of rotations?
199  inline bool valid() const;
200 
201  //- If no points were specified
202  inline bool empty() const;
203 
204  //- The number of points
205  inline label size() const;
206 
207  //- The points corresponding to mass centres
208  inline const pointField& points() const;
209 
210  //- The orientation of the points (mass centres)
211  inline const vectorField& angles() const;
212 
213  //- The local-to-global transformation for each point
214  inline const tensorField& rotations() const;
215 
216  //- Scale points by given factor.
217  // Zero and negative values are ignored.
218  void scalePoints(const scalar scaleFactor);
219 
220  //- The Euler-angle rotation order
221  inline quaternion::eulerOrder rotationOrder() const;
222 
223  //- Rotation angles in degrees
224  inline bool degrees() const;
225 
226  //- Relax the state
227  // alpha = 1 : no relaxation
228  // alpha < 1 : relaxation
229  // alpha = 0 : do nothing
230  void relax(const scalar alpha, const lumpedPointState& prev);
231 
232  //- Read input as dictionary content
233  bool readData
234  (
235  Istream& is,
236  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
237  const bool degrees = false
238  );
239 
240  //- Output as dictionary content
241  bool writeData(Ostream& os) const;
242 
243  //- Output as dictionary content
244  void writeDict(Ostream& os) const;
245 
246  //- Read input as plain content
247  bool readPlain
248  (
249  Istream& is,
250  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
251  const bool degrees = false
252  );
253 
254  //- Output as plain content
255  void writePlain(Ostream& os) const;
256 
257  //- Read input from file (master process only) using specified format
258  bool readData
259  (
260  const inputFormatType fmt,
261  const fileName& file,
262  const quaternion::eulerOrder rotOrder = quaternion::eulerOrder::ZXZ,
263  const bool degrees = false
264  );
265 
266  //- Output points/rotations as VTK file for debugging/visualization
267  // The points are written as vertices, rotation as a triangle
268  void writeVTP
269  (
270  const fileName& outputFile,
271  const labelListList& lines = labelListList(),
272  const labelList& pointIds = labelList::null()
273  ) const;
274 
275 
276  // Member Operators
277 
278  //- Assignment operator
279  void operator=(const lumpedPointState& rhs);
280 
281  //- Shift points by specified origin
282  void operator+=(const point& origin);
283 };
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace Foam
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #include "lumpedPointStateI.H"
293 
294 #endif
295 
296 // ************************************************************************* //
Foam::Enum< inputFormatType >
Foam::List::null
static const List< T > & null()
Return a null List.
Definition: ListI.H:109
Foam::lumpedPointState::operator+=
void operator+=(const point &origin)
Shift points by specified origin.
Definition: lumpedPointState.C:227
Foam::lumpedPointState::writeDict
void writeDict(Ostream &os) const
Output as dictionary content.
Definition: lumpedPointState.C:340
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
scalarField.H
Foam::lumpedPointState::inputFormatType
inputFormatType
Input format types.
Definition: lumpedPointState.H:119
quaternion.H
Foam::lumpedPointState::operator=
void operator=(const lumpedPointState &rhs)
Assignment operator.
Definition: lumpedPointState.C:216
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::lumpedPointState::relax
void relax(const scalar alpha, const lumpedPointState &prev)
Relax the state.
Definition: lumpedPointState.C:246
Foam::lumpedPointState::visUnused
static bool visUnused
Enable/disable visualization of unused points.
Definition: lumpedPointState.H:168
Foam::lumpedPointState::readPlain
bool readPlain(Istream &is, const quaternion::eulerOrder rotOrder=quaternion::eulerOrder::ZXZ, const bool degrees=false)
Read input as plain content.
Definition: lumpedPointState.C:275
Foam::lumpedPointState::points
const pointField & points() const
The points corresponding to mass centres.
Definition: lumpedPointStateI.H:46
scalarList.H
Foam::lumpedPointState::readData
bool readData(Istream &is, const quaternion::eulerOrder rotOrder=quaternion::eulerOrder::ZXZ, const bool degrees=false)
Read input as dictionary content.
Definition: lumpedPointState.C:320
Foam::lumpedPointState::inputFormatType::PLAIN
"plain" is a simple ASCII format
Foam::lumpedPointState
The state of lumped points corresponds to positions and rotations.
Definition: lumpedPointState.H:112
Foam::lumpedPointState::scalePoints
void scalePoints(const scalar scaleFactor)
Scale points by given factor.
Definition: lumpedPointState.C:236
Foam::Field< vector >
Foam::lumpedPointState::formatNames
static const Enum< inputFormatType > formatNames
Names for the input format types.
Definition: lumpedPointState.H:129
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::lumpedPointState::angles
const vectorField & angles() const
The orientation of the points (mass centres)
Definition: lumpedPointStateI.H:52
Foam::lumpedPointState::writePlain
void writePlain(Ostream &os) const
Output as plain content.
Definition: lumpedPointState.C:355
Foam::lumpedPointState::rotationOrder
quaternion::eulerOrder rotationOrder() const
The Euler-angle rotation order.
Definition: lumpedPointStateI.H:70
Foam::lumpedPointState::writeVTP
void writeVTP(const fileName &outputFile, const labelListList &lines=labelListList(), const labelList &pointIds=labelList::null()) const
Output points/rotations as VTK file for debugging/visualization.
Definition: lumpedPointStateWriter.C:36
Foam::quaternion::eulerOrder
eulerOrder
Euler-angle rotation order.
Definition: quaternion.H:102
Foam::lumpedPointState::empty
bool empty() const
If no points were specified.
Definition: lumpedPointStateI.H:34
Foam::lumpedPointState::valid
bool valid() const
Has positions and consistent number of rotations?
Definition: lumpedPointStateI.H:28
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::lumpedPointState::writeData
bool writeData(Ostream &os) const
Output as dictionary content.
Definition: lumpedPointState.C:333
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
tensorField.H
pointField.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::Vector< scalar >
Foam::List< labelList >
vectorField.H
dictionary.H
Foam::lumpedPointState::rotations
const tensorField & rotations() const
The local-to-global transformation for each point.
Definition: lumpedPointStateI.H:58
Foam::lumpedPointState::~lumpedPointState
virtual ~lumpedPointState()=default
Destructor.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::lumpedPointState::lumpedPointState
lumpedPointState()
Default construct.
Definition: lumpedPointState.C:115
Foam::lumpedPointState::size
label size() const
The number of points.
Definition: lumpedPointStateI.H:40
Foam::lumpedPointState::degrees
bool degrees() const
Rotation angles in degrees.
Definition: lumpedPointStateI.H:76
Foam::lumpedPointState::visLength
static scalar visLength
The length for visualization triangles.
Definition: lumpedPointState.H:171
lumpedPointStateI.H
Foam::lumpedPointState::inputFormatType::DICTIONARY
"dictionary" is the OpenFOAM dictionary format
Enum.H