rigidBodyMeshMotionSolver.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 
27 #include "polyMesh.H"
28 #include "pointPatchDist.H"
29 #include "pointConstraints.H"
30 #include "timeIOdictionary.H"
32 #include "forces.H"
33 #include "transformField.H"
34 #include "OneConstant.H"
35 #include "mathematicalConstants.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(rigidBodyMeshMotionSolver, 0);
43 
45  (
46  motionSolver,
47  rigidBodyMeshMotionSolver,
48  dictionary
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 Foam::rigidBodyMeshMotionSolver::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_(wordReList(dict.lookup("patches"))),
66  patchSet_(mesh.boundaryMesh().patchSet(patches_))
67 {}
68 
69 
71 (
72  const word& name,
73  const polyMesh& mesh,
74  const dictionary& dict
75 )
76 :
77  motionSolver(name, mesh, dict, typeName),
79  (
80  coeffDict(),
82  (
83  "rigidBodyMotionState",
84  mesh.time().timeName(),
85  "uniform",
86  mesh
87  ).headerOk()
89  (
90  IOobject
91  (
92  "rigidBodyMotionState",
93  mesh.time().timeName(),
94  "uniform",
95  mesh,
98  false
99  )
100  )
101  : coeffDict()
102  ),
103  test_(coeffDict().lookupOrDefault<Switch>("test", false)),
104  rhoInf_(1.0),
105  rhoName_(coeffDict().lookupOrDefault<word>("rho", "rho")),
106  ramp_(nullptr),
107  curTimeIndex_(-1),
108  meshSolverPtr_
109  (
111  (
112  name,
113  mesh,
115  (
116  IOobject
117  (
118  "rigidBodyMotionSolver:meshSolver",
119  mesh.time().constant(),
120  mesh
121  ),
122  coeffDict().subDict("meshSolver")
123  )
124  )
125  ),
126  meshSolver_(refCast<displacementMotionSolver>(meshSolverPtr_()))
127 {
128  if (rhoName_ == "rhoInf")
129  {
130  rhoInf_ = coeffDict().lookup<scalar>("rhoInf");
131  }
132 
133  if (coeffDict().found("ramp"))
134  {
135  ramp_ = Function1<scalar>::New("ramp", coeffDict());
136  }
137  else
138  {
139  ramp_ = new Function1s::OneConstant<scalar>("ramp");
140  }
141 
142  const dictionary& bodiesDict = coeffDict().subDict("bodies");
143 
144  forAllConstIter(IDLList<entry>, bodiesDict, iter)
145  {
146  const dictionary& bodyDict = iter().dict();
147 
148  if (bodyDict.found("patches"))
149  {
150  const label bodyID = this->bodyID(iter().keyword());
151 
152  if (bodyID == -1)
153  {
155  << "Body " << iter().keyword()
156  << " has been merged with another body"
157  " and cannot be assigned a set of patches"
158  << exit(FatalError);
159  }
160 
161  bodyMeshes_.append
162  (
163  new bodyMesh
164  (
165  mesh,
166  iter().keyword(),
167  bodyID,
168  bodyDict
169  )
170  );
171  }
172  }
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
177 
179 {}
180 
181 
182 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183 
186 {
187  return meshSolverPtr_->curPoints();
188 }
189 
190 
192 {
193  const Time& t = mesh().time();
194 
195  if (mesh().nPoints() != meshSolver_.points0().size())
196  {
198  << "The number of points in the mesh seems to have changed." << endl
199  << "In constant/polyMesh there are " << meshSolver_.points0().size()
200  << " points; in the current mesh there are " << mesh().nPoints()
201  << " points." << exit(FatalError);
202  }
203 
204  // Store the motion state at the beginning of the time-step
205  if (curTimeIndex_ != t.timeIndex())
206  {
207  newTime();
208  curTimeIndex_ = t.timeIndex();
209  }
210 
211  const scalar ramp = ramp_->value(t.value());
212 
214  {
215  g() =
216  ramp
218  }
219 
220  if (test_)
221  {
222  label nIter(coeffDict().lookup<label>("nIter"));
223 
224  for (label i=0; i<nIter; i++)
225  {
227  (
228  t.value(),
229  t.deltaTValue(),
230  scalarField(nDoF(), Zero),
232  );
233  }
234  }
235  else
236  {
238 
239  forAll(bodyMeshes_, bi)
240  {
241  const label bodyID = bodyMeshes_[bi].bodyID_;
242 
243  dictionary forcesDict;
244  forcesDict.add("type", functionObjects::forces::typeName);
245  forcesDict.add("patches", bodyMeshes_[bi].patches_);
246  forcesDict.add("rhoInf", rhoInf_);
247  forcesDict.add("rho", rhoName_);
248  forcesDict.add("CofR", vector::zero);
249 
250  functionObjects::forces f("forces", t, forcesDict);
251  f.calcForcesMoment();
252 
253  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
254  }
255 
257  (
258  t.value(),
259  t.deltaTValue(),
260  scalarField(nDoF(), Zero),
261  fx
262  );
263  }
264 
265  if (Pstream::master() && report())
266  {
267  forAll(bodyMeshes_, bi)
268  {
269  status(bodyMeshes_[bi].bodyID_);
270  }
271  }
272 
273  // Update the displacements
274  forAll(bodyMeshes_, bi)
275  {
276  forAllConstIter(labelHashSet, bodyMeshes_[bi].patchSet_, iter)
277  {
278  const label patchi = iter.key();
279 
280  const pointField patchPoints0
281  (
282  meshSolver_.pointDisplacement().boundaryField()[patchi]
283  .patchInternalField(meshSolver_.points0())
284  );
285 
286  meshSolver_.pointDisplacement().boundaryFieldRef()[patchi] ==
287  (
289  (
290  transform0(bodyMeshes_[bi].bodyID_),
291  patchPoints0
292  ) - patchPoints0
293  )();
294  }
295  }
296 
297  meshSolverPtr_->solve();
298 }
299 
300 
302 {
303  meshSolverPtr_->movePoints(points);
304 }
305 
306 
308 {
309  meshSolverPtr_->topoChange(map);
310 }
311 
312 
314 {
315  meshSolverPtr_->mapMesh(map);
316 }
317 
318 
320 (
321  const polyDistributionMap& map
322 )
323 {
324  meshSolverPtr_->distribute(map);
325 }
326 
327 
329 {
331  (
332  IOobject
333  (
334  "rigidBodyMotionState",
335  mesh().time().timeName(),
336  "uniform",
337  mesh(),
340  false
341  )
342  );
343 
344  state().write(dict);
345 
346  return
347  dict.regIOobject::writeObject
348  (
351  mesh().time().writeCompression(),
352  true
353  )
354  && motionSolver::write();
355 }
356 
357 
358 // ************************************************************************* //
void newTime()
Store the motion state at the beginning of the time-step.
Template class for intrusive linked lists.
Definition: ILList.H:50
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:1025
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
Virtual base class for displacement motion solver.
virtual bool write() const
Write motion state information for restart.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void status(const label bodyID) const
Report the status of the motion of the given body.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
bool foundObject(const word &name) const
Is the named Type found?
virtual bool write() const
Optionally write motion state information for restart.
Definition: motionSolver.C:128
Virtual base class for mesh motion solver.
Definition: motionSolver.H:56
const word & keyword() const
Return keyword.
Definition: motionSolver.H:131
virtual vector forceEff() const
Return the total force.
Definition: forces.C:1019
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
label nBodies() const
Return the number of bodies in the model (bodies().size())
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:143
const rigidBodyModelState & state() const
Return the motion state.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Macros for easy insertion into run-time selection tables.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1153
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
rigidBodyMeshMotionSolver(const word &name, const polyMesh &, const dictionary &dict)
Construct from polyMesh and dictionary.
Templated function that returns the corresponding 1 (one).
Definition: OneConstant.H:56
Spatial transformation functions for primitive fields.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
pointVectorField & pointDisplacement()
Return reference to the point motion displacement field.
Six degree of freedom motion for a rigid body.
stressControl lookup("compactNormalStress") >> compactNormalStress
const pointField & points
Pre-declare SubField and related Field type.
Definition: Field.H:56
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:47
virtual void topoChange(const polyTopoChangeMap &)
Update local data for topology changes.
void write(dictionary &dict) const
Write to dictionary.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
const Type & value() const
Return const reference to value.
word timeName
Definition: getTimeIndex.H:3
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:884
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
void solve(const scalar t, const scalar deltaT, const scalarField &tau, const Field< spatialVector > &fx)
Integrate velocities, orientation and position.
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
label nDoF() const
Return the number of degrees of freedom of the model.
const vector & g() const
Return the acceleration due to gravity.
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
motionSolver(const word &name, const polyMesh &mesh, const dictionary &, const word &type)
Construct from polyMesh and dictionary and type.
Definition: motionSolver.C:43
label patchi
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:203
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
void transformPoints(vectorField &, const spatialTransform &, const vectorField &)
Transform given vectorField of coordinates with the given spatialTransform.
label nPoints() const
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:76
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Calculates the forces and moments by integrating the pressure and skin-friction forces over a given l...
Definition: forces.H:204
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:137
A class for managing temporary objects.
Definition: PtrList.H:53
spatialTransform transform0(const label bodyID) const
Return the transformation of bodyID relative to the initial time.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
pointField & points0()
Return reference to the reference field.
virtual void solve()
Solve for motion.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
static autoPtr< motionSolver > New(const word &name, const polyMesh &, const dictionary &)
Select constructed from polyMesh and dictionary.
Definition: motionSolver.C:66
bool report() const
Return the report Switch.
Namespace for OpenFOAM.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864