sixDoFAccelerationSource.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-2023 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 "fvMatrices.H"
28 #include "geometricOneField.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39  (
40  fvModel,
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 #include "None.H"
51 #include "Constant.H"
52 #include "Uniform.H"
53 #include "ZeroConstant.H"
54 #include "OneConstant.H"
55 #include "Polynomial1.H"
56 #include "Sine.H"
57 #include "Square.H"
58 #include "Table.H"
59 #include "UniformTable1.H"
60 #include "NonUniformTable1.H"
61 #include "EmbeddedTableReader.H"
62 #include "FoamTableReader.H"
63 #include "Scale.H"
64 #include "CodedFunction1.H"
65 
67 
68 template<>
69 const char* const avType::vsType::typeName = "vectorVector";
70 
71 template<>
72 const char* const avType::vsType::componentNames[] = {"x", "y", "z"};
73 
74 template<>
75 const avType avType::vsType::vsType::zero(avType::uniform(vector::uniform(0)));
76 
77 template<>
78 const avType avType::vsType::one(avType::uniform(vector::uniform(1)));
79 
80 template<>
81 const avType avType::vsType::max(avType::uniform(vector::uniform(vGreat)));
82 
83 template<>
84 const avType avType::vsType::min(avType::uniform(vector::uniform(-vGreat)));
85 
86 template<>
88 (
89  avType::uniform(vector::uniform(rootVGreat))
90 );
91 
92 template<>
94 (
95  avType::uniform(vector::uniform(-rootVGreat))
96 );
97 
98 namespace Foam
99 {
100 
102 
104  namespace TableReaders
105  {
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
113 
114 void Foam::fv::sixDoFAccelerationSource::readCoeffs()
115 {
116  UName_ = coeffs().lookupOrDefault<word>("U", "U");
117 
118  accelerations_.reset
119  (
121  (
122  "accelerations",
123  coeffs()
124  ).ptr()
125  );
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
130 
132 (
133  const word& name,
134  const word& modelType,
135  const fvMesh& mesh,
136  const dictionary& dict
137 )
138 :
139  fvModel(name, modelType, mesh, dict),
140  UName_(coeffs().lookupOrDefault<word>("U", "U")),
141  accelerations_(nullptr),
142  g_
143  (
144  mesh.foundObject<uniformDimensionedVectorField>("g")
145  ? dimensionedVector(mesh.lookupObject<uniformDimensionedVectorField>("g"))
147  )
148 {
149  readCoeffs();
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
156 {
157  return wordList(1, UName_);
158 }
159 
160 
162 (
163  fvMatrix<vector>& eqn,
164  const word& fieldName
165 ) const
166 {
167  addForce(geometricOneField(), geometricOneField(), eqn, fieldName);
168 }
169 
170 
172 (
173  const volScalarField& rho,
174  fvMatrix<vector>& eqn,
175  const word& fieldName
176 ) const
177 {
178  addForce(geometricOneField(), rho, eqn, fieldName);
179 }
180 
181 
183 (
184  const volScalarField& alpha,
185  const volScalarField& rho,
186  fvMatrix<vector>& eqn,
187  const word& fieldName
188 ) const
189 {
190  addForce(alpha, rho, eqn, fieldName);
191 }
192 
193 
195 {
196  return true;
197 }
198 
199 
201 {}
202 
203 
205 {}
206 
207 
209 (
210  const polyDistributionMap&
211 )
212 {}
213 
214 
216 {
217  if (fvModel::read(dict))
218  {
219  readCoeffs();
220  return true;
221  }
222  else
223  {
224  return false;
225  }
226 }
227 
228 
229 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Generic GeometricField class.
Reads an interpolation table from the values entry in OpenFOAM-format.
static const Form zero
Definition: VectorSpace.H:113
static const char *const componentNames[]
Definition: VectorSpace.H:112
static const Form one
Definition: VectorSpace.H:114
static const Form rootMax
Definition: VectorSpace.H:117
static const Form max
Definition: VectorSpace.H:115
static const Form rootMin
Definition: VectorSpace.H:118
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
static const Form min
Definition: VectorSpace.H:116
static const char *const typeName
Definition: VectorSpace.H:111
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
Solid-body 6-DoF acceleration source.
virtual bool movePoints()
Update for mesh motion.
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.
virtual void addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
sixDoFAccelerationSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
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.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
makeTableReader(Embedded, trvType)
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
makeFunction1s(trvType, nullArg)
const dimensionSet dimAcceleration
defineTableReader(trvType)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict
Foam::fv::sixDoFAccelerationSource::accelerationVectors avType