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"
30 #include "makeFunction1s.H"
31 #include "makeTableReaders.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace fv
39 {
43  (
44  fvModel,
46  dictionary,
47  sixDoFAccelerationSource,
48  "sixDoFAccelerationSource"
49  );
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
57 
58 template<>
59 const char* const avType::vsType::typeName = "vectorVector";
60 
61 template<>
62 const char* const avType::vsType::componentNames[] = {"x", "y", "z"};
63 
64 template<>
65 const avType avType::vsType::vsType::zero(avType::uniform(vector::uniform(0)));
66 
67 template<>
68 const avType avType::vsType::one(avType::uniform(vector::uniform(1)));
69 
70 template<>
71 const avType avType::vsType::max(avType::uniform(vector::uniform(vGreat)));
72 
73 template<>
74 const avType avType::vsType::min(avType::uniform(vector::uniform(-vGreat)));
75 
76 template<>
78 (
79  avType::uniform(vector::uniform(rootVGreat))
80 );
81 
82 template<>
84 (
85  avType::uniform(vector::uniform(-rootVGreat))
86 );
87 
88 template<>
89 const avType avType::vsType::nan(avType::uniform(vector::uniform(NaN)));
90 
91 namespace Foam
92 {
95 }
96 
97 
98 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
99 
100 void Foam::fv::sixDoFAcceleration::readCoeffs()
101 {
102  UName_ = coeffs().lookupOrDefault<word>("U", "U");
103 
104  accelerations_.reset
105  (
107  (
108  "accelerations",
109  mesh().time().userUnits(),
110  unitNone,
111  coeffs()
112  ).ptr()
113  );
114 }
115 
116 
117 template<class AlphaFieldType, class RhoFieldType>
118 void Foam::fv::sixDoFAcceleration::addForce
119 (
120  const AlphaFieldType& alpha,
121  const RhoFieldType& rho,
122  const volVectorField& U,
123  fvMatrix<vector>& eqn
124 ) const
125 {
126  const Vector<vector> accelerations
127  (
128  accelerations_->value(mesh().time().value())
129  );
130 
131  // If gravitational force is present combine with the linear acceleration
132  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
133  {
135  mesh().lookupObjectRef<uniformDimensionedVectorField>("g");
136 
137  const uniformDimensionedScalarField& hRef =
138  mesh().lookupObject<uniformDimensionedScalarField>("hRef");
139 
140  g = g_ - dimensionedVector("a", dimAcceleration, accelerations.x());
141 
142  dimensionedScalar ghRef(- mag(g)*hRef);
143 
144  mesh().lookupObjectRef<volScalarField>("gh") = (g & mesh().C()) - ghRef;
145 
146  mesh().lookupObjectRef<surfaceScalarField>("ghf") =
147  (g & mesh().Cf()) - ghRef;
148  }
149  // ... otherwise include explicitly in the momentum equation
150  else
151  {
152  const dimensionedVector a("a", dimAcceleration, accelerations.x());
153  eqn -= alpha*rho*a;
154  }
155 
156  dimensionedVector Omega
157  (
158  "Omega",
159  dimensionSet(0, 0, -1, 0, 0),
160  accelerations.y()
161  );
162 
163  dimensionedVector dOmegaDT
164  (
165  "dOmegaDT",
166  dimensionSet(0, 0, -2, 0, 0),
167  accelerations.z()
168  );
169 
170  eqn -=
171  (
172  alpha*rho*(2*Omega ^ U) // Coriolis force
173  + alpha*rho*(Omega ^ (Omega ^ mesh().C())) // Centrifugal force
174  + alpha*rho*(dOmegaDT ^ mesh().C()) // Angular acceleration force
175  );
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
180 
182 (
183  const word& name,
184  const word& modelType,
185  const fvMesh& mesh,
186  const dictionary& dict
187 )
188 :
189  fvModel(name, modelType, mesh, dict),
190  UName_(coeffs().lookupOrDefault<word>("U", "U")),
191  accelerations_(nullptr),
192  g_
193  (
194  mesh.foundObject<uniformDimensionedVectorField>("g")
195  ? dimensionedVector(mesh.lookupObject<uniformDimensionedVectorField>("g"))
197  )
198 {
199  readCoeffs();
200 }
201 
202 
203 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204 
206 {
207  return wordList(1, UName_);
208 }
209 
210 
212 (
213  const volVectorField& U,
214  fvMatrix<vector>& eqn
215 ) const
216 {
217  addForce(geometricOneField(), geometricOneField(), U, eqn);
218 }
219 
220 
222 (
223  const volScalarField& rho,
224  const volVectorField& U,
225  fvMatrix<vector>& eqn
226 ) const
227 {
228  addForce(geometricOneField(), rho, U, eqn);
229 }
230 
231 
233 (
234  const volScalarField& alpha,
235  const volScalarField& rho,
236  const volVectorField& U,
237  fvMatrix<vector>& eqn
238 ) const
239 {
240  addForce(alpha, rho, U, eqn);
241 }
242 
243 
245 {
246  return true;
247 }
248 
249 
251 {}
252 
253 
255 {}
256 
257 
259 {}
260 
261 
263 {
264  if (fvModel::read(dict))
265  {
266  readCoeffs();
267  return true;
268  }
269  else
270  {
271  return false;
272  }
273 }
274 
275 
276 // ************************************************************************* //
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
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.
static const Form zero
Definition: VectorSpace.H:118
static const char *const componentNames[]
Definition: VectorSpace.H:117
static const Form one
Definition: VectorSpace.H:119
static const Form nan
Definition: VectorSpace.H:124
static const Form rootMax
Definition: VectorSpace.H:122
static const Form max
Definition: VectorSpace.H:120
static const Form rootMin
Definition: VectorSpace.H:123
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
static const Form min
Definition: VectorSpace.H:121
static const char *const typeName
Definition: VectorSpace.H:116
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
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:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
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...
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
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")
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:65
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
UniformDimensionedField< scalar > uniformDimensionedScalarField
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimAcceleration
const unitConversion unitNone
UniformDimensionedField< vector > uniformDimensionedVectorField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
makeFoamTableReaders(avType, nullArg)
makeFunction1s(avType, nullArg)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict
Foam::fv::sixDoFAcceleration::accelerationVectors avType