acceleration.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) 2018-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 "fvMesh.H"
27 #include "fvMatrix.H"
28 #include "geometricOneField.H"
29 #include "acceleration.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
41  (
42  fvModel,
44  dictionary,
45  accelerationSource,
46  "accelerationSource"
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 void Foam::fv::acceleration::readCoeffs()
55 {
56  UName_ = coeffs().lookupOrDefault<word>("U", "U");
57 
58  velocity_ =
60  (
61  "velocity",
62  mesh().time().userUnits(),
64  coeffs()
65  );
66 }
67 
68 
69 template<class AlphaRhoFieldType>
70 void Foam::fv::acceleration::add
71 (
72  const AlphaRhoFieldType& alphaRho,
73  fvMatrix<vector>& eqn
74 ) const
75 {
76  const DimensionedField<scalar, volMesh>& V = mesh().V();
77 
78  const scalar t = mesh().time().value();
79  const scalar dt = mesh().time().deltaTValue();
80  const vector dU = velocity_->value(t) - velocity_->value(t - dt);
81  const vector a = dU/mesh().time().deltaTValue();
82 
83  const labelUList cells = set_.cells();
84 
85  vectorField& eqnSource = eqn.source();
86  forAll(cells, i)
87  {
88  eqnSource[cells[i]] -= V[cells[i]]*alphaRho[cells[i]]*a;
89  }
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
94 
96 (
97  const word& name,
98  const word& modelType,
99  const fvMesh& mesh,
100  const dictionary& dict
101 )
102 :
103  fvModel(name, modelType, mesh, dict),
104  set_(mesh, coeffs()),
105  UName_(word::null),
106  velocity_(nullptr)
107 {
108  readCoeffs();
109 }
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  return wordList(1, UName_);
117 }
118 
119 
121 (
122  const volVectorField& U,
123  fvMatrix<vector>& eqn
124 ) const
125 {
126  add(geometricOneField(), eqn);
127 }
128 
129 
131 (
132  const volScalarField& rho,
133  const volVectorField& U,
134  fvMatrix<vector>& eqn
135 ) const
136 {
137  add(rho, eqn);
138 }
139 
140 
142 (
143  const volScalarField& alpha,
144  const volScalarField& rho,
145  const volVectorField& U,
146  fvMatrix<vector>& eqn
147 ) const
148 {
149  add((alpha*rho)(), eqn);
150 }
151 
152 
154 {
155  set_.movePoints();
156  return true;
157 }
158 
159 
161 {
162  set_.topoChange(map);
163 }
164 
165 
167 {
168  set_.mapMesh(map);
169 }
170 
171 
173 {
174  set_.distribute(map);
175 }
176 
177 
179 {
180  if (fvModel::read(dict))
181  {
182  set_.read(coeffs());
183  readCoeffs();
184  return true;
185  }
186  else
187  {
188  return false;
189  }
190 }
191 
192 
193 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 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
This fvModel applies an explicit acceleration force to components of the velocity field.
Definition: acceleration.H:78
virtual bool movePoints()
Update for mesh motion.
Definition: acceleration.C:153
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
Definition: acceleration.C:121
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: acceleration.C:114
acceleration(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: acceleration.C:96
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: acceleration.C:160
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: acceleration.C:172
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: acceleration.C:178
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: acceleration.C:166
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
const cellShapeList & cells
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.
List< word > wordList
A List of words.
Definition: fileName.H:54
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Field< vector > vectorField
Specialisation of Field<T> for vector.
const dimensionSet dimVelocity
UList< label > labelUList
Definition: UList.H:65
labelList fv(nPoints)
dictionary dict