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-2020 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  false
79  ),
80  pointMesh::New(mesh),
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,
114  false
115  )
116  )
117  : coeffDict()
118  ),
119  test_(coeffDict().getOrDefault("test", false)),
120  rhoInf_(1.0),
121  rhoName_(coeffDict().getOrDefault<word>("rho", "rho")),
122  ramp_(Function1<scalar>::NewIfPresent("ramp", coeffDict())),
123  curTimeIndex_(-1)
124 {
125  if (rhoName_ == "rhoInf")
126  {
127  readEntry("rhoInf", rhoInf_);
128  }
129 
130  const dictionary& bodiesDict = coeffDict().subDict("bodies");
131 
132  for (const entry& dEntry : bodiesDict)
133  {
134  const keyType& bodyName = dEntry.keyword();
135  const dictionary& bodyDict = dEntry.dict();
136 
137  if (bodyDict.found("patches"))
138  {
139  const label bodyID = model_.bodyID(bodyName);
140 
141  if (bodyID == -1)
142  {
144  << "Body " << bodyName
145  << " has been merged with another body"
146  " and cannot be assigned a set of patches"
147  << exit(FatalError);
148  }
149 
150  bodyMeshes_.append
151  (
152  new bodyMesh
153  (
154  mesh,
155  bodyName,
156  bodyID,
157  bodyDict
158  )
159  );
160  }
161  }
162 
163  // Calculate scaling factor everywhere for each meshed body
164  forAll(bodyMeshes_, bi)
165  {
166  const pointMesh& pMesh = pointMesh::New(mesh);
167 
168  pointPatchDist pDist(pMesh, bodyMeshes_[bi].patchSet_, points0());
169 
170  pointScalarField& scale = bodyMeshes_[bi].weight_;
171 
172  // Scaling: 1 up to di then linear down to 0 at do away from patches
173  scale.primitiveFieldRef() =
174  min
175  (
176  max
177  (
178  (bodyMeshes_[bi].do_ - pDist.primitiveField())
179  /(bodyMeshes_[bi].do_ - bodyMeshes_[bi].di_),
180  scalar(0)
181  ),
182  scalar(1)
183  );
184 
185  // Convert the scale function to a cosine
186  scale.primitiveFieldRef() =
187  min
188  (
189  max
190  (
191  0.5
192  - 0.5
193  *cos(scale.primitiveField()
195  scalar(0)
196  ),
197  scalar(1)
198  );
199 
200  pointConstraints::New(pMesh).constrain(scale);
201  //scale.write();
202  }
203 }
204 
205 
206 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
207 
210 {
211  return points0() + pointDisplacement_.primitiveField();
212 }
213 
214 
216 {
217  const Time& t = mesh().time();
218 
219  if (mesh().nPoints() != points0().size())
220  {
222  << "The number of points in the mesh seems to have changed." << endl
223  << "In constant/polyMesh there are " << points0().size()
224  << " points; in the current mesh there are " << mesh().nPoints()
225  << " points." << exit(FatalError);
226  }
227 
228  // Store the motion state at the beginning of the time-step
229  if (curTimeIndex_ != this->db().time().timeIndex())
230  {
231  model_.newTime();
232  curTimeIndex_ = this->db().time().timeIndex();
233  }
234 
235  const scalar ramp = (ramp_ ? ramp_->value(t.value()) : 1.0);
236 
238  {
239  model_.g() =
240  ramp*t.lookupObject<uniformDimensionedVectorField>("g").value();
241  }
242 
243  if (test_)
244  {
245  const label nIter(coeffDict().get<label>("nIter"));
246 
247  for (label i=0; i<nIter; i++)
248  {
249  model_.solve
250  (
251  t.value(),
252  t.deltaTValue(),
253  scalarField(model_.nDoF(), Zero),
254  Field<spatialVector>(model_.nBodies(), Zero)
255  );
256  }
257  }
258  else
259  {
260  const label nIter(coeffDict().getOrDefault<label>("nIter", 1));
261 
262  for (label i=0; i<nIter; i++)
263  {
264  Field<spatialVector> fx(model_.nBodies(), Zero);
265 
266  forAll(bodyMeshes_, bi)
267  {
268  const label bodyID = bodyMeshes_[bi].bodyID_;
269 
270  dictionary forcesDict;
271  forcesDict.add("type", functionObjects::forces::typeName);
272  forcesDict.add("patches", bodyMeshes_[bi].patches_);
273  forcesDict.add("rhoInf", rhoInf_);
274  forcesDict.add("rho", rhoName_);
275  forcesDict.add("CofR", vector::zero);
276 
277  functionObjects::forces f("forces", db(), forcesDict);
278  f.calcForcesMoment();
279 
280  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
281  }
282 
283  model_.solve
284  (
285  t.value(),
286  t.deltaTValue(),
287  scalarField(model_.nDoF(), Zero),
288  fx
289  );
290  }
291  }
292 
293  if (Pstream::master() && model_.report())
294  {
295  forAll(bodyMeshes_, bi)
296  {
297  model_.status(bodyMeshes_[bi].bodyID_);
298  }
299  }
300 
301  // Update the displacements
302  if (bodyMeshes_.size() == 1)
303  {
304  pointDisplacement_.primitiveFieldRef() = model_.transformPoints
305  (
306  bodyMeshes_[0].bodyID_,
307  bodyMeshes_[0].weight_,
308  points0()
309  ) - points0();
310  }
311  else
312  {
313  labelList bodyIDs(bodyMeshes_.size());
314  List<const scalarField*> weights(bodyMeshes_.size());
315  forAll(bodyIDs, bi)
316  {
317  bodyIDs[bi] = bodyMeshes_[bi].bodyID_;
318  weights[bi] = &bodyMeshes_[bi].weight_;
319  }
320 
321  pointDisplacement_.primitiveFieldRef() =
322  model_.transformPoints(bodyIDs, weights, points0()) - points0();
323  }
324 
325  // Displacement has changed. Update boundary conditions
327  (
328  pointDisplacement_.mesh()
329  ).constrainDisplacement(pointDisplacement_);
330 }
331 
332 
334 (
335  IOstreamOption streamOpt,
336  const bool valid
337 ) const
338 {
339  // Force ASCII writing
340  streamOpt.format(IOstream::ASCII);
341 
343  (
344  IOobject
345  (
346  "rigidBodyMotionState",
347  mesh().time().timeName(),
348  "uniform",
349  mesh(),
352  false
353  )
354  );
355 
356  model_.state().write(dict);
357  return dict.regIOobject::writeObject(streamOpt, valid);
358 }
359 
360 
362 {
364  {
365  model_.read(coeffDict());
366 
367  return true;
368  }
369 
370  return false;
371 }
372 
373 
374 // ************************************************************************* //
Foam::rigidBodyMeshMotion::read
virtual bool read()
Read dynamicMeshDict dictionary.
Definition: rigidBodyMeshMotion.C:361
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
mathematicalConstants.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::functionObjects::forces
Calculates the forces and moments by integrating the pressure and skin-friction forces over a given l...
Definition: forces.H:236
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
rigidBodyMeshMotion.H
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:364
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::IOstreamOption::format
streamFormat format() const noexcept
Get the current stream format.
Definition: IOstreamOption.H:289
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
Foam::rigidBodyMeshMotion::solve
virtual void solve()
Solve for motion.
Definition: rigidBodyMeshMotion.C:215
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
polyMesh.H
Foam::rigidBodyMeshMotion::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: rigidBodyMeshMotion.C:209
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:86
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::keyType
A class for handling keywords in dictionaries.
Definition: keyType.H:60
Foam::UniformDimensionedField< vector >
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
pointConstraints.H
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::motionSolver::read
virtual bool read()
Read dynamicMeshDict dictionary.
Definition: motionSolver.C:227
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:53
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:50
pointPatchDist.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::pointPatchDist
Calculation of distance to nearest patch for all points.
Definition: pointPatchDist.H:54
uniformDimensionedFields.H
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
Foam::IOstreamOption::ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
forces.H
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
f
labelList f(nPoints)
Foam::TimeState::deltaTValue
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:42
Foam::List< label >
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::zero
static const Vector< Cmpt > zero
Definition: VectorSpace.H:115
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:708
Foam::spatialVector
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:50
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::pointConstraints::constrain
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
Definition: pointConstraintsTemplates.C:131
Foam::GeometricField< scalar, pointPatchField, pointMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::rigidBodyMeshMotion::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write state using stream options.
Definition: rigidBodyMeshMotion.C:334
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265