uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField.C
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-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2019 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 
30 #include "pointPatchFields.H"
32 #include "Time.H"
33 #include "fvMesh.H"
34 #include "volFields.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
46 (
47  const pointPatch& p,
49 )
50 :
52  motion_(db().time()),
53  initialPoints_(p.localPoints()),
54  curTimeIndex_(-1)
55 {}
56 
57 
60 (
61  const pointPatch& p,
63  const dictionary& dict
64 )
65 :
67  motion_(dict, dict, db().time()),
68  curTimeIndex_(-1)
69 {
70  if (!dict.found("value"))
71  {
72  updateCoeffs();
73  }
74 
75  if (dict.found("initialPoints"))
76  {
77  initialPoints_ = vectorField("initialPoints", dict , p.size());
78  }
79  else
80  {
81  initialPoints_ = p.localPoints();
82  }
83 }
84 
85 
88 (
89  const uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField& ptf,
90  const pointPatch& p,
91  const DimensionedField<vector, pointMesh>& iF,
92  const pointPatchFieldMapper& mapper
93 )
94 :
95  fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
96  motion_(ptf.motion_),
97  initialPoints_(ptf.initialPoints_, mapper),
98  curTimeIndex_(-1)
99 {}
100 
101 
104 (
107 )
108 :
110  motion_(ptf.motion_),
111  initialPoints_(ptf.initialPoints_),
112  curTimeIndex_(-1)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 (
120  const pointPatchFieldMapper& m
121 )
122 {
124 
125  initialPoints_.autoMap(m);
126 }
127 
128 
130 (
131  const pointPatchField<vector>& ptf,
132  const labelList& addr
133 )
134 {
136  refCast
137  <
139  >(ptf);
142 
143  initialPoints_.rmap(uSDoFptf.initialPoints_, addr);
144 }
145 
146 
148 {
149  if (this->updated())
150  {
151  return;
152  }
153 
154  const polyMesh& mesh = this->internalField().mesh()();
155  const Time& t = mesh.time();
156 
157  // Store the motion state at the beginning of the time-step
158  bool firstIter = false;
159  if (curTimeIndex_ != t.timeIndex())
160  {
161  motion_.newTime();
162  curTimeIndex_ = t.timeIndex();
163  firstIter = true;
164  }
165 
166  vector gravity = Zero;
167 
168  if (db().time().foundObject<uniformDimensionedVectorField>("g"))
169  {
172 
173  gravity = g.value();
174  }
175 
176  // Do not modify the accelerations, except with gravity, where available
177  motion_.update
178  (
179  firstIter,
180  gravity*motion_.mass(),
181  Zero,
182  t.deltaTValue(),
183  t.deltaT0Value()
184  );
185 
187  (
188  motion_.transform(initialPoints_) - initialPoints_
189  );
190 
192 }
193 
194 
196 (
197  Ostream& os
198 ) const
199 {
201  motion_.write(os);
202  initialPoints_.writeEntry("initialPoints", os);
203  this->writeValueEntry(os);
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
210 (
212  uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField
213 );
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 } // End namespace Foam
218 
219 // ************************************************************************* //
dictionary dict
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
void write(Ostream &) const
Write.
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition: typeInfo.H:172
Various UniformDimensionedField types.
Foam::pointPatchFieldMapper.
virtual void write(Ostream &os) const
Write.
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:720
point transform(const point &initialPoints) const
Transform the given initial state point by the current motion.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
friend Ostream & operator(Ostream &, const Field< vector > &)
dynamicFvMesh & mesh
const Time & time() const noexcept
Return time registry.
Vector< scalar > vector
Definition: vector.H:57
bool updated() const noexcept
True if the boundary condition has already been updated.
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const objectRegistry & db() const
The associated objectRegistry.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
void newTime()
Store the motion state at the beginning of the time-step.
pointPatchField< vector > pointPatchVectorField
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:466
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:529
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
scalar mass() const
Return the mass.
makePointPatchTypeField(pointPatchVectorField, solidBodyMotionDisplacementPointPatchVectorField)
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
uncoupledSixDoFRigidBodyDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void rmap(const pointPatchField< Type > &, const labelList &)
Reverse map the given PointPatchField onto.
void update(bool firstIter, const vector &fGlobal, const vector &tauGlobal, scalar deltaT, scalar deltaT0)
Symplectic integration of velocities, orientation and position.
const DimensionedField< vector, pointMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127