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 "OneConstant.H"
37 #include "mathematicalConstants.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(rigidBodyMeshMotion, 0);
44 
46  (
47  motionSolver,
48  rigidBodyMeshMotion,
49  dictionary
50  );
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 Foam::rigidBodyMeshMotion::bodyMesh::bodyMesh
57 (
58  const polyMesh& mesh,
59  const word& name,
60  const label bodyID,
61  const dictionary& dict
62 )
63 :
64  name_(name),
65  bodyID_(bodyID),
66  patches_(dict.get<wordRes>("patches")),
67  patchSet_(mesh.boundaryMesh().patchSet(patches_)),
68  di_(dict.get<scalar>("innerDistance")),
69  do_(dict.get<scalar>("outerDistance")),
70  weight_
71  (
72  IOobject
73  (
74  name_ + ".motionScale",
75  mesh.time().timeName(),
76  mesh,
77  IOobject::NO_READ,
78  IOobject::NO_WRITE,
79  false
80  ),
81  pointMesh::New(mesh),
83  )
84 {}
85 
86 
87 Foam::rigidBodyMeshMotion::rigidBodyMeshMotion
88 (
89  const polyMesh& mesh,
90  const IOdictionary& dict
91 )
92 :
94  model_
95  (
96  mesh.time(),
97  coeffDict(),
98  IOobject
99  (
100  "rigidBodyMotionState",
101  mesh.time().timeName(),
102  "uniform",
103  mesh
104  ).typeHeaderOk<IOdictionary>(true)
105  ? IOdictionary
106  (
107  IOobject
108  (
109  "rigidBodyMotionState",
110  mesh.time().timeName(),
111  "uniform",
112  mesh,
115  false
116  )
117  )
118  : coeffDict()
119  ),
120  test_(coeffDict().getOrDefault("test", false)),
121  rhoInf_(1.0),
122  rhoName_(coeffDict().getOrDefault<word>("rho", "rho")),
123  ramp_(nullptr),
124  curTimeIndex_(-1)
125 {
126  if (rhoName_ == "rhoInf")
127  {
128  readEntry("rhoInf", rhoInf_);
129  }
130 
131  if (coeffDict().found("ramp"))
132  {
133  ramp_ = Function1<scalar>::New("ramp", coeffDict());
134  }
135  else
136  {
137  ramp_.reset(new Function1Types::OneConstant<scalar>("ramp"));
138  }
139 
140  const dictionary& bodiesDict = coeffDict().subDict("bodies");
141 
142  for (const entry& dEntry : bodiesDict)
143  {
144  const keyType& bodyName = dEntry.keyword();
145  const dictionary& bodyDict = dEntry.dict();
146 
147  if (bodyDict.found("patches"))
148  {
149  const label bodyID = model_.bodyID(bodyName);
150 
151  if (bodyID == -1)
152  {
154  << "Body " << bodyName
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  bodyName,
166  bodyID,
167  bodyDict
168  )
169  );
170  }
171  }
172 
173  // Calculate scaling factor everywhere for each meshed body
174  forAll(bodyMeshes_, bi)
175  {
176  const pointMesh& pMesh = pointMesh::New(mesh);
177 
178  pointPatchDist pDist(pMesh, bodyMeshes_[bi].patchSet_, points0());
179 
180  pointScalarField& scale = bodyMeshes_[bi].weight_;
181 
182  // Scaling: 1 up to di then linear down to 0 at do away from patches
183  scale.primitiveFieldRef() =
184  min
185  (
186  max
187  (
188  (bodyMeshes_[bi].do_ - pDist.primitiveField())
189  /(bodyMeshes_[bi].do_ - bodyMeshes_[bi].di_),
190  scalar(0)
191  ),
192  scalar(1)
193  );
194 
195  // Convert the scale function to a cosine
196  scale.primitiveFieldRef() =
197  min
198  (
199  max
200  (
201  0.5
202  - 0.5
203  *cos(scale.primitiveField()
205  scalar(0)
206  ),
207  scalar(1)
208  );
209 
210  pointConstraints::New(pMesh).constrain(scale);
211  //scale.write();
212  }
213 }
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
220 {
221  return points0() + pointDisplacement_.primitiveField();
222 }
223 
224 
226 {
227  const Time& t = mesh().time();
228 
229  if (mesh().nPoints() != points0().size())
230  {
232  << "The number of points in the mesh seems to have changed." << endl
233  << "In constant/polyMesh there are " << points0().size()
234  << " points; in the current mesh there are " << mesh().nPoints()
235  << " points." << exit(FatalError);
236  }
237 
238  // Store the motion state at the beginning of the time-step
239  if (curTimeIndex_ != this->db().time().timeIndex())
240  {
241  model_.newTime();
242  curTimeIndex_ = this->db().time().timeIndex();
243  }
244 
245  const scalar ramp = ramp_->value(t.value());
246 
248  {
249  model_.g() =
250  ramp*t.lookupObject<uniformDimensionedVectorField>("g").value();
251  }
252 
253  if (test_)
254  {
255  const label nIter(coeffDict().get<label>("nIter"));
256 
257  for (label i=0; i<nIter; i++)
258  {
259  model_.solve
260  (
261  t.value(),
262  t.deltaTValue(),
263  scalarField(model_.nDoF(), Zero),
264  Field<spatialVector>(model_.nBodies(), Zero)
265  );
266  }
267  }
268  else
269  {
270  const label nIter(coeffDict().getOrDefault<label>("nIter", 1));
271 
272  for (label i=0; i<nIter; i++)
273  {
274  Field<spatialVector> fx(model_.nBodies(), Zero);
275 
276  forAll(bodyMeshes_, bi)
277  {
278  const label bodyID = bodyMeshes_[bi].bodyID_;
279 
280  dictionary forcesDict;
281  forcesDict.add("type", functionObjects::forces::typeName);
282  forcesDict.add("patches", bodyMeshes_[bi].patches_);
283  forcesDict.add("rhoInf", rhoInf_);
284  forcesDict.add("rho", rhoName_);
285  forcesDict.add("CofR", vector::zero);
286 
287  functionObjects::forces f("forces", db(), forcesDict);
288  f.calcForcesMoment();
289 
290  fx[bodyID] = ramp*spatialVector(f.momentEff(), f.forceEff());
291  }
292 
293  model_.solve
294  (
295  t.value(),
296  t.deltaTValue(),
297  scalarField(model_.nDoF(), Zero),
298  fx
299  );
300  }
301  }
302 
303  if (Pstream::master() && model_.report())
304  {
305  forAll(bodyMeshes_, bi)
306  {
307  model_.status(bodyMeshes_[bi].bodyID_);
308  }
309  }
310 
311  // Update the displacements
312  if (bodyMeshes_.size() == 1)
313  {
314  pointDisplacement_.primitiveFieldRef() = model_.transformPoints
315  (
316  bodyMeshes_[0].bodyID_,
317  bodyMeshes_[0].weight_,
318  points0()
319  ) - points0();
320  }
321  else
322  {
323  labelList bodyIDs(bodyMeshes_.size());
324  List<const scalarField*> weights(bodyMeshes_.size());
325  forAll(bodyIDs, bi)
326  {
327  bodyIDs[bi] = bodyMeshes_[bi].bodyID_;
328  weights[bi] = &bodyMeshes_[bi].weight_;
329  }
330 
331  pointDisplacement_.primitiveFieldRef() =
332  model_.transformPoints(bodyIDs, weights, points0()) - points0();
333  }
334 
335  // Displacement has changed. Update boundary conditions
337  (
338  pointDisplacement_.mesh()
339  ).constrainDisplacement(pointDisplacement_);
340 }
341 
342 
344 (
345  IOstreamOption streamOpt,
346  const bool valid
347 ) const
348 {
349  // Force ASCII writing
350  streamOpt.format(IOstream::ASCII);
351 
353  (
354  IOobject
355  (
356  "rigidBodyMotionState",
357  mesh().time().timeName(),
358  "uniform",
359  mesh(),
362  false
363  )
364  );
365 
366  model_.state().write(dict);
367  return dict.regIOobject::writeObject(streamOpt, valid);
368 }
369 
370 
372 {
374  {
375  model_.read(coeffDict());
376 
377  return true;
378  }
379 
380  return false;
381 }
382 
383 
384 // ************************************************************************* //
Foam::rigidBodyMeshMotion::read
virtual bool read()
Read dynamicMeshDict dictionary.
Definition: rigidBodyMeshMotion.C:371
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:59
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:785
Foam::IOstreamOption::format
streamFormat format() const noexcept
Get the current stream format.
Definition: IOstreamOption.H:289
Foam::rigidBodyMeshMotion::solve
virtual void solve()
Solve for motion.
Definition: rigidBodyMeshMotion.C:225
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:219
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::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::Function1Types::OneConstant
Templated function that returns the corresponding 1 (one).
Definition: OneConstant.H:59
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::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:439
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
found
bool found
Definition: TABSMDCalcMethod2.H:32
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:372
forces.H
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
f
labelList f(nPoints)
OneConstant.H
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:248
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > 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::Function1::New
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict, const word &redirectType=word::null)
Selector.
Definition: Function1New.C:35
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:344
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265