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-2021 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"
37
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 false
79 ),
80 pointMesh::New(mesh),
81 dimensionedScalar(dimless, Zero)
82 )
83{}
84
85
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 false
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();
210 }
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 tmp<pointField> ttransformedPts(new pointField(mesh().points()));
228 pointField& transformedPts = ttransformedPts.ref();
229
230 UIndirectList<point>(transformedPts, pointIDs()) =
231 pointField(newPoints.ref(), 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
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()
388 ).constrainDisplacement(pointDisplacement_);
389}
390
391
393(
394 IOstreamOption streamOpt,
395 const bool valid
396) const
397{
398 // Force ASCII writing
399 streamOpt.format(IOstream::ASCII);
400
402 (
404 (
405 "rigidBodyMotionState",
406 mesh().time().timeName(),
407 "uniform",
408 mesh(),
411 false
412 )
413 );
414
415 model_.state().write(dict);
416 return dict.regIOobject::writeObject(streamOpt, valid);
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.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
The IOstreamOption is a simple container for options an IOstream can normally have.
streamFormat format() const noexcept
Get the current stream format.
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:113
virtual bool read()
Re-read model coefficients if they have changed.
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:43
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
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
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.
Virtual base class for displacement motion solver.
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
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
A class for handling keywords in dictionaries.
Definition: keyType.H:71
Virtual base class for mesh motion solver.
Definition: motionSolver.H:61
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:144
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:150
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
const Type & lookupObject(const word &name, const bool recursive=false) const
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
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:81
label nPoints() const noexcept
Number of mesh points.
Rigid-body mesh motion solver for fvMesh.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write state using stream options.
virtual void solve()
Solve for motion.
virtual bool read()
Read dynamicMeshDict dictionary.
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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:47
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:50
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
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
dimensionedScalar cos(const dimensionedScalar &ds)
label timeIndex
Definition: getTimeIndex.H:30
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.
A non-counting (dummy) refCount.
Definition: refCount.H:59
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))