rigidBodyMeshMotion.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "rigidBodyMeshMotion.H"
27 #include "polyMesh.H"
28 #include "polyTopoChangeMap.H"
29 #include "pointDist.H"
30 #include "pointConstraints.H"
31 #include "timeIOdictionary.H"
33 #include "forces.H"
34 #include "OneConstant.H"
35 #include "mathematicalConstants.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
43 
45  (
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 Foam::List<Foam::septernion> Foam::rigidBodyMeshMotion::transforms0() const
56 {
57  List<septernion> transforms0(bodyMeshes_.size() + 1);
58 
59  forAll(bodyMeshes_, bi)
60  {
61  // Calculate the septernion equivalent of the transformation
62  transforms0[bi] = septernion(transform0(bodyMeshes_[bi].bodyIndex_));
63  }
64 
65  transforms0[bodyMeshes_.size()] = septernion::I;
66 
67  return transforms0;
68 }
69 
70 
71 Foam::List<Foam::scalar>& Foam::rigidBodyMeshMotion::weights
72 (
73  const label pointi,
74  List<scalar>& w
75 ) const
76 {
77  // Initialise to 1 for the far-field weight
78  scalar sum1mw = 1;
79 
80  forAll(bodyMeshes_, bi)
81  {
82  w[bi] = bodyMeshes_[bi].weight_[pointi];
83  sum1mw += w[bi]/(1 + small - w[bi]);
84  }
85 
86  // Calculate the limiter for wi/(1 - wi) to ensure the sum(wi) = 1
87  scalar lambda = 1/sum1mw;
88 
89  // Limit wi/(1 - wi) and sum the resulting wi
90  scalar sumw = 0;
91  forAll(bodyMeshes_, bi)
92  {
93  w[bi] = lambda*w[bi]/(1 + small - w[bi]);
94  sumw += w[bi];
95  }
96 
97  // Calculate the weight for the stationary far-field
98  w[bodyMeshes_.size()] = 1 - sumw;
99 
100  return w;
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
106 Foam::rigidBodyMeshMotion::bodyMesh::bodyMesh
107 (
108  const polyMesh& mesh,
109  const word& name,
110  const label bodyID,
111  const dictionary& dict
112 )
113 :
114  name_(name),
115  bodyIndex_(bodyID),
116  patches_(wordReList(dict.lookup("patches"))),
117  patchSet_(mesh.boundaryMesh().patchSet(patches_)),
118  pointZones_(dict.lookupOrDefault("pointZones", wordReList::null())),
119  pointZoneSet_(mesh.pointZones().zoneSet(pointZones_)),
120  di_(dict.lookup<scalar>("innerDistance")),
121  do_(dict.lookup<scalar>("outerDistance")),
122  weight_
123  (
124  IOobject
125  (
126  name_ + ".motionScale",
127  mesh.time().name(),
128  mesh,
129  IOobject::NO_READ,
130  IOobject::NO_WRITE
131  ),
132  pointMesh::New(mesh),
134  )
135 {}
136 
137 
139 (
140  const word& name,
141  const polyMesh& mesh,
142  const dictionary& dict
143 )
144 :
145  displacementMotionSolver(name, mesh, dict, typeName),
146  RBD::rigidBodyMotion
147  (
148  dict,
150  (
151  "rigidBodyMotionState",
152  mesh.time().name(),
153  "uniform",
154  mesh
155  ).headerOk()
157  (
158  IOobject
159  (
160  "rigidBodyMotionState",
161  mesh.time().name(),
162  "uniform",
163  mesh,
164  IOobject::READ_IF_PRESENT,
165  IOobject::NO_WRITE,
166  false
167  )
168  )
169  : dict
170  ),
171  test_(dict.lookupOrDefault<Switch>("test", false)),
172  nIter_(test_ ? dict.lookup<label>("nIter") : 0),
173  rhoInf_(1.0),
174  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
175  ramp_(nullptr),
176  curTimeIndex_(-1)
177 {
178  if (rhoName_ == "rhoInf")
179  {
180  rhoInf_ = dict.lookup<scalar>("rhoInf");
181  }
182 
183  if (dict.found("ramp"))
184  {
185  ramp_ = Function1<scalar>::New("ramp", dimTime, dimless, dict);
186  }
187  else
188  {
189  ramp_ = new Function1s::OneConstant<scalar>("ramp");
190  }
191 
192  const dictionary& bodiesDict = dict.subDict("bodies");
193 
194  forAllConstIter(IDLList<entry>, bodiesDict, iter)
195  {
196  const dictionary& bodyDict = iter().dict();
197 
198  if (bodyDict.found("patches"))
199  {
200  const label bodyID = this->bodyIndex(iter().keyword());
201 
202  if (bodyID == -1)
203  {
205  << "Body " << iter().keyword()
206  << " has been merged with another body"
207  " and cannot be assigned a set of patches"
208  << exit(FatalError);
209  }
210 
211  bodyMeshes_.append
212  (
213  new bodyMesh
214  (
215  mesh,
216  iter().keyword(),
217  bodyID,
218  bodyDict
219  )
220  );
221  }
222  }
223 
224  const pointMesh& pMesh = pointMesh::New(mesh);
225 
226  // Calculate scaling factor everywhere for each meshed body
227  forAll(bodyMeshes_, bi)
228  {
229  const pointDist pDist
230  (
231  pMesh,
232  bodyMeshes_[bi].patchSet_,
233  bodyMeshes_[bi].pointZoneSet_,
234  points0(),
235  bodyMeshes_[bi].do_
236  );
237 
238  bodyMeshes_[bi].weight_.primitiveFieldRef() =
239  bodyMeshes_[bi].weight(pDist.primitiveField());
240 
241  pointConstraints::New(pMesh).constrain(bodyMeshes_[bi].weight_);
242  }
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
247 
249 {}
250 
251 
252 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
253 
254 template<class Type>
255 Type Foam::rigidBodyMeshMotion::bodyMesh::weight
256 (
257  const Type& pDist
258 ) const
259 {
260  // Scaling: 1 up to di then linear down to 0 at do away from patches
261  Type weight
262  (
263  min(max((do_ - pDist)/(do_ - di_), scalar(0)), scalar(1))
264  );
265 
266  // Convert the weight function to a cosine
267  weight =
268  min
269  (
270  max
271  (
272  0.5 - 0.5*cos(weight*Foam::constant::mathematical::pi),
273  scalar(0)
274  ),
275  scalar(1)
276  );
277 
278  return weight;
279 }
280 
281 
284 {
285  return points0() + pointDisplacement_.primitiveField();
286 }
287 
288 
290 {
291  const Time& t = mesh().time();
292 
293  if (mesh().nPoints() != points0().size())
294  {
296  << "The number of points in the mesh seems to have changed." << endl
297  << "In constant/polyMesh there are " << points0().size()
298  << " points; in the current mesh there are " << mesh().nPoints()
299  << " points." << exit(FatalError);
300  }
301 
302  // Store the motion state at the beginning of the time-step
303  if (curTimeIndex_ != t.timeIndex())
304  {
305  newTime();
306  curTimeIndex_ = t.timeIndex();
307  }
308 
309  const scalar ramp = ramp_->value(t.value());
310 
311  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
312  {
313  g() =
314  ramp
316  }
317 
318  if (test_)
319  {
320  for (label i=0; i<nIter_; i++)
321  {
323  (
324  t.value(),
325  t.deltaTValue(),
326  scalarField(nDoF(), Zero),
327  Field<spatialVector>(nBodies(), Zero)
328  );
329  }
330  }
331  else
332  {
333  Field<spatialVector> fx(nBodies(), Zero);
334 
335  forAll(bodyMeshes_, bi)
336  {
337  const label bodyID = bodyMeshes_[bi].bodyIndex_;
338 
340  (
341  functionObjects::forces::typeName,
342  t,
344  (
345  "type", functionObjects::forces::typeName,
346  "patches", bodyMeshes_[bi].patches_,
347  "rhoInf", rhoInf_,
348  "rho", rhoName_,
349  "CofR", vector::zero
350  )
351  );
352 
353  f.calcForcesMoments();
354 
355  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
356  }
357 
359  (
360  t.value(),
361  t.deltaTValue(),
362  scalarField(nDoF(), Zero),
363  fx
364  );
365  }
366 
367  if (Pstream::master() && report())
368  {
369  forAll(bodyMeshes_, bi)
370  {
371  status(bodyMeshes_[bi].bodyIndex_);
372  }
373  }
374 
375  vectorField& pointDisplacement = pointDisplacement_.primitiveFieldRef();
376  const pointField& points0 = this->points0();
377 
378  // Update the displacements
379  if (bodyMeshes_.size() == 1)
380  {
381  const label bodyID = bodyMeshes_[0].bodyIndex_;
382  const septernion transform0(this->transform0(bodyID));
383  const scalarField& weight = bodyMeshes_[0].weight_;
384 
385  forAll(points0, pointi)
386  {
387  // Don't move where weight ~= 0
388  if (weight[pointi] <= small)
389  {}
390  // Use solid-body motion where weight ~= 1
391  else if (weight[pointi] > 1 - small)
392  {
393  pointDisplacement[pointi] =
394  transform0.transformPoint(points0[pointi])
395  - points0[pointi];
396  }
397  // Slerp septernion interpolation
398  else
399  {
400  pointDisplacement[pointi] =
401  slerp(septernion::I, transform0, weight[pointi])
402  .transformPoint(points0[pointi])
403  - points0[pointi];
404  }
405  }
406  }
407  else
408  {
409  const List<septernion> transforms0(this->transforms0());
410  List<scalar> w(transforms0.size());
411 
412  forAll(points0, pointi)
413  {
414  pointDisplacement[pointi] =
415  average(transforms0, weights(pointi, w))
416  .transformPoint(points0[pointi])
417  - points0[pointi];
418  }
419  }
420 
421  // Displacement has changed. Update boundary conditions
423  (
424  pointDisplacement_.mesh()
425  ).constrainDisplacement(pointDisplacement_);
426 }
427 
428 
430 {
431  // pointMesh already updates pointFields
432 
433  // Get the new points either from the map or the mesh
434  const pointField& points = mesh().points();
435 
436  const pointMesh& pMesh = pointMesh::New(mesh());
437 
438  pointField newPoints0(mesh().points());
439 
440  forAll(newPoints0, pointi)
441  {
442  if (map.pointMap()[pointi] < 0)
443  {
445  << "Cannot determine co-ordinates of introduced vertices."
446  << " New vertex " << pointi << " at co-ordinate "
447  << points[pointi] << exit(FatalError);
448  }
449  }
450 
451  // Iterate to update the transformation of the new points to the
452  // corresponding points0, required because the body-point weights are
453  // calculated for points0
454  for (int iter=0; iter<3; iter++)
455  {
456  // Calculate scaling factor everywhere for each meshed body
457  forAll(bodyMeshes_, bi)
458  {
459  const pointDist pDist
460  (
461  pMesh,
462  bodyMeshes_[bi].patchSet_,
463  bodyMeshes_[bi].pointZoneSet_,
464  newPoints0,
465  bodyMeshes_[bi].do_
466  );
467 
468  pointScalarField& weight = bodyMeshes_[bi].weight_;
469 
470  forAll(newPoints0, pointi)
471  {
472  const label oldPointi = map.pointMap()[pointi];
473 
474  if (map.reversePointMap()[oldPointi] != pointi)
475  {
476  weight[pointi] = bodyMeshes_[bi].weight(pDist[pointi]);
477  }
478  }
479 
480  pointConstraints::New(pMesh).constrain(weight);
481  }
482 
483  // Set directly mapped points
484  forAll(newPoints0, pointi)
485  {
486  const label oldPointi = map.pointMap()[pointi];
487 
488  if (map.reversePointMap()[oldPointi] == pointi)
489  {
490  newPoints0[pointi] = points0_[oldPointi];
491  }
492  }
493 
494  // Interpolate indirectly mapped points
495  if (bodyMeshes_.size() == 1)
496  {
497  const label bodyID = bodyMeshes_[0].bodyIndex_;
498  const septernion transform0(this->transform0(bodyID));
499  const scalarField& weight = bodyMeshes_[0].weight_;
500 
501  forAll(newPoints0, pointi)
502  {
503  const label oldPointi = map.pointMap()[pointi];
504 
505  if (map.reversePointMap()[oldPointi] != pointi)
506  {
507  // Don't move where weight ~= 0
508  if (weight[pointi] <= small)
509  {
510  newPoints0[pointi] = points[pointi];
511  }
512  // Use solid-body motion where weight ~= 1
513  else if (weight[pointi] > 1 - small)
514  {
515  newPoints0[pointi] =
516  transform0.invTransformPoint(points[pointi]);
517  }
518  // Slerp septernion interpolation
519  else
520  {
521  newPoints0[pointi] =
522  slerp(septernion::I, transform0, weight[pointi])
523  .invTransformPoint(points[pointi]);
524  }
525  }
526  }
527  }
528  else
529  {
530  const List<septernion> transforms0(this->transforms0());
531  List<scalar> w(transforms0.size());
532 
533  forAll(newPoints0, pointi)
534  {
535  const label oldPointi = map.pointMap()[pointi];
536 
537  if (map.reversePointMap()[oldPointi] != pointi)
538  {
539  newPoints0[pointi] =
540  average(transforms0, weights(pointi, w))
541  .invTransformPoint(newPoints0[pointi]);
542  }
543  }
544  }
545  }
546 
547  // Move into base class storage and mark as to-be-written
548  points0_.primitiveFieldRef() = newPoints0;
549  points0_.writeOpt() = IOobject::AUTO_WRITE;
550  points0_.instance() = mesh().time().name();
551 }
552 
553 
555 {
557  (
558  IOobject
559  (
560  "rigidBodyMotionState",
561  mesh().time().name(),
562  "uniform",
563  mesh(),
566  false
567  )
568  );
569 
570  state().write(dict);
571 
572  return
574  && dict.regIOobject::writeObject
575  (
578  mesh().time().writeCompression(),
579  true
580  );
581 }
582 
583 
584 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Templated function that returns the corresponding 1 (one).
Definition: OneConstant.H:59
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
Template class for intrusive linked lists.
Definition: ILList.H:67
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:203
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
label bodyIndex(const word &name) const
Return the ID of the body with the given name.
spatialTransform transform0(const label bodyID) const
Return the transformation of bodyID relative to the initial time.
void solve(const scalar t, const scalar deltaT, const scalarField &tau, const Field< spatialVector > &fx)
Integrate velocities, orientation and position.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const Form zero
Definition: VectorSpace.H:118
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
static std::tuple< const Entries &... > entries(const Entries &...)
Construct an entries tuple from which to make a dictionary.
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
Virtual base class for displacement motion solver.
Calculates the forces and moments by integrating the pressure and skin-friction forces over a given l...
Definition: forces.H:188
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
const word & keyword() const
Return keyword.
Definition: motionSolver.H:127
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:133
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
void constrain(PointField< Type > &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),.
Calculates the distance to the specified sets of patch and pointZone points or for all points.
Definition: pointDist.H:54
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
virtual bool write() const
Write points0 if the mesh topology changed.
pointField & points0()
Return reference to the reference field.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const
Reverse point map.
const labelList & pointMap() const
Old point map.
label nPoints() const
Rigid-body mesh motion solver for fvMesh.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
rigidBodyMeshMotion(const word &name, const polyMesh &, const dictionary &dict)
Construct from polyMesh and dictionary.
virtual bool write() const
Write motion state information for restart.
virtual void solve()
Solve for motion.
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:66
vector invTransformPoint(const vector &v) const
Inverse Transform the given coordinate point.
Definition: septernionI.H:88
vector transformPoint(const vector &v) const
Transform the given coordinate point.
Definition: septernionI.H:82
static const septernion I
Definition: septernion.H:83
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
label nPoints
dimensionedScalar lambda(viscosity->lookup("lambda"))
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
dimensioned< Type > average(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:47
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
quaternion slerp(const quaternion &qa, const quaternion &qb, const scalar t)
Spherical linear interpolation of quaternions.
Definition: quaternion.C:55
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList f(nPoints)
dictionary dict