Loading...
Searching...
No Matches
sixDoFRigidBodyMotionSolver.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) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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::sixDoFRigidBodyMotionSolver::sixDoFRigidBodyMotionSolver
56(
57 const polyMesh& mesh,
58 const IOdictionary& dict
59)
60:
62 motion_
63 (
64 coeffDict(),
66 (
67 "sixDoFRigidBodyMotionState",
68 mesh.time().timeName(),
69 "uniform",
70 mesh
71 ).typeHeaderOk<IOdictionary>(true)
73 (
75 (
76 "sixDoFRigidBodyMotionState",
77 mesh.time().timeName(),
78 "uniform",
79 mesh,
80 IOobject::READ_IF_PRESENT,
81 IOobject::NO_WRITE,
82 IOobject::NO_REGISTER
83 )
84 )
85 : coeffDict(),
86 mesh.time()
87 ),
88 patches_(coeffDict().get<wordRes>("patches")),
89 patchSet_(mesh.boundaryMesh().patchSet(patches_)),
90 di_(coeffDict().get<scalar>("innerDistance")),
91 do_(coeffDict().get<scalar>("outerDistance")),
92 test_(coeffDict().getOrDefault("test", false)),
93 rhoInf_(1.0),
94 rhoName_(coeffDict().getOrDefault<word>("rho", "rho")),
95 scale_
96 (
98 (
99 "motionScale",
100 mesh.time().timeName(),
101 mesh,
102 IOobject::NO_READ,
103 IOobject::NO_WRITE,
104 IOobject::NO_REGISTER
105 ),
108 ),
109 curTimeIndex_(-1),
110 cOfGdisplacement_(coeffDict().getOrDefault<word>("cOfGdisplacement", "none"))
111{
112 if (rhoName_ == "rhoInf")
113 {
114 coeffDict().readEntry("rhoInf", rhoInf_);
115 }
116
117 // Calculate scaling factor everywhere
118 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
119
120 {
121 const pointMesh& pMesh = pointMesh::New(mesh);
122
123 pointPatchDist pDist(pMesh, patchSet_, points0());
124
125 // Scaling: 1 up to di then linear down to 0 at do away from patches
126 scale_.primitiveFieldRef() =
127 min
128 (
129 max
130 (
131 (do_ - pDist.primitiveField())/(do_ - di_),
132 scalar(0)
133 ),
134 scalar(1)
135 );
136
137 // Convert the scale function to a cosine
138 scale_.primitiveFieldRef() =
139 min
140 (
141 max
142 (
143 0.5
144 - 0.5
145 *cos(scale_.primitiveField()
147 scalar(0)
148 ),
149 scalar(1)
150 );
151
152 pointConstraints::New(pMesh).constrain(scale_);
153 scale_.write();
154 }
155}
156
157// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158
161{
162 return motion_;
163}
164
165
168{
169 tmp<pointField> newPoints
170 (
171 points0() + pointDisplacement_.primitiveField()
172 );
173
174 if (!moveAllCells())
175 {
176 auto ttransformedPts = tmp<pointField>::New(mesh().points());
177 auto& transformedPts = ttransformedPts.ref();
178
179 UIndirectList<point>(transformedPts, pointIDs()) =
180 pointField(newPoints.ref(), pointIDs());
181
182 return ttransformedPts;
183 }
184
185 return newPoints;
186}
187
188
190{
191 const Time& t = mesh().time();
192
193 if (mesh().nPoints() != points0().size())
194 {
196 << "The number of points in the mesh seems to have changed." << endl
197 << "In constant/polyMesh there are " << points0().size()
198 << " points; in the current mesh there are " << mesh().nPoints()
199 << " points." << exit(FatalError);
200 }
201
202 // Store the motion state at the beginning of the time-step
203 bool firstIter = false;
204 if (curTimeIndex_ != this->db().time().timeIndex())
205 {
206 motion_.newTime();
207 curTimeIndex_ = this->db().time().timeIndex();
208 firstIter = true;
209 }
210
212
213 if (db().time().foundObject<uniformDimensionedVectorField>("g"))
214 {
216 }
217 else
218 {
219 coeffDict().readIfPresent("g", g);
220 }
221
222 // const scalar ramp = clamp((this->db().time().value() - 5)/10, 0, 1);
223 const scalar ramp = 1.0;
224
225 if (test_)
226 {
227 motion_.update
228 (
229 firstIter,
230 ramp*(motion_.mass()*g.value()),
231 ramp*(motion_.mass()*(motion_.momentArm() ^ g.value())),
232 t.deltaTValue(),
233 t.deltaT0Value()
234 );
235 }
236 else
237 {
238 dictionary forcesDict;
239
240 forcesDict.add("type", functionObjects::forces::typeName);
241 forcesDict.add("patches", patches_);
242 forcesDict.add("rhoInf", rhoInf_);
243 forcesDict.add("rho", rhoName_);
244 forcesDict.add("CofR", motion_.centreOfRotation());
245
246 vector oldPos = motion_.centreOfRotation();
247
248 functionObjects::forces f("forces", db(), forcesDict);
249
250 f.calcForcesMoments();
251
252 motion_.update
253 (
254 firstIter,
255 ramp*(f.forceEff() + motion_.mass()*g.value()),
256 ramp
257 *(
258 f.momentEff()
259 + motion_.mass()*(motion_.momentArm() ^ g.value())
260 ),
261 t.deltaTValue(),
262 t.deltaT0Value()
263 );
264
265 if (cOfGdisplacement_ != "none")
266 {
267 if
268 (
269 db().time().foundObject<uniformDimensionedVectorField>
270 (
271 cOfGdisplacement_
272 )
273 )
274 {
275 auto& disp =
276 db().time().lookupObjectRef<uniformDimensionedVectorField>
277 (
278 cOfGdisplacement_
279 );
280
281 disp += (motion_.centreOfRotation() - oldPos);
282 }
283 }
284 }
285
286 // Update the displacements
287 pointDisplacement_.primitiveFieldRef() =
288 motion_.transform(points0(), scale_) - points0();
289
290 // Displacement has changed. Update boundary conditions
292 (
293 pointDisplacement_.mesh()
295}
296
297
299(
300 IOstreamOption streamOpt,
301 const bool writeOnProc
302) const
303{
305 (
307 (
308 "sixDoFRigidBodyMotionState",
309 mesh().time().timeName(),
310 "uniform",
311 mesh(),
315 )
316 );
317
318 motion_.state().write(dict);
319 return dict.regIOobject::writeObject(streamOpt, writeOnProc);
320}
321
322
324{
326 {
327 motion_.read(coeffDict());
328
329 return true;
330 }
331
332 return false;
333}
334
335
336// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
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.
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
scalar deltaT0Value() const noexcept
Return old time step value.
Definition TimeStateI.H:55
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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.
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,...
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...
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.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
pointVectorField pointDisplacement_
Point motion field.
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
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
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
virtual tmp< pointField > newPoints()
Provide new points for motion. Solves for motion.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
virtual bool read()
Read dynamicMeshDict dictionary.
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...
void constrain(GeometricField< Type, pointPatchField, pointMesh > &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),.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Calculation of distance to nearest patch for all points.
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
label nPoints() const noexcept
Number of mesh points.
6-DoF solid-body mesh motion solver for an 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.
const sixDoFRigidBodyMotion & motion() const
Return the six DoF motion object.
virtual bool read()
Read dynamicMeshDict dictionary.
Six degree of freedom motion for a rigid body.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
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
const labelList & pointIDs() const noexcept
The point ids (for cell set/zone subset).
Definition zoneMotion.H:102
bool moveAllCells() const noexcept
Move all cells?
Definition zoneMotion.H:110
#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
const pointField & points
label nPoints
word timeName
Definition getTimeIndex.H:3
constexpr scalar pi(M_PI)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
const dimensionSet dimless
Dimensionless.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
UniformDimensionedField< vector > uniformDimensionedVectorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimAcceleration
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
label timeIndex
labelList f(nPoints)
dictionary dict
Various UniformDimensionedField types.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))