Loading...
Searching...
No Matches
rigidBodyMeshMotionSolver.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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
31#include "polyMesh.H"
32#include "pointPatchDist.H"
33#include "pointConstraints.H"
35#include "forces.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
43
45 (
49 );
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55Foam::rigidBodyMeshMotionSolver::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{}
68
69
70Foam::rigidBodyMeshMotionSolver::rigidBodyMeshMotionSolver
71(
72 const polyMesh& mesh,
73 const IOdictionary& dict
74)
75:
77 model_
78 (
79 mesh.time(),
80 coeffDict(),
82 (
83 "rigidBodyMotionState",
84 mesh.time().timeName(),
85 "uniform",
86 mesh
87 ).typeHeaderOk<IOdictionary>(true)
89 (
91 (
92 "rigidBodyMotionState",
93 mesh.time().timeName(),
94 "uniform",
95 mesh,
96 IOobject::READ_IF_PRESENT,
97 IOobject::NO_WRITE,
98 IOobject::NO_REGISTER
99 )
100 )
101 : coeffDict()
102 ),
103 test_(coeffDict().getOrDefault("test", false)),
104 rhoInf_(1.0),
105 rhoName_(coeffDict().getOrDefault<word>("rho", "rho")),
106 curTimeIndex_(-1),
107 meshSolverPtr_
108 (
110 (
111 mesh,
113 (
115 (
116 "rigidBodyMotionSolver:meshSolver",
117 mesh.time().constant(),
118 mesh
119 ),
120 coeffDict().subDict("meshSolver")
121 )
122 )
123 ),
124 meshSolver_(refCast<displacementMotionSolver>(meshSolverPtr_()))
125{
126 if (rhoName_ == "rhoInf")
127 {
128 coeffDict().readEntry("rhoInf", rhoInf_);
129 }
130
131 const dictionary& bodiesDict = coeffDict().subDict("bodies");
132
133 for (const entry& dEntry : bodiesDict)
134 {
135 const keyType& bodyName = dEntry.keyword();
136 const dictionary& bodyDict = dEntry.dict();
137
138 if (bodyDict.found("patches"))
139 {
140 const label bodyID = model_.bodyID(bodyName);
141
142 if (bodyID == -1)
143 {
145 << "Body " << bodyName
146 << " has been merged with another body"
147 " and cannot be assigned a set of patches"
148 << exit(FatalError);
149 }
150
151 bodyMeshes_.append
152 (
153 new bodyMesh
154 (
155 mesh,
156 bodyName,
157 bodyID,
158 bodyDict
159 )
160 );
161 }
163}
164
165
166// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167
170{
171 return meshSolverPtr_->curPoints();
172}
173
174
176{
177 const Time& t = mesh().time();
178
179 if (mesh().nPoints() != meshSolver_.points0().size())
180 {
182 << "The number of points in the mesh seems to have changed." << endl
183 << "In constant/polyMesh there are " << meshSolver_.points0().size()
184 << " points; in the current mesh there are " << mesh().nPoints()
185 << " points." << exit(FatalError);
186 }
187
188 // Store the motion state at the beginning of the time-step
189 if (curTimeIndex_ != this->db().time().timeIndex())
190 {
191 model_.newTime();
192 curTimeIndex_ = this->db().time().timeIndex();
193 }
194
196 {
197 model_.g() = t.lookupObject<uniformDimensionedVectorField>("g").value();
198 }
199
200 if (test_)
201 {
202 const label nIter(coeffDict().get<label>("nIter"));
203
204 for (label i=0; i<nIter; i++)
205 {
206 model_.solve
207 (
208 t.value(),
209 t.deltaTValue(),
210 scalarField(model_.nDoF(), Zero),
211 Field<spatialVector>(model_.nBodies(), Zero)
212 );
213 }
214 }
215 else
216 {
217 Field<spatialVector> fx(model_.nBodies(), Zero);
218
219 forAll(bodyMeshes_, bi)
220 {
221 const label bodyID = bodyMeshes_[bi].bodyID_;
222
223 dictionary forcesDict;
224 forcesDict.add("type", functionObjects::forces::typeName);
225 forcesDict.add("patches", bodyMeshes_[bi].patches_);
226 forcesDict.add("rhoInf", rhoInf_);
227 forcesDict.add("rho", rhoName_);
228 forcesDict.add("CofR", vector::zero);
229
230 functionObjects::forces f("forces", db(), forcesDict);
231 f.calcForcesMoments();
232
233 fx[bodyID] = spatialVector(f.momentEff(), f.forceEff());
234 }
235
236 model_.solve
237 (
238 t.value(),
239 t.deltaTValue(),
240 scalarField(model_.nDoF(), Zero),
241 fx
242 );
243 }
244
245 if (Pstream::master() && model_.report())
246 {
247 forAll(bodyMeshes_, bi)
248 {
249 model_.status(bodyMeshes_[bi].bodyID_);
250 }
251 }
252
253 // Update the displacements
254 forAll(bodyMeshes_, bi)
255 {
256 for (const label patchi : bodyMeshes_[bi].patchSet_)
257 {
258 pointField patchPoints0
259 (
260 meshSolver_.pointDisplacement().boundaryField()[patchi]
261 .patchInternalField(meshSolver_.points0())
262 );
263
264 meshSolver_.pointDisplacement().boundaryFieldRef()[patchi] ==
265 (
266 model_.transformPoints
267 (
268 bodyMeshes_[bi].bodyID_,
269 patchPoints0
270 ) - patchPoints0
271 )();
273 }
274
275 meshSolverPtr_->solve();
276}
277
278
280(
281 IOstreamOption streamOpt,
282 const bool writeOnProc
283) const
284{
285 // Force ASCII writing
286 streamOpt.format(IOstreamOption::ASCII);
287
289 (
291 (
292 "rigidBodyMotionState",
293 mesh().time().timeName(),
294 "uniform",
295 mesh(),
299 )
300 );
301
302 model_.state().write(dict);
303 return dict.regIOobject::writeObject(streamOpt, writeOnProc);
304}
305
306
308{
309 if (motionSolver::read())
310 {
311 model_.read(coeffDict());
312
313 return true;
314 }
315
316 return false;
317}
318
321{
322 meshSolverPtr_->movePoints(points);
323}
324
325
327{
328 meshSolverPtr_->updateMesh(mpm);
329}
330
331
332// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
A simple container for options an IOstream can normally have.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
label bodyID(const word &name) const
Return the ID of the body with the given name.
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static const Form zero
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
friend class entry
Declare friendship with the entry class for IO.
Definition dictionary.H:338
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
const Type & value() const noexcept
Return const reference to value.
Virtual base class for displacement motion solver.
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition forces.H:321
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
A class for handling keywords in dictionaries.
Definition keyType.H:69
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
static autoPtr< motionSolver > New(const polyMesh &)
Select constructed from polyMesh.
motionSolver(const polyMesh &mesh)
Construct from polyMesh.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
virtual bool read()
Read dynamicMeshDict dictionary.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
label nPoints() const noexcept
Number of mesh points.
Rigid-body mesh motion solver for fvMesh.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write state using stream options.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
virtual void solve()
Solve for motion.
virtual bool read()
Read dynamicMeshDict dictionary.
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const pointField & points
label nPoints
word timeName
Definition getTimeIndex.H:3
Different types of constants.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
SpatialVector< scalar > spatialVector
SpatialVector of scalars.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
UniformDimensionedField< vector > uniformDimensionedVectorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
label timeIndex
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Various UniformDimensionedField types.