momentumTransportModel.H
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) 2013-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 Class
25  Foam::momentumTransportModel
26 
27 Description
28  Abstract base class for turbulence models (RAS, LES and laminar).
29 
30 SourceFiles
31  momentumTransportModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef momentumTransportModel_H
36 #define momentumTransportModel_H
37 
38 #include "IOdictionary.H"
39 #include "primitiveFieldsFwd.H"
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 #include "fvMatricesFwd.H"
43 #include "geometricOneField.H"
44 #include "viscosity.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declarations
52 class fvMesh;
53 
54 /*---------------------------------------------------------------------------*\
55  Class momentumTransportModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 :
60  public IOdictionary
61 {
62 protected:
63 
64  // Protected data
65 
66  const Time& runTime_;
67  const fvMesh& mesh_;
68 
71  const surfaceScalarField& phi_;
72 
73  const viscosity& viscosity_;
74 
75 
76  // Protected member functions
77 
79  (
80  const objectRegistry& obr,
81  const word& group,
82  bool registerObject = false
83  );
84 
85 
86  template<class MomentumTransportModel>
88  (
89  const typename MomentumTransportModel::alphaField& alpha,
90  const typename MomentumTransportModel::rhoField& rho,
91  const volVectorField& U,
93  const surfaceScalarField& phi,
94  const viscosity& viscosity
95  );
96 
97 
98 public:
99 
100  //- Runtime type information
101  TypeName("momentumTransport");
102 
103 
104  // Constructors
105 
106  //- Construct from components
108  (
109  const volVectorField& U,
111  const surfaceScalarField& phi,
112  const viscosity& viscosity
113  );
114 
115  //- Disallow default bitwise copy construction
117 
118 
119  //- Destructor
120  virtual ~momentumTransportModel()
121  {}
122 
123 
124  // Member Functions
125 
126  //- Read model coefficients if they have changed
127  virtual bool read() = 0;
128 
129  const Time& time() const
130  {
131  return runTime_;
132  }
133 
134  const fvMesh& mesh() const
135  {
136  return mesh_;
137  }
138 
139  //- Const access to the coefficients dictionary
140  virtual const dictionary& coeffDict() const = 0;
141 
142  inline word groupName(const word& name) const
143  {
145  }
146 
147  //- Helper function to return the name of the turbulence G field
148  inline word GName() const
149  {
150  return typedName("G");
151  }
152 
153  //- Access function to velocity field
154  inline const volVectorField& U() const
155  {
156  return U_;
157  }
158 
159  //- Access function to phase flux field
160  inline const surfaceScalarField& alphaRhoPhi() const
161  {
162  return alphaRhoPhi_;
163  }
164 
165  //- Return the volumetric flux field
166  virtual tmp<surfaceScalarField> phi() const;
167 
168  //- Access function to fluid properties
169  const class Foam::viscosity& properties() const
170  {
171  return viscosity_;
172  }
173 
174  //- Return the near wall distance
175  const volScalarField::Boundary& y() const;
176 
177  //- Return the laminar viscosity
178  virtual tmp<volScalarField> nu() const
179  {
180  return this->viscosity_.nu();
181  }
182 
183  //- Return the laminar viscosity on patchi
184  virtual tmp<scalarField> nu(const label patchi) const
185  {
186  return this->viscosity_.nu(patchi);
187  }
188 
189  //- Return the turbulence viscosity
190  virtual tmp<volScalarField> nut() const = 0;
191 
192  //- Return the turbulence viscosity on patch
193  virtual tmp<scalarField> nut(const label patchi) const = 0;
194 
195  //- Return the effective viscosity
196  virtual tmp<volScalarField> nuEff() const = 0;
197 
198  //- Return the effective viscosity on patch
199  virtual tmp<scalarField> nuEff(const label patchi) const = 0;
200 
201  //- Return the turbulence kinetic energy
202  virtual tmp<volScalarField> k() const = 0;
203 
204  //- Return the turbulence kinetic energy dissipation rate
205  virtual tmp<volScalarField> epsilon() const = 0;
206 
207  //- Return the turbulence specific dissipation rate
208  virtual tmp<volScalarField> omega() const = 0;
209 
210  //- Return the stress tensor [m^2/s^2]
211  virtual tmp<volSymmTensorField> sigma() const = 0;
212 
213  //- Return the explicit stress correction matrix
214  // for the momentum equation.
215  // Combined with the implicit part by divDevTau
216  // Also creates and registers the devTauCorrFlux field
217  // which is cached in the matrix if fluxRequired.
219  (
220  const tmp<volTensorField>& devTauCorr,
222  ) const;
223 
224  //- Validate the fields after construction
225  // Update derived fields as required
226  virtual void validate();
227 
228  //- Predict the momentum transport coefficients if possible
229  // without solving momentum transport model equations
230  virtual void predict() = 0;
231 
232  //- Solve the momentum transport model equations
233  // and correct the momentum transport coefficients
234  virtual void correct() = 0;
235 
236 
237  // Member Operators
238 
239  //- Disallow default bitwise assignment
240  void operator=(const momentumTransportModel&) = delete;
241 };
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace Foam
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 #ifdef NoRepository
252 #endif
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #endif
257 
258 // ************************************************************************* //
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:346
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
static word groupName(Name name, const word &group)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
word groupName(const word &name) const
virtual void validate()
Validate the fields after construction.
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
word GName() const
Helper function to return the name of the turbulence G field.
const volVectorField & U() const
Access function to velocity field.
const volScalarField::Boundary & y() const
Return the near wall distance.
const surfaceScalarField & phi_
const class Foam::viscosity & properties() const
Access function to fluid properties.
virtual tmp< volScalarField > omega() const =0
Return the turbulence specific dissipation rate.
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
virtual bool read()=0
Read model coefficients if they have changed.
momentumTransportModel(const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
void operator=(const momentumTransportModel &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual tmp< fvVectorMatrix > divDevTauCorr(const tmp< volTensorField > &devTauCorr, volVectorField &U) const
Return the explicit stress correction matrix.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual ~momentumTransportModel()
Destructor.
virtual void correct()=0
Solve the momentum transport model equations.
virtual void predict()=0
Predict the momentum transport coefficients if possible.
const surfaceScalarField & alphaRhoPhi_
static autoPtr< MomentumTransportModel > New(const typename MomentumTransportModel::alphaField &alpha, const typename MomentumTransportModel::rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
virtual tmp< volSymmTensorField > sigma() const =0
Return the stress tensor [m^2/s^2].
TypeName("momentumTransport")
Runtime type information.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
Registry of regIOobjects.
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
A class for handling words, derived from string.
Definition: word.H:62
Forward declarations of fvMatrix specialisations.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
Foam::surfaceFields.