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 | Copyright (C) 2016 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 
26 #include "rigidBodyMeshMotion.H"
28 #include "polyMesh.H"
29 #include "pointPatchDist.H"
30 #include "pointConstraints.H"
32 #include "forces.H"
33 #include "mathematicalConstants.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(rigidBodyMeshMotion, 0);
40 
42  (
43  motionSolver,
44  rigidBodyMeshMotion,
45  dictionary
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::rigidBodyMeshMotion::bodyMesh::bodyMesh
53 (
54  const polyMesh& mesh,
55  const word& name,
56  const label bodyID,
57  const dictionary& dict
58 )
59 :
60  name_(name),
61  bodyID_(bodyID),
62  patches_(wordReList(dict.lookup("patches"))),
63  patchSet_(mesh.boundaryMesh().patchSet(patches_)),
64  di_(readScalar(dict.lookup("innerDistance"))),
65  do_(readScalar(dict.lookup("outerDistance"))),
66  weight_
67  (
68  IOobject
69  (
70  name_ + ".motionScale",
71  mesh.time().timeName(),
72  mesh,
73  IOobject::NO_READ,
74  IOobject::NO_WRITE,
75  false
76  ),
77  pointMesh::New(mesh),
78  dimensionedScalar("zero", dimless, 0.0)
79  )
80 {}
81 
82 
83 Foam::rigidBodyMeshMotion::rigidBodyMeshMotion
84 (
85  const polyMesh& mesh,
86  const IOdictionary& dict
87 )
88 :
89  displacementMotionSolver(mesh, dict, typeName),
90  model_
91  (
92  coeffDict(),
93  IOobject
94  (
95  "rigidBodyMotionState",
96  mesh.time().timeName(),
97  "uniform",
98  mesh
99  ).headerOk()
100  ? IOdictionary
101  (
102  IOobject
103  (
104  "rigidBodyMotionState",
105  mesh.time().timeName(),
106  "uniform",
107  mesh,
110  false
111  )
112  )
113  : coeffDict()
114  ),
115  test_(coeffDict().lookupOrDefault<Switch>("test", false)),
116  rhoInf_(1.0),
117  rhoName_(coeffDict().lookupOrDefault<word>("rho", "rho")),
118  curTimeIndex_(-1)
119 {
120  if (rhoName_ == "rhoInf")
121  {
122  rhoInf_ = readScalar(coeffDict().lookup("rhoInf"));
123  }
124 
125  const dictionary& bodiesDict = coeffDict().subDict("bodies");
126 
127  forAllConstIter(IDLList<entry>, bodiesDict, iter)
128  {
129  const dictionary& bodyDict = iter().dict();
130 
131  if (bodyDict.found("patches"))
132  {
133  bodyMeshes_.append
134  (
135  new bodyMesh
136  (
137  mesh,
138  iter().keyword(),
139  model_.bodyID(iter().keyword()),
140  bodyDict
141  )
142  );
143  }
144  }
145 
146  // Calculate scaling factor everywhere for each meshed body
147  forAll(bodyMeshes_, bi)
148  {
149  const pointMesh& pMesh = pointMesh::New(mesh);
150 
151  pointPatchDist pDist(pMesh, bodyMeshes_[bi].patchSet_, points0());
152 
153  pointScalarField& scale = bodyMeshes_[bi].weight_;
154 
155  // Scaling: 1 up to di then linear down to 0 at do away from patches
156  scale.primitiveFieldRef() =
157  min
158  (
159  max
160  (
161  (bodyMeshes_[bi].do_ - pDist.primitiveField())
162  /(bodyMeshes_[bi].do_ - bodyMeshes_[bi].di_),
163  scalar(0)
164  ),
165  scalar(1)
166  );
167 
168  // Convert the scale function to a cosine
169  scale.primitiveFieldRef() =
170  min
171  (
172  max
173  (
174  0.5
175  - 0.5
176  *cos(scale.primitiveField()
178  scalar(0)
179  ),
180  scalar(1)
181  );
182 
183  pointConstraints::New(pMesh).constrain(scale);
184  //scale.write();
185  }
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
190 
192 {}
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
199 {
201 }
202 
203 
205 {
206  const Time& t = mesh().time();
207 
208  if (mesh().nPoints() != points0().size())
209  {
211  << "The number of points in the mesh seems to have changed." << endl
212  << "In constant/polyMesh there are " << points0().size()
213  << " points; in the current mesh there are " << mesh().nPoints()
214  << " points." << exit(FatalError);
215  }
216 
217  // Store the motion state at the beginning of the time-step
218  if (curTimeIndex_ != this->db().time().timeIndex())
219  {
220  model_.newTime();
221  curTimeIndex_ = this->db().time().timeIndex();
222  }
223 
224  if (db().foundObject<uniformDimensionedVectorField>("g"))
225  {
226  model_.g() =
228  }
229 
230  if (test_)
231  {
232  label nIter(readLabel(coeffDict().lookup("nIter")));
233 
234  for (label i=0; i<nIter; i++)
235  {
236  model_.solve
237  (
238  t.deltaTValue(),
239  scalarField(model_.nDoF(), Zero),
241  );
242  }
243  }
244  else
245  {
246  Field<spatialVector> fx(model_.nBodies(), Zero);
247 
248  forAll(bodyMeshes_, bi)
249  {
250  const label bodyID = bodyMeshes_[bi].bodyID_;
251 
252  dictionary forcesDict;
253  forcesDict.add("type", functionObjects::forces::typeName);
254  forcesDict.add("patches", bodyMeshes_[bi].patches_);
255  forcesDict.add("rhoInf", rhoInf_);
256  forcesDict.add("rho", rhoName_);
257  forcesDict.add("CofR", vector::zero);
258 
259  functionObjects::forces f("forces", db(), forcesDict);
260  f.calcForcesMoment();
261 
262  fx[bodyID] = spatialVector(f.momentEff(), f.forceEff());
263  }
264 
265  model_.solve
266  (
267  t.deltaTValue(),
268  scalarField(model_.nDoF(), Zero),
269  fx
270  );
271  }
272 
273  if (Pstream::master() && model_.report())
274  {
275  forAll(bodyMeshes_, bi)
276  {
277  model_.status(bodyMeshes_[bi].bodyID_);
278  }
279  }
280 
281  // Update the displacements
282  if (bodyMeshes_.size() == 1)
283  {
285  (
286  bodyMeshes_[0].bodyID_,
287  bodyMeshes_[0].weight_,
288  points0()
289  ) - points0();
290  }
291  else
292  {
293  labelList bodyIDs(bodyMeshes_.size());
294  List<const scalarField*> weights(bodyMeshes_.size());
295  forAll(bodyIDs, bi)
296  {
297  bodyIDs[bi] = bodyMeshes_[bi].bodyID_;
298  weights[bi] = &bodyMeshes_[bi].weight_;
299  }
300 
302  model_.transformPoints(bodyIDs, weights, points0()) - points0();
303  }
304 
305  // Displacement has changed. Update boundary conditions
307  (
309  ).constrainDisplacement(pointDisplacement_);
310 }
311 
312 
314 (
318 ) const
319 {
321  (
322  IOobject
323  (
324  "rigidBodyMotionState",
325  mesh().time().timeName(),
326  "uniform",
327  mesh(),
330  false
331  )
332  );
333 
334  model_.state().write(dict);
335  return dict.regIOobject::write();
336 }
337 
338 
340 {
342  {
343  model_.read(coeffDict());
344 
345  return true;
346  }
347  else
348  {
349  return false;
350  }
351 }
352 
353 
354 // ************************************************************************* //
void newTime()
Store the motion state at the beginning of the time-step.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
const Time & time() const
Return time.
const rigidBodyModelState & state() const
Return the motion state.
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
dictionary dict
Virtual base class for displacement motion solver.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
pointVectorField pointDisplacement_
Point motion field.
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
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:129
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
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:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label bodyID(const word &name) const
Return the ID of the body with the given name.
void status(const label bodyID) const
Report the status of the motion of the given body.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
virtual vector forceEff() const
Return the total force.
Definition: forces.C:887
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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.
Definition: Switch.H:60
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write state using given format, version and compression.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
bool read(const dictionary &dict)
Read coefficients dictionary and update system parameters,.
label size() const
Return number of elements in list.
Definition: DLListBaseI.H:77
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
label nBodies() const
Return the number of bodies in the model (bodies().size())
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:737
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
static const pointMesh & New(const polyMesh &mesh)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
virtual bool read()
Read dynamicMeshDict dictionary.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
label nDoF() const
Return the number of degrees of freedom of the model.
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar cos(const dimensionedScalar &ds)
virtual void solve()
Solve for motion.
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
const Time & time() const
Return time.
Definition: IOobject.C:227
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Definition: spatialVector.H:47
Calculation of distance to nearest patch for all points.
void solve(scalar deltaT, const scalarField &tau, const Field< spatialVector > &fx)
Integrate velocities, orientation and position.
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:123
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
void write(dictionary &dict) const
Write to dictionary.
const vector & g() const
Return the acceleration due to gravity.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
word timeName
Definition: getTimeIndex.H:3
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:747
static const zero Zero
Definition: zero.H:91
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
virtual bool read()
Read dynamicMeshDict dictionary.
Definition: motionSolver.C:189
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
pointField & points0()
Return reference to the reference field.
tmp< pointField > transformPoints(const label bodyID, const pointField &initialPoints) const
Transform the given initial pointField of the specified body.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:893
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
IOdictionary(const IOobject &)
Construct given an IOobject.
Definition: IOdictionary.C:45
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
label timeIndex
Definition: getTimeIndex.H:4
Version number type.
Definition: IOstream.H:96
label nPoints() const
bool headerOk()
Read and check header info.
Definition: IOobject.C:400
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
This function object calculates the forces and moments by integrating the pressure and skin-friction ...
Definition: forces.H:196
A class for managing temporary objects.
Definition: PtrList.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
bool report() const
Return the report Switch.
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451