sixDoFRigidBodyMotionI.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-2015 OpenFOAM Foundation
9  Copyright (C) 2016-2023 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 inline Foam::tensor Foam::sixDoFRigidBodyMotion::rotationTensorX
32 (
33  scalar phi
34 ) const
35 {
36  return tensor
37  (
38  1, 0, 0,
39  0, Foam::cos(phi), -Foam::sin(phi),
41  );
42 }
43 
44 
45 inline Foam::tensor Foam::sixDoFRigidBodyMotion::rotationTensorY
46 (
47  scalar phi
48 ) const
49 {
50  return tensor
51  (
52  Foam::cos(phi), 0, Foam::sin(phi),
53  0, 1, 0,
54  -Foam::sin(phi), 0, Foam::cos(phi)
55  );
56 }
57 
58 
59 inline Foam::tensor Foam::sixDoFRigidBodyMotion::rotationTensorZ
60 (
61  scalar phi
62 ) const
63 {
64  return tensor
65  (
66  Foam::cos(phi), -Foam::sin(phi), 0,
67  Foam::sin(phi), Foam::cos(phi), 0,
68  0, 0, 1
69  );
70 }
71 
72 
74 Foam::sixDoFRigidBodyMotion::rotate
75 (
76  const tensor& Q0,
77  const vector& pi0,
78  const scalar deltaT
79 ) const
80 {
81  Tuple2<tensor, vector> Qpi(Q0, pi0);
82  tensor& Q = Qpi.first();
83  vector& pi = Qpi.second();
84 
85  tensor R = rotationTensorX(0.5*deltaT*pi.x()/momentOfInertia_.xx());
86  pi = pi & R;
87  Q = Q & R;
88 
89  R = rotationTensorY(0.5*deltaT*pi.y()/momentOfInertia_.yy());
90  pi = pi & R;
91  Q = Q & R;
92 
93  R = rotationTensorZ(deltaT*pi.z()/momentOfInertia_.zz());
94  pi = pi & R;
95  Q = Q & R;
96 
97  R = rotationTensorY(0.5*deltaT*pi.y()/momentOfInertia_.yy());
98  pi = pi & R;
99  Q = Q & R;
100 
101  R = rotationTensorX(0.5*deltaT*pi.x()/momentOfInertia_.xx());
102  pi = pi & R;
103  Q = Q & R;
104 
105  return Qpi;
106 }
107 
108 
110 Foam::sixDoFRigidBodyMotion::restraints() const
111 {
112  return restraints_;
113 }
114 
115 
117 Foam::sixDoFRigidBodyMotion::constraints() const
118 {
119  return constraints_;
120 }
121 
122 
123 inline const Foam::point&
124 Foam::sixDoFRigidBodyMotion::initialCentreOfRotation() const
125 {
126  return initialCentreOfRotation_;
127 }
128 
129 
130 inline const Foam::tensor&
131 Foam::sixDoFRigidBodyMotion::initialQ() const
132 {
133  return initialQ_;
134 }
135 
137 inline const Foam::tensor& Foam::sixDoFRigidBodyMotion::Q() const
138 {
139  return motionState_.Q();
140 }
141 
142 
143 inline const Foam::vector& Foam::sixDoFRigidBodyMotion::v() const
144 {
145  return motionState_.v();
146 }
147 
148 
149 inline const Foam::vector& Foam::sixDoFRigidBodyMotion::a() const
150 {
151  return motionState_.a();
152 }
153 
154 
155 inline const Foam::vector& Foam::sixDoFRigidBodyMotion::pi() const
156 {
157  return motionState_.pi();
158 }
159 
160 
161 inline const Foam::vector& Foam::sixDoFRigidBodyMotion::tau() const
162 {
163  return motionState_.tau();
164 }
165 
166 
167 inline Foam::point& Foam::sixDoFRigidBodyMotion::initialCentreOfRotation()
168 {
169  return initialCentreOfRotation_;
170 }
171 
172 
173 inline Foam::tensor& Foam::sixDoFRigidBodyMotion::initialQ()
174 {
175  return initialQ_;
176 }
177 
179 inline Foam::tensor& Foam::sixDoFRigidBodyMotion::Q()
180 {
181  return motionState_.Q();
182 }
183 
184 
186 {
187  return motionState_.v();
188 }
189 
190 
191 inline Foam::vector& Foam::sixDoFRigidBodyMotion::a()
192 {
193  return motionState_.a();
194 }
195 
196 
197 inline Foam::vector& Foam::sixDoFRigidBodyMotion::pi()
198 {
199  return motionState_.pi();
200 }
201 
202 
203 inline Foam::vector& Foam::sixDoFRigidBodyMotion::tau()
204 {
205  return motionState_.tau();
206 }
207 
208 
209 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
210 
211 inline Foam::scalar Foam::sixDoFRigidBodyMotion::mass() const
212 {
213  return mass_;
214 }
215 
216 
217 inline const Foam::diagTensor&
219 {
220  return momentOfInertia_;
221 }
222 
223 
226 {
227  return motionState_;
228 }
229 
230 
232 {
233  return motionState_.centreOfRotation();
234 }
235 
236 
237 inline const Foam::point&
239 {
240  return initialCentreOfMass_;
241 }
242 
245 {
246  return transform(initialCentreOfMass_);
247 }
248 
249 
251 {
252  return centreOfMass() - motionState_.centreOfRotation();
253 }
254 
255 
256 inline const Foam::tensor&
258 {
259  return Q();
260 }
262 
264 {
265  return Q() & (inv(momentOfInertia_) & pi());
266 }
267 
268 inline const Foam::Time& Foam::sixDoFRigidBodyMotion::time() const
269 {
270  return time_;
271 }
272 
273 inline bool Foam::sixDoFRigidBodyMotion::report() const
274 {
275  return report_;
276 }
278 inline bool Foam::sixDoFRigidBodyMotion::updateConstraints() const
279 {
280  return updateConstraints_;
281 }
282 
285 {
286  motionState0_ = motionState_;
287 }
288 
289 
291 {
292  return motionState_.centreOfRotation();
293 }
294 
295 
297 (
298  const point& pt
299 ) const
300 {
301  return (omega() ^ (pt - centreOfRotation())) + v();
302 }
303 
304 
306 (
307  const point& initialPoint
308 ) const
309 {
310  return
311  (
312  centreOfRotation()
313  + (Q() & initialQ_.T() & (initialPoint - initialCentreOfRotation_))
314  );
315 }
316 
317 // ************************************************************************* //
bool report() const
Return the report Switch.
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
point centreOfMass() const
Return the current centre of mass.
const point & centreOfRotation() const
Return the current centre of rotation.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Tensor< scalar > tensor
Definition: symmTensor.H:57
point transform(const point &initialPoints) const
Transform the given initial state point by the current motion.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const sixDoFRigidBodyMotionState & state() const
Return the motion state.
vector momentArm() const
Return the current momentArm.
dimensionedScalar cos(const dimensionedScalar &ds)
const diagTensor & momentOfInertia() const
Return the inertia tensor.
vector omega() const
Return the angular velocity in the global frame.
Holds the motion state of sixDoF object. Wrapped up together to allow rapid scatter to other processo...
const tensor & orientation() const
Return the orientation tensor, Q.
constexpr scalar pi(M_PI)
Vector< scalar > vector
Definition: vector.H:57
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void newTime()
Store the motion state at the beginning of the time-step.
point velocity(const point &pt) const
Return the velocity of a position.
const Time & time() const
Return time.
#define R(A, B, C, D, E, F, K, M)
scalar mass() const
Return the mass.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Tensor of scalars, i.e. Tensor<scalar>.
const point & initialCentreOfMass() const
Return the initial centre of mass.
const vector & v() const
Return the current velocity.