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 -------------------------------------------------------------------------------
11 License
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"
36 #include "mathematicalConstants.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(rigidBodyMeshMotion, 0);
43 
45  (
46  motionSolver,
47  rigidBodyMeshMotion,
48  dictionary
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 Foam::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 
86 Foam::rigidBodyMeshMotion::rigidBodyMeshMotion
87 (
88  const polyMesh& mesh,
89  const IOdictionary& dict
90 )
91 :
93  model_
94  (
95  mesh.time(),
96  coeffDict(),
97  IOobject
98  (
99  "rigidBodyMotionState",
100  mesh.time().timeName(),
101  "uniform",
102  mesh
103  ).typeHeaderOk<IOdictionary>(true)
104  ? IOdictionary
105  (
106  IOobject
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();
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 
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()
388  ).constrainDisplacement(pointDisplacement_);
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  (
403  IOobject
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 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Virtual base class for displacement motion solver.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
UniformDimensionedField< vector > uniformDimensionedVectorField
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write state using stream options.
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
"ascii" (normal default)
Various UniformDimensionedField types.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
A simple container for options an IOstream can normally have.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:173
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Macros for easy insertion into run-time selection tables.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual bool read()
Read dynamicMeshDict dictionary.
dictionary()
Default construct, a top-level empty dictionary.
Definition: dictionary.C:68
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...
word timeName
Definition: getTimeIndex.H:3
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
virtual void solve()
Solve for motion.
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:46
constexpr scalar pi(M_PI)
Vector< scalar > vector
Definition: vector.H:57
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual bool read()
Read dynamicMeshDict dictionary.
Definition: motionSolver.C:220
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
Definition: PtrList.H:344
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
Nothing to be read.
pointField & points0() noexcept
Return reference to the reference (&#39;0&#39;) pointField.
label bodyID(const word &name) const
Return the ID of the body with the given name.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
friend class entry
Declare friendship with the entry class for IO.
Definition: dictionary.H:325
List< label > labelList
A List of labels.
Definition: List.H:62
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:165
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
streamFormat format() const noexcept
Get the current stream format.
Do not request registration (bool: false)
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127