sixDoFRigidBodyDisplacementPointPatchVectorField.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) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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#include "forces.H"
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40namespace Foam
41{
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45sixDoFRigidBodyDisplacementPointPatchVectorField::
46sixDoFRigidBodyDisplacementPointPatchVectorField
47(
48 const pointPatch& p,
50)
51:
53 motion_(db().time()),
54 initialPoints_(p.localPoints()),
55 rhoInf_(1.0),
56 rhoName_("rho"),
57 lookupGravity_(-1),
58 g_(Zero),
59 curTimeIndex_(-1)
60{}
61
62
63sixDoFRigidBodyDisplacementPointPatchVectorField::
64sixDoFRigidBodyDisplacementPointPatchVectorField
65(
66 const pointPatch& p,
68 const dictionary& dict
69)
70:
72 motion_(dict, dict, db().time()),
73 rhoInf_(1.0),
74 rhoName_(dict.getOrDefault<word>("rho", "rho")),
75 lookupGravity_(-1),
76 g_(Zero),
77 curTimeIndex_(-1)
78{
79 if (rhoName_ == "rhoInf")
80 {
81 dict.readEntry("rhoInf", rhoInf_);
82 }
83
84 if (dict.readIfPresent("g", g_))
85 {
86 lookupGravity_ = -2;
87 }
88
89 if (!dict.found("value"))
90 {
92 }
93
94 if (dict.found("initialPoints"))
95 {
96 initialPoints_ = vectorField("initialPoints", dict , p.size());
97 }
98 else
99 {
100 initialPoints_ = p.localPoints();
101 }
102}
103
104
105sixDoFRigidBodyDisplacementPointPatchVectorField::
106sixDoFRigidBodyDisplacementPointPatchVectorField
107(
109 const pointPatch& p,
111 const pointPatchFieldMapper& mapper
112)
113:
114 fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
115 motion_(ptf.motion_),
116 initialPoints_(ptf.initialPoints_, mapper),
117 rhoInf_(ptf.rhoInf_),
118 rhoName_(ptf.rhoName_),
119 lookupGravity_(ptf.lookupGravity_),
120 g_(ptf.g_),
121 curTimeIndex_(-1)
122{}
123
124
125sixDoFRigidBodyDisplacementPointPatchVectorField::
126sixDoFRigidBodyDisplacementPointPatchVectorField
127(
130)
131:
133 motion_(ptf.motion_),
134 initialPoints_(ptf.initialPoints_),
135 rhoInf_(ptf.rhoInf_),
136 rhoName_(ptf.rhoName_),
137 lookupGravity_(ptf.lookupGravity_),
138 g_(ptf.g_),
139 curTimeIndex_(-1)
140{}
141
142
143// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144
146(
147 const pointPatchFieldMapper& m
148)
149{
151
152 initialPoints_.autoMap(m);
153}
154
155
157(
158 const pointPatchField<vector>& ptf,
159 const labelList& addr
160)
161{
163 refCast<const sixDoFRigidBodyDisplacementPointPatchVectorField>(ptf);
164
166
167 initialPoints_.rmap(sDoFptf.initialPoints_, addr);
168}
169
170
172{
173 if (this->updated())
174 {
175 return;
176 }
177
178 if (lookupGravity_ < 0)
179 {
180 if (db().time().foundObject<uniformDimensionedVectorField>("g"))
181 {
182 if (lookupGravity_ == -2)
183 {
185 << "Specifying the value of g in this boundary condition "
186 << "when g is available from the database is considered "
187 << "a fatal error to avoid the possibility of inconsistency"
188 << exit(FatalError);
189 }
190 else
191 {
192 lookupGravity_ = 1;
193 }
194 }
195 else
196 {
197 lookupGravity_ = 0;
198 }
199 }
200
201 const polyMesh& mesh = this->internalField().mesh()();
202 const Time& t = mesh.time();
203 const pointPatch& ptPatch = this->patch();
204
205 // Store the motion state at the beginning of the time-step
206 bool firstIter = false;
207 if (curTimeIndex_ != t.timeIndex())
208 {
209 motion_.newTime();
210 curTimeIndex_ = t.timeIndex();
211 firstIter = true;
212 }
213
214 dictionary forcesDict;
215
216 forcesDict.add("type", functionObjects::forces::typeName);
217 forcesDict.add("patches", wordList(1, ptPatch.name()));
218 forcesDict.add("rhoInf", rhoInf_);
219 forcesDict.add("rho", rhoName_);
220 forcesDict.add("CofR", motion_.centreOfRotation());
221
222 functionObjects::forces f("forces", db(), forcesDict);
223
224 f.calcForcesMoments();
225
226 // Get the forces on the patch faces at the current positions
227
228 if (lookupGravity_ == 1)
229 {
232
233 g_ = g.value();
234 }
235
236 // scalar ramp = min(max((t.value() - 5)/10, 0), 1);
237 scalar ramp = 1.0;
238
239 motion_.update
240 (
241 firstIter,
242 ramp*(f.forceEff() + motion_.mass()*g_),
243 ramp*(f.momentEff() + motion_.mass()*(motion_.momentArm() ^ g_)),
244 t.deltaTValue(),
245 t.deltaT0Value()
246 );
247
249 (
250 motion_.transform(initialPoints_) - initialPoints_
251 );
252
254}
255
256
258{
260
261 os.writeEntry("rho", rhoName_);
262
263 if (rhoName_ == "rhoInf")
264 {
265 os.writeEntry("rhoInf", rhoInf_);
266 }
267
268 if (lookupGravity_ == 0 || lookupGravity_ == -2)
269 {
270 os.writeEntry("g", g_);
271 }
272
273 motion_.write(os);
274
275 initialPoints_.writeEntry("initialPoints", os);
276
277 writeEntry("value", os);
278}
279
280
281// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282
284(
287);
288
289// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290
291} // End namespace Foam
292
293// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
Generic templated field type.
Definition: Field.H:82
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:403
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:614
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:466
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
scalar deltaT0Value() const noexcept
Return old time step value.
Definition: TimeStateI.H:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const Type & value() const
Return const reference to value.
A FixedValue boundary condition for pointField.
virtual bool write()
Write the output fields.
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition: forces.H:325
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const Time & time() const noexcept
Return time registry.
const Type & lookupObject(const word &name, const bool recursive=false) const
Foam::pointPatchFieldMapper.
Abstract base class for point-mesh patch fields.
const pointPatch & patch() const
Return patch.
const objectRegistry & db() const
Return local objectRegistry.
const DimensionedField< Type, pointMesh > & internalField() const
Return dimensioned internal field reference.
bool updated() const
Return true if the boundary condition has already been updated.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
virtual const word & name() const =0
Return name.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual void autoMap(const pointPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const pointPatchField< vector > &, const labelList &)
Reverse map the given pointPatchField onto this pointPatchField.
void newTime()
Store the motion state at the beginning of the time-step.
void update(bool firstIter, const vector &fGlobal, const vector &tauGlobal, scalar deltaT, scalar deltaT0)
Symplectic integration of velocities, orientation and position.
void write(Ostream &) const
Write.
point transform(const point &initialPoints) const
Transform the given initial state point by the current motion.
vector momentArm() const
Return the current momentArm.
scalar mass() const
Return the mass.
const point & centreOfRotation() const
Return the current centre of rotation.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
labelList f(nPoints)
dictionary dict
static const char *const typeName
The type name used in ensight case files.