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-2024 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  di_(dict.lookup<scalar>("innerDistance")),
119  do_(dict.lookup<scalar>("outerDistance")),
120  weight_
121  (
122  IOobject
123  (
124  name_ + ".motionScale",
125  mesh.time().name(),
126  mesh,
127  IOobject::NO_READ,
128  IOobject::NO_WRITE
129  ),
130  pointMesh::New(mesh),
132  )
133 {}
134 
135 
137 (
138  const word& name,
139  const polyMesh& mesh,
140  const dictionary& dict
141 )
142 :
143  displacementMotionSolver(name, mesh, dict, typeName),
144  RBD::rigidBodyMotion
145  (
146  coeffDict(),
148  (
149  "rigidBodyMotionState",
150  mesh.time().name(),
151  "uniform",
152  mesh
153  ).headerOk()
155  (
156  IOobject
157  (
158  "rigidBodyMotionState",
159  mesh.time().name(),
160  "uniform",
161  mesh,
162  IOobject::READ_IF_PRESENT,
163  IOobject::NO_WRITE,
164  false
165  )
166  )
167  : coeffDict()
168  ),
169  test_(coeffDict().lookupOrDefault<Switch>("test", false)),
170  rhoInf_(1.0),
171  rhoName_(coeffDict().lookupOrDefault<word>("rho", "rho")),
172  ramp_(nullptr),
173  curTimeIndex_(-1)
174 {
175  if (rhoName_ == "rhoInf")
176  {
177  rhoInf_ = coeffDict().lookup<scalar>("rhoInf");
178  }
179 
180  if (coeffDict().found("ramp"))
181  {
182  ramp_ = Function1<scalar>::New("ramp", dimTime, dimless, coeffDict());
183  }
184  else
185  {
186  ramp_ = new Function1s::OneConstant<scalar>("ramp");
187  }
188 
189  const dictionary& bodiesDict = coeffDict().subDict("bodies");
190 
191  forAllConstIter(IDLList<entry>, bodiesDict, iter)
192  {
193  const dictionary& bodyDict = iter().dict();
194 
195  if (bodyDict.found("patches"))
196  {
197  const label bodyID = this->bodyIndex(iter().keyword());
198 
199  if (bodyID == -1)
200  {
202  << "Body " << iter().keyword()
203  << " has been merged with another body"
204  " and cannot be assigned a set of patches"
205  << exit(FatalError);
206  }
207 
208  bodyMeshes_.append
209  (
210  new bodyMesh
211  (
212  mesh,
213  iter().keyword(),
214  bodyID,
215  bodyDict
216  )
217  );
218  }
219  }
220 
221  const pointMesh& pMesh = pointMesh::New(mesh);
222 
223  // Calculate scaling factor everywhere for each meshed body
224  forAll(bodyMeshes_, bi)
225  {
226  const pointDist pDist
227  (
228  pMesh,
229  bodyMeshes_[bi].patchSet_,
230  points0(),
231  bodyMeshes_[bi].do_
232  );
233 
234  bodyMeshes_[bi].weight_.primitiveFieldRef() =
235  bodyMeshes_[bi].weight(pDist.primitiveField());
236 
237  pointConstraints::New(pMesh).constrain(bodyMeshes_[bi].weight_);
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
243 
245 {}
246 
247 
248 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
249 
250 template<class Type>
251 Type Foam::rigidBodyMeshMotion::bodyMesh::weight
252 (
253  const Type& pDist
254 ) const
255 {
256  // Scaling: 1 up to di then linear down to 0 at do away from patches
257  Type weight
258  (
259  min(max((do_ - pDist)/(do_ - di_), scalar(0)), scalar(1))
260  );
261 
262  // Convert the weight function to a cosine
263  weight =
264  min
265  (
266  max
267  (
268  0.5 - 0.5*cos(weight*Foam::constant::mathematical::pi),
269  scalar(0)
270  ),
271  scalar(1)
272  );
273 
274  return weight;
275 }
276 
277 
280 {
281  return points0() + pointDisplacement_.primitiveField();
282 }
283 
284 
286 {
287  const Time& t = mesh().time();
288 
289  if (mesh().nPoints() != points0().size())
290  {
292  << "The number of points in the mesh seems to have changed." << endl
293  << "In constant/polyMesh there are " << points0().size()
294  << " points; in the current mesh there are " << mesh().nPoints()
295  << " points." << exit(FatalError);
296  }
297 
298  // Store the motion state at the beginning of the time-step
299  if (curTimeIndex_ != t.timeIndex())
300  {
301  newTime();
302  curTimeIndex_ = t.timeIndex();
303  }
304 
305  const scalar ramp = ramp_->value(t.value());
306 
307  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
308  {
309  g() =
310  ramp
311  *mesh().lookupObject<uniformDimensionedVectorField>("g").value();
312  }
313 
314  if (test_)
315  {
316  label nIter(coeffDict().lookup<label>("nIter"));
317 
318  for (label i=0; i<nIter; i++)
319  {
321  (
322  t.value(),
323  t.deltaTValue(),
324  scalarField(nDoF(), Zero),
325  Field<spatialVector>(nBodies(), Zero)
326  );
327  }
328  }
329  else
330  {
331  Field<spatialVector> fx(nBodies(), Zero);
332 
333  forAll(bodyMeshes_, bi)
334  {
335  const label bodyID = bodyMeshes_[bi].bodyIndex_;
336 
338  (
339  functionObjects::forces::typeName,
340  t,
341  dictionary
342  (
343  "type", functionObjects::forces::typeName,
344  "patches", bodyMeshes_[bi].patches_,
345  "rhoInf", rhoInf_,
346  "rho", rhoName_,
347  "CofR", vector::zero
348  )
349  );
350 
351  f.calcForcesMoments();
352 
353  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
354  }
355 
357  (
358  t.value(),
359  t.deltaTValue(),
360  scalarField(nDoF(), Zero),
361  fx
362  );
363  }
364 
365  if (Pstream::master() && report())
366  {
367  forAll(bodyMeshes_, bi)
368  {
369  status(bodyMeshes_[bi].bodyIndex_);
370  }
371  }
372 
373  vectorField& pointDisplacement = pointDisplacement_.primitiveFieldRef();
374  const pointField& points0 = this->points0();
375 
376  // Update the displacements
377  if (bodyMeshes_.size() == 1)
378  {
379  const label bodyID = bodyMeshes_[0].bodyIndex_;
380  const septernion transform0(this->transform0(bodyID));
381  const scalarField& weight = bodyMeshes_[0].weight_;
382 
383  forAll(points0, pointi)
384  {
385  // Don't move where weight ~= 0
386  if (weight[pointi] <= small)
387  {}
388  // Use solid-body motion where weight ~= 1
389  else if (weight[pointi] > 1 - small)
390  {
391  pointDisplacement[pointi] =
392  transform0.transformPoint(points0[pointi])
393  - points0[pointi];
394  }
395  // Slerp septernion interpolation
396  else
397  {
398  pointDisplacement[pointi] =
399  slerp(septernion::I, transform0, weight[pointi])
400  .transformPoint(points0[pointi])
401  - points0[pointi];
402  }
403  }
404  }
405  else
406  {
407  const List<septernion> transforms0(this->transforms0());
408  List<scalar> w(transforms0.size());
409 
410  forAll(points0, pointi)
411  {
412  pointDisplacement[pointi] =
413  average(transforms0, weights(pointi, w))
414  .transformPoint(points0[pointi])
415  - points0[pointi];
416  }
417  }
418 
419  // Displacement has changed. Update boundary conditions
421  (
422  pointDisplacement_.mesh()
423  ).constrainDisplacement(pointDisplacement_);
424 }
425 
426 
428 {
429  // pointMesh already updates pointFields
430 
431  // Get the new points either from the map or the mesh
432  const pointField& points = mesh().points();
433 
434  const pointMesh& pMesh = pointMesh::New(mesh());
435 
436  pointField newPoints0(mesh().points());
437 
438  forAll(newPoints0, pointi)
439  {
440  if (map.pointMap()[pointi] < 0)
441  {
443  << "Cannot determine co-ordinates of introduced vertices."
444  << " New vertex " << pointi << " at co-ordinate "
445  << points[pointi] << exit(FatalError);
446  }
447  }
448 
449  // Iterate to update the transformation of the new points to the
450  // corresponding points0, required because the body-point weights are
451  // calculated for points0
452  for (int iter=0; iter<3; iter++)
453  {
454  // Calculate scaling factor everywhere for each meshed body
455  forAll(bodyMeshes_, bi)
456  {
457  const pointDist pDist
458  (
459  pMesh,
460  bodyMeshes_[bi].patchSet_,
461  newPoints0,
462  bodyMeshes_[bi].do_
463  );
464 
465  pointScalarField& weight = bodyMeshes_[bi].weight_;
466 
467  forAll(newPoints0, pointi)
468  {
469  const label oldPointi = map.pointMap()[pointi];
470 
471  if (map.reversePointMap()[oldPointi] != pointi)
472  {
473  weight[pointi] = bodyMeshes_[bi].weight(pDist[pointi]);
474  }
475  }
476 
477  pointConstraints::New(pMesh).constrain(weight);
478  }
479 
480  // Set directly mapped points
481  forAll(newPoints0, pointi)
482  {
483  const label oldPointi = map.pointMap()[pointi];
484 
485  if (map.reversePointMap()[oldPointi] == pointi)
486  {
487  newPoints0[pointi] = points0_[oldPointi];
488  }
489  }
490 
491  // Interpolate indirectly mapped points
492  if (bodyMeshes_.size() == 1)
493  {
494  const label bodyID = bodyMeshes_[0].bodyIndex_;
495  const septernion transform0(this->transform0(bodyID));
496  const scalarField& weight = bodyMeshes_[0].weight_;
497 
498  forAll(newPoints0, pointi)
499  {
500  const label oldPointi = map.pointMap()[pointi];
501 
502  if (map.reversePointMap()[oldPointi] != pointi)
503  {
504  // Don't move where weight ~= 0
505  if (weight[pointi] <= small)
506  {
507  newPoints0[pointi] = points[pointi];
508  }
509  // Use solid-body motion where weight ~= 1
510  else if (weight[pointi] > 1 - small)
511  {
512  newPoints0[pointi] =
513  transform0.invTransformPoint(points[pointi]);
514  }
515  // Slerp septernion interpolation
516  else
517  {
518  newPoints0[pointi] =
519  slerp(septernion::I, transform0, weight[pointi])
520  .invTransformPoint(points[pointi]);
521  }
522  }
523  }
524  }
525  else
526  {
527  const List<septernion> transforms0(this->transforms0());
528  List<scalar> w(transforms0.size());
529 
530  forAll(newPoints0, pointi)
531  {
532  const label oldPointi = map.pointMap()[pointi];
533 
534  if (map.reversePointMap()[oldPointi] != pointi)
535  {
536  newPoints0[pointi] =
537  average(transforms0, weights(pointi, w))
538  .invTransformPoint(newPoints0[pointi]);
539  }
540  }
541  }
542  }
543 
544  // Move into base class storage and mark as to-be-written
545  points0_.primitiveFieldRef() = newPoints0;
546  points0_.writeOpt() = IOobject::AUTO_WRITE;
547  points0_.instance() = mesh().time().name();
548 }
549 
550 
552 {
554  (
555  IOobject
556  (
557  "rigidBodyMotionState",
558  mesh().time().name(),
559  "uniform",
560  mesh(),
563  false
564  )
565  );
566 
567  state().write(dict);
568 
569  return
571  && dict.regIOobject::writeObject
572  (
575  mesh().time().writeCompression(),
576  true
577  );
578 }
579 
580 
581 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
const Type & value() const
Return const reference to value.
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
Virtual base class for mesh motion solver.
Definition: motionSolver.H:57
const word & keyword() const
Return keyword.
Definition: motionSolver.H:131
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:143
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:137
const Time & time() const
Return time.
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
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.
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:531
A class for handling words, derived from string.
Definition: word.H:62
#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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &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)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
error FatalError
quaternion slerp(const quaternion &qa, const quaternion &qb, const scalar t)
Spherical linear interpolation of quaternions.
Definition: quaternion.C:55
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList f(nPoints)
dictionary dict