Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::lumpedPointState
28
29Description
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
54SeeAlso
55 Foam::coordinateRotations::euler, Foam::quaternion
56
57SourceFiles
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
77namespace Foam
78{
79
80/*---------------------------------------------------------------------------*\
81 Class lumpedPointState Declaration
82\*---------------------------------------------------------------------------*/
83
85{
86public:
87
88 // Data Types
89
90 //- Input format types
91 enum class inputFormatType
92 {
93 PLAIN,
95 };
96
97
98 // Static Data
99
100 //- Names for the input format types
101 static const Enum<inputFormatType> formatNames;
102
103
104private:
105
106 // Private Data
107
108 //- Positions of lumped points
109 pointField points_;
110
111 //- Orientation of lumped points (as Euler angles)
112 vectorField angles_;
114 //- The Euler-angle rotation order (default: zxz)
117 //- Euler angles measured in degrees
118 bool degrees_;
119
120 //- Tensor rotation of lumped points
121 mutable std::unique_ptr<tensorField> rotationPtr_;
122
123
124 // Private Member Functions
126 void calcRotations() const;
127
128 void readDict
129 (
130 const dictionary& dict,
132 const bool degrees = false
133 );
134
135public:
136
137 // Public Data
138
139 //- Enable/disable visualization of unused points
140 static bool visUnused;
141
142 //- The length for visualization triangles
143 static scalar visLength;
144
145
146 // Constructors
147
148 //- Default construct
150
151 //- Copy construct
153
154 //- Copy construct from points and angles
156 (
157 const pointField& pts,
158 const vectorField& ang,
160 const bool degrees = false
161 );
162
163 //- Copy construct from points with zero-rotation
164 explicit lumpedPointState
165 (
166 const pointField& pts,
168 const bool degrees = false
169 );
170
171 //- Construct from points with zero-rotation
172 explicit lumpedPointState
173 (
176 const bool degrees = false
177 );
178
179 //- Construct from dictionary
180 explicit lumpedPointState
182 const dictionary& dict,
184 const bool degrees = false
185 );
186
187
188 //- Destructor
189 virtual ~lumpedPointState() = default;
190
191
192 // Member Functions
193
194 //- Has positions and consistent number of rotations?
195 inline bool good() const;
196
197 //- If no points were specified
198 inline bool empty() const;
199
200 //- The number of points
201 inline label size() const;
202
203 //- Same as good()
204 bool valid() const { return good(); }
205
206 //- The points corresponding to mass centres
207 inline const pointField& points() const;
208
209 //- The orientation of the points (mass centres)
210 inline const vectorField& angles() const;
211
212 //- The local-to-global transformation for each point
213 inline const tensorField& rotations() const;
214
215 //- Scale points by given factor.
216 // Zero and negative values are ignored.
217 void scalePoints(const scalar scaleFactor);
218
219 //- The Euler-angle rotation order
221
222 //- Rotation angles in degrees
223 inline bool degrees() const;
224
225 //- Relax the state
226 // alpha = 1 : no relaxation
227 // alpha < 1 : relaxation
228 // alpha = 0 : do nothing
229 void relax(const scalar alpha, const lumpedPointState& prev);
230
231 //- Read input as dictionary content
232 bool readData
233 (
234 Istream& is,
236 const bool degrees = false
237 );
238
239 //- Output as dictionary content
240 bool writeData(Ostream& os) const;
242 //- Output as dictionary content
243 void writeDict(Ostream& os) const;
244
245 //- Read input as plain content
246 bool readPlain
247 (
248 Istream& is,
250 const bool degrees = false
251 );
252
253 //- Output as plain content
254 void writePlain(Ostream& os) const;
255
256 //- Read input from file (master process only) using specified format
257 bool readData
258 (
259 const inputFormatType fmt,
260 const fileName& file,
262 const bool degrees = false
263 );
265 //- Output points/rotations as VTK file for debugging/visualization
266 // The points are written as vertices, rotation as a triangle
267 void writeVTP
268 (
269 const fileName& outputFile,
270 const labelListList& lines = labelListList(),
271 const labelList& pointIds = labelList::null()
272 ) const;
273
274
275 // Member Operators
276
277 //- Assignment operator
278 void operator=(const lumpedPointState& rhs);
279
280 //- Shift points by specified origin
281 void operator+=(const point& origin);
282};
283
284
285// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286
287} // End namespace Foam
288
289// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290
291#include "lumpedPointStateI.H"
292
293#endif
294
295// ************************************************************************* //
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
static const List< label > & null() noexcept
Definition List.H:138
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
The state of lumped points corresponds to positions and rotations.
const tensorField & rotations() const
The local-to-global transformation for each point.
lumpedPointState()
Default construct.
void writePlain(Ostream &os) const
Output as plain content.
inputFormatType
Input format types.
@ DICTIONARY
"dictionary" is the OpenFOAM dictionary format
@ PLAIN
"plain" is a simple ASCII format
bool degrees() const
Rotation angles in degrees.
bool valid() const
Same as good().
label size() const
The number of points.
bool writeData(Ostream &os) const
Output as dictionary content.
static bool visUnused
Enable/disable visualization of unused points.
bool empty() const
If no points were specified.
quaternion::eulerOrder rotationOrder() const
The Euler-angle rotation order.
static scalar visLength
The length for visualization triangles.
static const Enum< inputFormatType > formatNames
Names for the input format types.
bool readData(Istream &is, const quaternion::eulerOrder rotOrder=quaternion::eulerOrder::ZXZ, const bool degrees=false)
Read input as dictionary content.
void operator=(const lumpedPointState &rhs)
Assignment operator.
void writeDict(Ostream &os) const
Output as dictionary content.
bool good() const
Has positions and consistent number of rotations?
void writeVTP(const fileName &outputFile, const labelListList &lines=labelListList(), const labelList &pointIds=labelList::null()) const
Output points/rotations as VTK file for debugging/visualization.
const pointField & points() const
The points corresponding to mass centres.
bool readPlain(Istream &is, const quaternion::eulerOrder rotOrder=quaternion::eulerOrder::ZXZ, const bool degrees=false)
Read input as plain content.
void operator+=(const point &origin)
Shift points by specified origin.
virtual ~lumpedPointState()=default
Destructor.
void scalePoints(const scalar scaleFactor)
Scale points by given factor.
const vectorField & angles() const
The orientation of the points (mass centres).
eulerOrder
Euler-angle rotation order.
Definition quaternion.H:116
A class for managing temporary objects.
Definition tmp.H:75
UEqn relax()
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
vectorField pointField
pointField is a vectorField.
volScalarField & alpha
dictionary dict
const pointField & pts