Loading...
Searching...
No Matches
rigidBodyMeshMotion.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) 2016-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2022 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
29#include "rigidBodyMeshMotion.H"
31#include "polyMesh.H"
32#include "pointPatchDist.H"
33#include "pointConstraints.H"
35#include "forces.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
43
45 (
49 );
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55Foam::rigidBodyMeshMotion::bodyMesh::bodyMesh
56(
57 const polyMesh& mesh,
58 const word& name,
59 const label bodyID,
60 const dictionary& dict
61)
62:
63 name_(name),
64 bodyID_(bodyID),
65 patches_(dict.get<wordRes>("patches")),
66 patchSet_(mesh.boundaryMesh().patchSet(patches_)),
67 di_(dict.get<scalar>("innerDistance")),
68 do_(dict.get<scalar>("outerDistance")),
69 weight_
70 (
71 IOobject
72 (
73 name_ + ".motionScale",
74 mesh.time().timeName(),
75 mesh,
76 IOobject::NO_READ,
77 IOobject::NO_WRITE,
78 IOobject::NO_REGISTER
79 ),
82 )
83{}
84
85
86Foam::rigidBodyMeshMotion::rigidBodyMeshMotion
87(
88 const polyMesh& mesh,
89 const IOdictionary& dict
90)
91:
93 model_
94 (
95 mesh.time(),
96 coeffDict(),
98 (
99 "rigidBodyMotionState",
100 mesh.time().timeName(),
101 "uniform",
102 mesh
103 ).typeHeaderOk<IOdictionary>(true)
105 (
107 (
108 "rigidBodyMotionState",
109 mesh.time().timeName(),
110 "uniform",
111 mesh,
112 IOobject::READ_IF_PRESENT,
113 IOobject::NO_WRITE,
114 IOobject::NO_REGISTER
115 )
116 )
117 : coeffDict()
118 ),
119 test_(coeffDict().getOrDefault("test", false)),
120 rhoInf_(1.0),
121 rhoName_(coeffDict().getOrDefault<word>("rho", "rho")),
122 ramp_
123 (
124 Function1<scalar>::NewIfPresent("ramp", coeffDict(), word::null, &mesh)
125 ),
126 curTimeIndex_(-1),
127 cOfGdisplacement_
128 (
129 coeffDict().getOrDefault<word>("cOfGdisplacement", "none")
130 ),
131 bodyIdCofG_(coeffDict().getOrDefault<label>("bodyIdCofG", -1))
132{
133 if (rhoName_ == "rhoInf")
134 {
135 readEntry("rhoInf", rhoInf_);
136 }
137
138 const dictionary& bodiesDict = coeffDict().subDict("bodies");
139
140 for (const entry& dEntry : bodiesDict)
141 {
142 const keyType& bodyName = dEntry.keyword();
143 const dictionary& bodyDict = dEntry.dict();
144
145 if (bodyDict.found("patches"))
146 {
147 const label bodyID = model_.bodyID(bodyName);
148
149 if (bodyID == -1)
150 {
152 << "Body " << bodyName
153 << " has been merged with another body"
154 " and cannot be assigned a set of patches"
155 << exit(FatalError);
156 }
157
158 bodyMeshes_.append
159 (
160 new bodyMesh
161 (
162 mesh,
163 bodyName,
164 bodyID,
165 bodyDict
166 )
167 );
168 }
169 }
170
171 // Calculate scaling factor everywhere for each meshed body
172 forAll(bodyMeshes_, bi)
173 {
174 const pointMesh& pMesh = pointMesh::New(mesh);
175
176 pointPatchDist pDist(pMesh, bodyMeshes_[bi].patchSet_, points0());
177
178 pointScalarField& scale = bodyMeshes_[bi].weight_;
179
180 // Scaling: 1 up to di then linear down to 0 at do away from patches
181 scale.primitiveFieldRef() =
182 min
183 (
184 max
185 (
186 (bodyMeshes_[bi].do_ - pDist.primitiveField())
187 /(bodyMeshes_[bi].do_ - bodyMeshes_[bi].di_),
188 scalar(0)
189 ),
190 scalar(1)
191 );
192
193 // Convert the scale function to a cosine
194 scale.primitiveFieldRef() =
195 min
196 (
197 max
198 (
199 0.5
200 - 0.5
201 *cos(scale.primitiveField()
203 scalar(0)
204 ),
205 scalar(1)
206 );
207
208 pointConstraints::New(pMesh).constrain(scale);
209 //scale.write();
211}
212
213
214// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215
218{
219 tmp<pointField> newPoints(points0() + pointDisplacement_.primitiveField());
220
221 if (moveAllCells())
222 {
223 return newPoints;
224 }
225 else
226 {
227 auto ttransformedPts = tmp<pointField>::New(mesh().points());
228 auto& transformedPts = ttransformedPts.ref();
229
230 UIndirectList<point>(transformedPts, pointIDs()) =
232
233 return ttransformedPts;
234 }
235}
236
237
239{
240 const Time& t = mesh().time();
241
242 if (mesh().nPoints() != points0().size())
243 {
245 << "The number of points in the mesh seems to have changed." << endl
246 << "In constant/polyMesh there are " << points0().size()
247 << " points; in the current mesh there are " << mesh().nPoints()
248 << " points." << exit(FatalError);
249 }
250
251 // Store the motion state at the beginning of the time-step
252 if (curTimeIndex_ != this->db().time().timeIndex())
253 {
254 model_.newTime();
255 curTimeIndex_ = this->db().time().timeIndex();
256 }
257
258 const scalar ramp = (ramp_ ? ramp_->value(t.value()) : 1.0);
259
260 if (t.foundObject<uniformDimensionedVectorField>("g"))
261 {
262 model_.g() =
263 ramp*t.lookupObject<uniformDimensionedVectorField>("g").value();
264 }
265
266 vector oldPos(vector::uniform(GREAT));
267 if (bodyIdCofG_ != -1)
268 {
269 oldPos = model_.cCofR(bodyIdCofG_);
270 }
271
272 if (test_)
273 {
274 const label nIter(coeffDict().get<label>("nIter"));
275
276 for (label i=0; i<nIter; i++)
277 {
278 model_.solve
279 (
280 t.value(),
281 t.deltaTValue(),
282 scalarField(model_.nDoF(), Zero),
283 Field<spatialVector>(model_.nBodies(), Zero)
284 );
285 }
286 }
287 else
288 {
289 const label nIter(coeffDict().getOrDefault<label>("nIter", 1));
290
291 for (label i=0; i<nIter; i++)
292 {
293 Field<spatialVector> fx(model_.nBodies(), Zero);
294
295 forAll(bodyMeshes_, bi)
296 {
297 const label bodyID = bodyMeshes_[bi].bodyID_;
298
299 dictionary forcesDict;
300 forcesDict.add("type", functionObjects::forces::typeName);
301 forcesDict.add("patches", bodyMeshes_[bi].patches_);
302 forcesDict.add("rhoInf", rhoInf_);
303 forcesDict.add("rho", rhoName_);
304 forcesDict.add("CofR", vector::zero);
305
306 functionObjects::forces f("forces", db(), forcesDict);
307 f.calcForcesMoments();
308
309 fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
310 }
311
312 model_.solve
313 (
314 t.value(),
315 t.deltaTValue(),
316 scalarField(model_.nDoF(), Zero),
317 fx
318 );
319 }
320
321 if (cOfGdisplacement_ != "none")
322 {
323 if (bodyIdCofG_ != -1)
324 {
325 if
326 (
327 db().time().foundObject<uniformDimensionedVectorField>
328 (
329 cOfGdisplacement_
330 )
331 )
332 {
333 auto& disp =
334 db().time().lookupObjectRef<uniformDimensionedVectorField>
335 (
336 cOfGdisplacement_
337 );
338
339 disp.value() += model_.cCofR(bodyIdCofG_) - oldPos;
340 }
341 }
342 else
343 {
345 << "CofGdisplacement is different to none." << endl
346 << "The model needs the entry body reference Id: bodyIdCofG."
347 << exit(FatalError);
348 }
349 }
350 }
351
352 if (Pstream::master() && model_.report())
353 {
354 forAll(bodyMeshes_, bi)
355 {
356 model_.status(bodyMeshes_[bi].bodyID_);
357 }
358 }
359
360 // Update the displacements
361 if (bodyMeshes_.size() == 1)
362 {
363 pointDisplacement_.primitiveFieldRef() = model_.transformPoints
364 (
365 bodyMeshes_[0].bodyID_,
366 bodyMeshes_[0].weight_,
367 points0()
368 ) - points0();
369 }
370 else
371 {
372 labelList bodyIDs(bodyMeshes_.size());
373 List<const scalarField*> weights(bodyMeshes_.size());
374 forAll(bodyIDs, bi)
375 {
376 bodyIDs[bi] = bodyMeshes_[bi].bodyID_;
377 weights[bi] = &bodyMeshes_[bi].weight_;
378 }
379
380 pointDisplacement_.primitiveFieldRef() =
381 model_.transformPoints(bodyIDs, weights, points0()) - points0();
382
383 }
384 // Displacement has changed. Update boundary conditions
386 (
387 pointDisplacement_.mesh()
389}
390
391
393(
394 IOstreamOption streamOpt,
395 const bool writeOnProc
396) const
397{
398 // Force ASCII writing
399 streamOpt.format(IOstreamOption::ASCII);
400
402 (
404 (
405 "rigidBodyMotionState",
406 mesh().time().timeName(),
407 "uniform",
408 mesh(),
412 )
413 );
414
415 model_.state().write(dict);
416 return dict.regIOobject::writeObject(streamOpt, writeOnProc);
417}
418
419
421{
423 {
424 model_.read(coeffDict());
425
426 return true;
427 }
428
429 return false;
430}
431
432
433// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
A simple container for options an IOstream can normally have.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
label bodyID(const word &name) const
Return the ID of the body with the given name.
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static const Form zero
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
friend class entry
Declare friendship with the entry class for IO.
Definition dictionary.H:338
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
const Type & value() const noexcept
Return const reference to value.
Virtual base class for displacement motion solver.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
pointVectorField pointDisplacement_
Point motion field.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition forces.H:321
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
A class for handling keywords in dictionaries.
Definition keyType.H:69
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
virtual tmp< pointField > newPoints()
Provide new points for motion. Solves for motion.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
virtual bool read()
Read dynamicMeshDict dictionary.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
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 constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Calculation of distance to nearest patch for all points.
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
label nPoints() const noexcept
Number of mesh points.
Rigid-body mesh motion solver for fvMesh.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write state using stream options.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void solve()
Solve for motion.
virtual bool read()
Read dynamicMeshDict dictionary.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
const labelList & pointIDs() const noexcept
The point ids (for cell set/zone subset).
Definition zoneMotion.H:102
bool moveAllCells() const noexcept
Move all cells?
Definition zoneMotion.H:110
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const pointField & points
label nPoints
word timeName
Definition getTimeIndex.H:3
constexpr scalar pi(M_PI)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
UniformDimensionedField< vector > uniformDimensionedVectorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
label timeIndex
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Various UniformDimensionedField types.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))