sixDoFAcceleration.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) 2015-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 
26 #include "sixDoFAcceleration.H"
27 #include "fvMatrices.H"
28 #include "geometricOneField.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
41  (
42  fvModel,
44  dictionary,
45  sixDoFAccelerationSource,
46  "sixDoFAccelerationSource"
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 void Foam::fv::sixDoFAcceleration::readCoeffs(const dictionary& dict)
55 {
56  UName_ = dict.lookupOrDefault<word>("U", "U");
57 
58  acceleration_.reset
59  (
61  (
62  "acceleration",
63  mesh().time().userUnits(),
65  dict
66  ).ptr()
67  );
68 
69  angularVelocity_.reset
70  (
72  (
73  "angularVelocity",
74  mesh().time().userUnits(),
76  dict
77  ).ptr()
78  );
79 
80  angularAcceleration_.reset
81  (
83  (
84  "angularAcceleration",
85  mesh().time().userUnits(),
87  dict
88  ).ptr()
89  );
90 }
91 
92 
93 template<class AlphaFieldType, class RhoFieldType>
94 void Foam::fv::sixDoFAcceleration::addForce
95 (
96  const AlphaFieldType& alpha,
97  const RhoFieldType& rho,
98  const volVectorField& U,
99  fvMatrix<vector>& eqn
100 ) const
101 {
102  const dimensionedVector a
103  (
104  "a",
106  acceleration_->value(mesh().time().value())
107  );
108 
109  // If gravitational force is present combine with the linear acceleration
110  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
111  {
114 
115  const uniformDimensionedScalarField& hRef =
117 
118  g = g_ - a;
119 
120  dimensionedScalar ghRef(- mag(g)*hRef);
121 
122  mesh().lookupObjectRef<volScalarField>("gh") = (g & mesh().C()) - ghRef;
123 
125  (g & mesh().Cf()) - ghRef;
126  }
127  // ... otherwise include explicitly in the momentum equation
128  else
129  {
130  eqn -= alpha*rho*a;
131  }
132 
133  const dimensionedVector Omega
134  (
135  "Omega",
137  angularVelocity_->value(mesh().time().value())
138  );
139 
140  const dimensionedVector dOmegaDT
141  (
142  "dOmegaDT",
144  angularAcceleration_->value(mesh().time().value())
145  );
146 
147  eqn -=
148  (
149  alpha*rho*(2*Omega ^ U) // Coriolis force
150  + alpha*rho*(Omega ^ (Omega ^ mesh().C())) // Centrifugal force
151  + alpha*rho*(dOmegaDT ^ mesh().C()) // Angular acceleration force
152  );
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157 
159 (
160  const word& name,
161  const word& modelType,
162  const fvMesh& mesh,
163  const dictionary& dict
164 )
165 :
166  fvModel(name, modelType, mesh, dict),
167  UName_(dict.lookupOrDefault<word>("U", "U")),
168  acceleration_(nullptr),
169  angularVelocity_(nullptr),
170  angularAcceleration_(nullptr),
171  g_
172  (
173  mesh.foundObject<uniformDimensionedVectorField>("g")
176  )
177 {
178  readCoeffs(coeffs(dict));
179 }
180 
181 
182 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183 
185 {
186  return wordList(1, UName_);
187 }
188 
189 
191 (
192  const volVectorField& U,
193  fvMatrix<vector>& eqn
194 ) const
195 {
196  addForce(geometricOneField(), geometricOneField(), U, eqn);
197 }
198 
199 
201 (
202  const volScalarField& rho,
203  const volVectorField& U,
204  fvMatrix<vector>& eqn
205 ) const
206 {
207  addForce(geometricOneField(), rho, U, eqn);
208 }
209 
210 
212 (
213  const volScalarField& alpha,
214  const volScalarField& rho,
215  const volVectorField& U,
216  fvMatrix<vector>& eqn
217 ) const
218 {
219  addForce(alpha, rho, U, eqn);
220 }
221 
222 
224 {
225  return true;
226 }
227 
228 
230 {}
231 
232 
234 {}
235 
236 
238 {}
239 
240 
242 {
243  if (fvModel::read(dict))
244  {
245  readCoeffs(coeffs(dict));
246  return true;
247  }
248  else
249  {
250  return false;
251  }
252 }
253 
254 
255 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const volVectorField & C() const
Return cell centres.
const surfaceVectorField & Cf() const
Return face centres.
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
Solid-body 6-DoF acceleration source.
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
sixDoFAcceleration(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Type & lookupObjectRef(const word &name) const
Lookup and return the object reference of the given Type.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
static const coefficient C("C", dimTemperature, 234.5)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
UniformDimensionedField< scalar > uniformDimensionedScalarField
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimAcceleration
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
UniformDimensionedField< vector > uniformDimensionedVectorField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
const unitConversion unitRadians
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict