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-2022 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 "pointPatchDist.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].bodyID_));
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  bodyID_(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", 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->bodyID(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 pointPatchDist pDist(pMesh, bodyMeshes_[bi].patchSet_, points0());
227 
228  bodyMeshes_[bi].weight_.primitiveFieldRef() =
229  bodyMeshes_[bi].weight(pDist.primitiveField());
230 
231  pointConstraints::New(pMesh).constrain(bodyMeshes_[bi].weight_);
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
237 
239 {}
240 
241 
242 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
243 
244 template<class Type>
245 Type Foam::rigidBodyMeshMotion::bodyMesh::weight
246 (
247  const Type& pDist
248 ) const
249 {
250  // Scaling: 1 up to di then linear down to 0 at do away from patches
251  Type weight
252  (
253  min(max((do_ - pDist)/(do_ - di_), scalar(0)), scalar(1))
254  );
255 
256  // Convert the weight function to a cosine
257  weight =
258  min
259  (
260  max
261  (
262  0.5 - 0.5*cos(weight*Foam::constant::mathematical::pi),
263  scalar(0)
264  ),
265  scalar(1)
266  );
267 
268  return weight;
269 }
270 
271 
274 {
275  return points0() + pointDisplacement_.primitiveField();
276 }
277 
278 
280 {
281  const Time& t = mesh().time();
282 
283  if (mesh().nPoints() != points0().size())
284  {
286  << "The number of points in the mesh seems to have changed." << endl
287  << "In constant/polyMesh there are " << points0().size()
288  << " points; in the current mesh there are " << mesh().nPoints()
289  << " points." << exit(FatalError);
290  }
291 
292  // Store the motion state at the beginning of the time-step
293  if (curTimeIndex_ != t.timeIndex())
294  {
295  newTime();
296  curTimeIndex_ = t.timeIndex();
297  }
298 
299  const scalar ramp = ramp_->value(t.value());
300 
301  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
302  {
303  g() =
304  ramp
305  *mesh().lookupObject<uniformDimensionedVectorField>("g").value();
306  }
307 
308  if (test_)
309  {
310  label nIter(coeffDict().lookup<label>("nIter"));
311 
312  for (label i=0; i<nIter; i++)
313  {
315  (
316  t.value(),
317  t.deltaTValue(),
318  scalarField(nDoF(), Zero),
319  Field<spatialVector>(nBodies(), Zero)
320  );
321  }
322  }
323  else
324  {
325  Field<spatialVector> fx(nBodies(), Zero);
326 
327  forAll(bodyMeshes_, bi)
328  {
329  const label bodyID = bodyMeshes_[bi].bodyID_;
330 
331  dictionary forcesDict;
332  forcesDict.add("type", functionObjects::forces::typeName);
333  forcesDict.add("patches", bodyMeshes_[bi].patches_);
334  forcesDict.add("rhoInf", rhoInf_);
335  forcesDict.add("rho", rhoName_);
336  forcesDict.add("CofR", vector::zero);
337 
338  functionObjects::forces f("forces", t, forcesDict);
339  f.calcForcesMoment();
340 
341  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
342  }
343 
345  (
346  t.value(),
347  t.deltaTValue(),
348  scalarField(nDoF(), Zero),
349  fx
350  );
351  }
352 
353  if (Pstream::master() && report())
354  {
355  forAll(bodyMeshes_, bi)
356  {
357  status(bodyMeshes_[bi].bodyID_);
358  }
359  }
360 
361  vectorField& pointDisplacement = pointDisplacement_.primitiveFieldRef();
362  const pointField& points0 = this->points0();
363 
364  // Update the displacements
365  if (bodyMeshes_.size() == 1)
366  {
367  const septernion transform0(this->transform0(bodyMeshes_[0].bodyID_));
368  const scalarField& weight = bodyMeshes_[0].weight_;
369 
370  forAll(points0, pointi)
371  {
372  // Move non-stationary points
373  if (weight[pointi] > small)
374  {
375  // Use solid-body motion where weight = 1
376  if (weight[pointi] > 1 - small)
377  {
378  pointDisplacement[pointi] =
379  transform0.transformPoint(points0[pointi])
380  - points0[pointi];
381  }
382  // Slerp septernion interpolation
383  else
384  {
385  pointDisplacement[pointi] =
386  slerp(septernion::I, transform0, weight[pointi])
387  .transformPoint(points0[pointi])
388  - points0[pointi];
389  }
390  }
391  }
392  }
393  else
394  {
395  const List<septernion> transforms0(this->transforms0());
396  List<scalar> w(transforms0.size());
397 
398  forAll(points0, pointi)
399  {
400  pointDisplacement[pointi] =
401  average(transforms0, weights(pointi, w))
402  .transformPoint(points0[pointi])
403  - points0[pointi];
404  }
405  }
406 
407  // Displacement has changed. Update boundary conditions
409  (
410  pointDisplacement_.mesh()
411  ).constrainDisplacement(pointDisplacement_);
412 }
413 
414 
416 {
417  // pointMesh already updates pointFields
418 
419  // Get the new points either from the map or the mesh
420  const pointField& points =
421  (
422  map.hasMotionPoints()
423  ? map.preMotionPoints()
424  : mesh().points()
425  );
426 
427  const pointMesh& pMesh = pointMesh::New(mesh());
428  pointField points0(mesh().points());
429 
430  // Iterate to update the transformation of the new points to the
431  // corresponding points0, required because the body-point weights are
432  // calculated for points0
433  for (int iter=0; iter<3; iter++)
434  {
435  // Calculate scaling factor everywhere for each meshed body
436  forAll(bodyMeshes_, bi)
437  {
438  const pointPatchDist pDist
439  (
440  pMesh,
441  bodyMeshes_[bi].patchSet_,
442  points0
443  );
444 
445  pointScalarField& weight = bodyMeshes_[bi].weight_;
446 
447  forAll(points0, pointi)
448  {
449  const label oldPointi = map.pointMap()[pointi];
450 
451  if (oldPointi >= 0)
452  {
453  if (map.reversePointMap()[oldPointi] != pointi)
454  {
455  weight[pointi] = bodyMeshes_[bi].weight(pDist[pointi]);
456  }
457  }
458  else
459  {
461  << "Cannot determine co-ordinates "
462  "of introduced vertices."
463  << " New vertex " << pointi << " at co-ordinate "
464  << points[pointi] << exit(FatalError);
465  }
466  }
467 
468  pointConstraints::New(pMesh).constrain(weight);
469  }
470 
471  forAll(points0, pointi)
472  {
473  const label oldPointi = map.pointMap()[pointi];
474 
475  if (oldPointi >= 0)
476  {
477  if (map.reversePointMap()[oldPointi] == pointi)
478  {
479  // points0[pointi] = points0_[oldPointi];
480  points0[pointi] = points0_[pointi];
481  }
482  else
483  {
484  if (bodyMeshes_.size() == 1)
485  {
486  // Use solid-body motion where weight = 1
487  if (bodyMeshes_[0].weight_[pointi] > 1 - small)
488  {
489  points0[pointi] =
490  transform0(bodyMeshes_[0].bodyID_).inv()
491  .transformPoint(points[pointi]);
492  }
493  // Slerp septernion interpolation
494  else
495  {
496  points0[pointi] =
497  slerp
498  (
500  septernion
501  (
502  transform0(bodyMeshes_[0].bodyID_)
503  ),
504  bodyMeshes_[0].weight_[pointi]
505  ).invTransformPoint(points[pointi]);
506  }
507  }
508  else
509  {
510  const List<septernion> transforms0(this->transforms0());
511  List<scalar> w(transforms0.size());
512 
513  forAll(points0, pointi)
514  {
515  points0[pointi] =
516  average(transforms0, weights(pointi, w))
517  .invTransformPoint(points0[pointi]);
518  }
519  }
520  }
521  }
522  }
523  }
524 
525  points0_.transfer(points0);
526 
527  // points0 changed - set to write and check-in to database
528  points0_.rename("points0");
529  points0_.writeOpt() = IOobject::AUTO_WRITE;
530  points0_.instance() = mesh().time().name();
531  points0_.checkIn();
532 }
533 
534 
536 {
538  (
539  IOobject
540  (
541  "rigidBodyMotionState",
542  mesh().time().name(),
543  "uniform",
544  mesh(),
547  false
548  )
549  );
550 
551  state().write(dict);
552 
553  return
554  dict.regIOobject::writeObject
555  (
558  mesh().time().writeCompression(),
559  true
560  )
562 }
563 
564 
565 // ************************************************************************* //
#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 autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
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 internal 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 bodyID(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:113
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1169
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
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:208
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),.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
Calculation of distance to nearest patch for all points.
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 pointField & preMotionPoints() const
Pre-motion point positions.
const labelList & reversePointMap() const
Reverse point map.
const labelList & pointMap() const
Old point map.
bool hasMotionPoints() const
Has valid preMotionPoints?
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 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:306
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar cos(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
labelList f(nPoints)
dictionary dict