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-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 
27 #include "polyMesh.H"
28 #include "polyTopoChangeMap.H"
29 #include "pointConstraints.H"
30 #include "timeIOdictionary.H"
32 #include "forces.H"
33 #include "OneConstant.H"
34 #include "mathematicalConstants.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
42 
44  (
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 Foam::rigidBodyMeshMotionSolver::bodyMesh::bodyMesh
55 (
56  const polyMesh& mesh,
57  const word& name,
58  const label bodyID,
59  const dictionary& dict
60 )
61 :
62  name_(name),
63  bodyIndex_(bodyID),
64  patches_(wordReList(dict.lookup("patches"))),
65  patchSet_(mesh.boundaryMesh().patchSet(patches_))
66 {}
67 
68 
70 (
71  const word& name,
72  const polyMesh& mesh,
73  const dictionary& dict
74 )
75 :
76  motionSolver(name, mesh, typeName),
77  RBD::rigidBodyMotion
78  (
79  dict,
81  (
82  "rigidBodyMotionState",
83  mesh.time().name(),
84  "uniform",
85  mesh
86  ).headerOk()
88  (
89  IOobject
90  (
91  "rigidBodyMotionState",
92  mesh.time().name(),
93  "uniform",
94  mesh,
95  IOobject::READ_IF_PRESENT,
96  IOobject::NO_WRITE,
97  false
98  )
99  )
100  : dict
101  ),
102  test_(dict.lookupOrDefault<Switch>("test", false)),
103  nIter_(test_ ? dict.lookup<label>("nIter") : 0),
104  rhoInf_(1.0),
105  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
106  ramp_(nullptr),
107  curTimeIndex_(-1),
108  meshSolverPtr_
109  (
111  (
112  name,
113  mesh,
115  (
116  IOobject
117  (
118  typedName("meshSolver"),
119  mesh.time().constant(),
120  mesh
121  ),
122  dict.subDict("meshSolver")
123  )
124  )
125  ),
126  meshSolver_(refCast<displacementMotionSolver>(meshSolverPtr_()))
127 {
128  if (rhoName_ == "rhoInf")
129  {
130  rhoInf_ = dict.lookup<scalar>("rhoInf");
131  }
132 
133  if (dict.found("ramp"))
134  {
135  ramp_ = Function1<scalar>::New("ramp", dimTime, dimless, dict);
136  }
137  else
138  {
139  ramp_ = new Function1s::OneConstant<scalar>("ramp");
140  }
141 
142  const dictionary& bodiesDict = dict.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->bodyIndex(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 
213  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
214  {
215  g() =
216  ramp
218  }
219 
220  if (test_)
221  {
222  for (label i=0; i<nIter_; i++)
223  {
225  (
226  t.value(),
227  t.deltaTValue(),
228  scalarField(nDoF(), Zero),
229  Field<spatialVector>(nBodies(), Zero)
230  );
231  }
232  }
233  else
234  {
235  Field<spatialVector> fx(nBodies(), Zero);
236 
237  forAll(bodyMeshes_, bi)
238  {
239  const label bodyID = bodyMeshes_[bi].bodyIndex_;
240 
242  (
243  functionObjects::forces::typeName,
244  t,
246  (
247  "type", functionObjects::forces::typeName,
248  "patches", bodyMeshes_[bi].patches_,
249  "rhoInf", rhoInf_,
250  "rho", rhoName_,
251  "CofR", vector::zero
252  )
253  );
254 
255  f.calcForcesMoments();
256 
257  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
258  }
259 
261  (
262  t.value(),
263  t.deltaTValue(),
264  scalarField(nDoF(), Zero),
265  fx
266  );
267  }
268 
269  if (Pstream::master() && report())
270  {
271  forAll(bodyMeshes_, bi)
272  {
273  status(bodyMeshes_[bi].bodyIndex_);
274  }
275  }
276 
277  // Update the displacements
278  forAll(bodyMeshes_, bi)
279  {
280  forAllConstIter(labelHashSet, bodyMeshes_[bi].patchSet_, iter)
281  {
282  const label patchi = iter.key();
283 
284  const pointField patchPoints0
285  (
286  meshSolver_.pointDisplacement().boundaryField()[patchi]
287  .patchInternalField(meshSolver_.points0())
288  );
289 
290  meshSolver_.pointDisplacement().boundaryFieldRef()[patchi] ==
291  (
293  (
294  transform0(bodyMeshes_[bi].bodyIndex_),
295  patchPoints0
296  ) - patchPoints0
297  )();
298  }
299  }
300 
301  meshSolverPtr_->solve();
302 }
303 
304 
306 {
307  meshSolverPtr_->movePoints(points);
308 }
309 
310 
312 {
313  meshSolverPtr_->topoChange(map);
314 }
315 
316 
318 {
319  meshSolverPtr_->mapMesh(map);
320 }
321 
322 
324 (
325  const polyDistributionMap& map
326 )
327 {
328  meshSolverPtr_->distribute(map);
329 }
330 
331 
333 {
335  (
336  IOobject
337  (
338  "rigidBodyMotionState",
339  mesh().time().name(),
340  "uniform",
341  mesh(),
344  false
345  )
346  );
347 
348  state().write(dict);
349 
350  return
352  && dict.regIOobject::writeObject
353  (
356  mesh().time().writeCompression(),
357  true
358  );
359 }
360 
361 
362 // ************************************************************************* //
#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 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
Template class for intrusive linked lists.
Definition: ILList.H:67
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
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
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.
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
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.
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
virtual bool write() const
Optionally write motion state information for restart.
Definition: motionSolver.C:134
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.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
rigidBodyMeshMotionSolver(const word &name, const polyMesh &, const dictionary &dict)
Construct from polyMesh and dictionary.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual bool write() const
Write motion state information for restart.
virtual void solve()
Solve for motion.
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
label patchi
const pointField & points
label nPoints
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
void transformPoints(vectorField &, const spatialTransform &, const vectorField &)
Transform given vectorField of coordinates with the given spatialTransform.
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)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
const dimensionSet dimTime
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:47
defineTypeNameAndDebug(combustionModel, 0)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
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
labelList f(nPoints)
dictionary dict