sixDoFRigidBodyMotionSolver.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) 2013-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 "pointDist.H"
29 #include "pointConstraints.H"
30 #include "timeIOdictionary.H"
32 #include "forces.H"
33 #include "mathematicalConstants.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 
43  (
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const word& name,
56  const polyMesh& mesh,
57  const dictionary& dict
58 )
59 :
60  displacementMotionSolver(name, mesh, dict, typeName),
62  (
63  coeffDict(),
65  (
66  "sixDoFRigidBodyMotionState",
67  mesh.time().name(),
68  "uniform",
69  mesh
70  ).headerOk()
72  (
73  IOobject
74  (
75  "sixDoFRigidBodyMotionState",
76  mesh.time().name(),
77  "uniform",
78  mesh,
79  IOobject::READ_IF_PRESENT,
80  IOobject::NO_WRITE,
81  false
82  )
83  )
84  : coeffDict()
85  ),
86  patches_(wordReList(coeffDict().lookup("patches"))),
87  patchSet_(mesh.boundaryMesh().patchSet(patches_)),
88  di_(coeffDict().lookup<scalar>("innerDistance")),
89  do_(coeffDict().lookup<scalar>("outerDistance")),
90  test_(coeffDict().lookupOrDefault<Switch>("test", false)),
91  rhoInf_(1.0),
92  rhoName_(coeffDict().lookupOrDefault<word>("rho", "rho")),
93  scale_
94  (
95  IOobject
96  (
97  "motionScale",
98  mesh.time().name(),
99  mesh,
100  IOobject::NO_READ,
101  IOobject::NO_WRITE,
102  false
103  ),
104  pointMesh::New(mesh),
106  ),
107  curTimeIndex_(-1)
108 {
109  if (rhoName_ == "rhoInf")
110  {
111  rhoInf_ = coeffDict().lookup<scalar>("rhoInf");
112  }
113 
114  // Calculate scaling factor everywhere
115  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
116 
117  {
118  const pointMesh& pMesh = pointMesh::New(mesh);
119 
120  pointDist pDist(pMesh, patchSet_, points0(), do_);
121 
122  // Scaling: 1 up to di then linear down to 0 at do away from patches
123  scale_.primitiveFieldRef() =
124  min
125  (
126  max
127  (
128  (do_ - pDist.primitiveField())/(do_ - di_),
129  scalar(0)
130  ),
131  scalar(1)
132  );
133 
134  // Convert the scale function to a cosine
135  scale_.primitiveFieldRef() =
136  min
137  (
138  max
139  (
140  0.5
141  - 0.5
142  *cos(scale_.primitiveField()
144  scalar(0)
145  ),
146  scalar(1)
147  );
148 
149  pointConstraints::New(pMesh).constrain(scale_);
150  scale_.write();
151  }
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
156 
158 {}
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
165 {
166  return points0() + pointDisplacement_.primitiveField();
167 }
168 
169 
171 {
172  const Time& t = mesh().time();
173 
174  if (mesh().nPoints() != points0().size())
175  {
177  << "The number of points in the mesh seems to have changed." << endl
178  << "In constant/polyMesh there are " << points0().size()
179  << " points; in the current mesh there are " << mesh().nPoints()
180  << " points." << exit(FatalError);
181  }
182 
183  // Store the motion state at the beginning of the time-stepbool
184  bool firstIter = false;
185  if (curTimeIndex_ != t.timeIndex())
186  {
187  newTime();
188  curTimeIndex_ = t.timeIndex();
189  firstIter = true;
190  }
191 
193 
194  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
195  {
196  g = mesh().lookupObject<uniformDimensionedVectorField>("g");
197  }
198  else if (coeffDict().found("g"))
199  {
200  coeffDict().lookup("g") >> g;
201  }
202 
203  // scalar ramp = min(max((t.value() - 5)/10, 0), 1);
204  scalar ramp = 1.0;
205 
206  if (test_)
207  {
208  update
209  (
210  firstIter,
211  ramp*(mass()*g.value()),
212  ramp*(mass()*(momentArm() ^ g.value())),
213  t.deltaTValue(),
214  t.deltaT0Value()
215  );
216  }
217  else
218  {
220  (
221  functionObjects::forces::typeName,
222  t,
223  dictionary
224  (
225  "type", functionObjects::forces::typeName,
226  "patches", patches_,
227  "rhoInf", rhoInf_,
228  "rho", rhoName_,
229  "CofR", centreOfRotation()
230  )
231  );
232 
233  f.calcForcesMoments();
234 
235  update
236  (
237  firstIter,
238  ramp*(f.forceEff() + mass()*g.value()),
239  ramp
240  *(
241  f.momentEff()
242  + mass()*(momentArm() ^ g.value())
243  ),
244  t.deltaTValue(),
245  t.deltaT0Value()
246  );
247  }
248 
249  // Update the displacements
250  pointDisplacement_.primitiveFieldRef() =
251  transform(points0(), scale_) - points0();
252 
253  // Displacement has changed. Update boundary conditions
255  (
256  pointDisplacement_.mesh()
257  ).constrainDisplacement(pointDisplacement_);
258 }
259 
260 
262 {
264  (
265  IOobject
266  (
267  "sixDoFRigidBodyMotionState",
268  mesh().time().name(),
269  "uniform",
270  mesh(),
273  false
274  )
275  );
276 
277  state().write(dict);
278 
279  return
280  dict.regIOobject::writeObject
281  (
284  mesh().time().writeCompression(),
285  true
286  )
288 }
289 
290 
291 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
scalar deltaT0Value() const
Return old time step value.
Definition: TimeStateI.H:40
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
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 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 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.
void constrain(PointField< Type > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Calculates the distance to the specified sets of patch and pointZone points or for all points.
Definition: pointDist.H:54
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
virtual bool write() const
Write points0 if the mesh topology changed.
pointField & points0()
Return reference to the reference field.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual bool write(const bool write=true) const
Write using setting from DB.
6-DoF solid-body mesh motion solver for an fvMesh.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
sixDoFRigidBodyMotionSolver(const word &name, const polyMesh &, const dictionary &dict)
Construct from polyMesh and dictionary.
virtual bool write() const
Write motion state information for restart.
virtual void solve()
Solve for motion.
Six degree of freedom motion for a rigid body.
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 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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
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 dimAcceleration
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dimensionedScalar cos(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict