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-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 
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, dict, typeName),
77  RBD::rigidBodyMotion
78  (
79  coeffDict(),
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  : coeffDict()
101  ),
102  test_(coeffDict().lookupOrDefault<Switch>("test", false)),
103  rhoInf_(1.0),
104  rhoName_(coeffDict().lookupOrDefault<word>("rho", "rho")),
105  ramp_(nullptr),
106  curTimeIndex_(-1),
107  meshSolverPtr_
108  (
110  (
111  name,
112  mesh,
114  (
115  IOobject
116  (
117  typedName("meshSolver"),
118  mesh.time().constant(),
119  mesh
120  ),
121  coeffDict().subDict("meshSolver")
122  )
123  )
124  ),
125  meshSolver_(refCast<displacementMotionSolver>(meshSolverPtr_()))
126 {
127  if (rhoName_ == "rhoInf")
128  {
129  rhoInf_ = coeffDict().lookup<scalar>("rhoInf");
130  }
131 
132  if (coeffDict().found("ramp"))
133  {
134  ramp_ = Function1<scalar>::New("ramp", dimTime, dimless, coeffDict());
135  }
136  else
137  {
138  ramp_ = new Function1s::OneConstant<scalar>("ramp");
139  }
140 
141  const dictionary& bodiesDict = coeffDict().subDict("bodies");
142 
143  forAllConstIter(IDLList<entry>, bodiesDict, iter)
144  {
145  const dictionary& bodyDict = iter().dict();
146 
147  if (bodyDict.found("patches"))
148  {
149  const label bodyID = this->bodyIndex(iter().keyword());
150 
151  if (bodyID == -1)
152  {
154  << "Body " << iter().keyword()
155  << " has been merged with another body"
156  " and cannot be assigned a set of patches"
157  << exit(FatalError);
158  }
159 
160  bodyMeshes_.append
161  (
162  new bodyMesh
163  (
164  mesh,
165  iter().keyword(),
166  bodyID,
167  bodyDict
168  )
169  );
170  }
171  }
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
176 
178 {}
179 
180 
181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 
185 {
186  return meshSolverPtr_->curPoints();
187 }
188 
189 
191 {
192  const Time& t = mesh().time();
193 
194  if (mesh().nPoints() != meshSolver_.points0().size())
195  {
197  << "The number of points in the mesh seems to have changed." << endl
198  << "In constant/polyMesh there are " << meshSolver_.points0().size()
199  << " points; in the current mesh there are " << mesh().nPoints()
200  << " points." << exit(FatalError);
201  }
202 
203  // Store the motion state at the beginning of the time-step
204  if (curTimeIndex_ != t.timeIndex())
205  {
206  newTime();
207  curTimeIndex_ = t.timeIndex();
208  }
209 
210  const scalar ramp = ramp_->value(t.value());
211 
212  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
213  {
214  g() =
215  ramp
216  *mesh().lookupObject<uniformDimensionedVectorField>("g").value();
217  }
218 
219  if (test_)
220  {
221  label nIter(coeffDict().lookup<label>("nIter"));
222 
223  for (label i=0; i<nIter; i++)
224  {
226  (
227  t.value(),
228  t.deltaTValue(),
229  scalarField(nDoF(), Zero),
230  Field<spatialVector>(nBodies(), Zero)
231  );
232  }
233  }
234  else
235  {
236  Field<spatialVector> fx(nBodies(), Zero);
237 
238  forAll(bodyMeshes_, bi)
239  {
240  const label bodyID = bodyMeshes_[bi].bodyIndex_;
241 
243  (
244  functionObjects::forces::typeName,
245  t,
246  dictionary
247  (
248  "type", functionObjects::forces::typeName,
249  "patches", bodyMeshes_[bi].patches_,
250  "rhoInf", rhoInf_,
251  "rho", rhoName_,
252  "CofR", vector::zero
253  )
254  );
255 
256  f.calcForcesMoments();
257 
258  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
259  }
260 
262  (
263  t.value(),
264  t.deltaTValue(),
265  scalarField(nDoF(), Zero),
266  fx
267  );
268  }
269 
270  if (Pstream::master() && report())
271  {
272  forAll(bodyMeshes_, bi)
273  {
274  status(bodyMeshes_[bi].bodyIndex_);
275  }
276  }
277 
278  // Update the displacements
279  forAll(bodyMeshes_, bi)
280  {
281  forAllConstIter(labelHashSet, bodyMeshes_[bi].patchSet_, iter)
282  {
283  const label patchi = iter.key();
284 
285  const pointField patchPoints0
286  (
287  meshSolver_.pointDisplacement().boundaryField()[patchi]
288  .patchInternalField(meshSolver_.points0())
289  );
290 
291  meshSolver_.pointDisplacement().boundaryFieldRef()[patchi] ==
292  (
294  (
295  transform0(bodyMeshes_[bi].bodyIndex_),
296  patchPoints0
297  ) - patchPoints0
298  )();
299  }
300  }
301 
302  meshSolverPtr_->solve();
303 }
304 
305 
307 {
308  meshSolverPtr_->movePoints(points);
309 }
310 
311 
313 {
314  meshSolverPtr_->topoChange(map);
315 }
316 
317 
319 {
320  meshSolverPtr_->mapMesh(map);
321 }
322 
323 
325 (
326  const polyDistributionMap& map
327 )
328 {
329  meshSolverPtr_->distribute(map);
330 }
331 
332 
334 {
336  (
337  IOobject
338  (
339  "rigidBodyMotionState",
340  mesh().time().name(),
341  "uniform",
342  mesh(),
345  false
346  )
347  );
348 
349  state().write(dict);
350 
351  return
353  && dict.regIOobject::writeObject
354  (
357  mesh().time().writeCompression(),
358  true
359  );
360 }
361 
362 
363 // ************************************************************************* //
#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 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 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
virtual bool write() const
Optionally write motion state information for restart.
Definition: motionSolver.C:128
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.
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.
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: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
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:129
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
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:176
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
error FatalError
labelList f(nPoints)
dictionary dict