turbulenceModel.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-2019 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::turbulenceModel
26 
27 Description
28  Abstract base class for turbulence models (RAS, LES and laminar).
29 
30 SourceFiles
31  turbulenceModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef turbulenceModel_H
36 #define turbulenceModel_H
37 
38 #include "IOdictionary.H"
39 #include "primitiveFieldsFwd.H"
40 #include "volFieldsFwd.H"
41 #include "surfaceFieldsFwd.H"
42 #include "fvMatricesFwd.H"
43 #include "nearWallDist.H"
44 #include "geometricOneField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declarations
52 class fvMesh;
53 
54 /*---------------------------------------------------------------------------*\
55  Class turbulenceModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class turbulenceModel
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  //- Near wall distance boundary field
75 
76 
77 public:
78 
79  //- Runtime type information
80  TypeName("turbulenceModel");
81 
82  //- Default name of the turbulence properties dictionary
83  static const word propertiesName;
84 
85 
86  // Constructors
87 
88  //- Construct from components
90  (
91  const volVectorField& U,
93  const surfaceScalarField& phi,
94  const word& propertiesName
95  );
96 
97  //- Disallow default bitwise copy construction
98  turbulenceModel(const turbulenceModel&) = delete;
99 
100 
101  //- Destructor
102  virtual ~turbulenceModel()
103  {}
104 
105 
106  // Member Functions
107 
108  //- Read model coefficients if they have changed
109  virtual bool read() = 0;
111  const Time& time() const
112  {
113  return runTime_;
114  }
116  const fvMesh& mesh() const
117  {
118  return mesh_;
119  }
120 
121  //- Const access to the coefficients dictionary
122  virtual const dictionary& coeffDict() const = 0;
123 
124  //- Helper function to return the name of the turbulence G field
125  inline word GName() const
126  {
127  return word(type() + ":G");
128  }
129 
130  //- Access function to velocity field
131  inline const volVectorField& U() const
132  {
133  return U_;
134  }
135 
136  //- Access function to phase flux field
137  inline const surfaceScalarField& alphaRhoPhi() const
138  {
139  return alphaRhoPhi_;
140  }
141 
142  //- Return the volumetric flux field
143  virtual tmp<surfaceScalarField> phi() const;
144 
145  //- Return the near wall distances
146  const nearWallDist& y() const
147  {
148  return y_;
149  }
150 
151  //- Return the laminar viscosity
152  virtual tmp<volScalarField> nu() const = 0;
153 
154  //- Return the laminar viscosity on patch
155  virtual tmp<scalarField> nu(const label patchi) const = 0;
156 
157  //- Return the turbulence viscosity
158  virtual tmp<volScalarField> nut() const = 0;
159 
160  //- Return the turbulence viscosity on patch
161  virtual tmp<scalarField> nut(const label patchi) const = 0;
162 
163  //- Return the effective viscosity
164  virtual tmp<volScalarField> nuEff() const = 0;
165 
166  //- Return the effective viscosity on patch
167  virtual tmp<scalarField> nuEff(const label patchi) const = 0;
168 
169  //- Return the laminar dynamic viscosity
170  virtual tmp<volScalarField> mu() const = 0;
171 
172  //- Return the laminar dynamic viscosity on patch
173  virtual tmp<scalarField> mu(const label patchi) const = 0;
174 
175  //- Return the turbulence dynamic viscosity
176  virtual tmp<volScalarField> mut() const = 0;
177 
178  //- Return the turbulence dynamic viscosity on patch
179  virtual tmp<scalarField> mut(const label patchi) const = 0;
180 
181  //- Return the effective dynamic viscosity
182  virtual tmp<volScalarField> muEff() const = 0;
183 
184  //- Return the effective dynamic viscosity on patch
185  virtual tmp<scalarField> muEff(const label patchi) const = 0;
186 
187  //- Return the turbulence kinetic energy
188  virtual tmp<volScalarField> k() const = 0;
189 
190  //- Return the turbulence kinetic energy dissipation rate
191  virtual tmp<volScalarField> epsilon() const = 0;
192 
193  //- Return the Reynolds stress tensor
194  virtual tmp<volSymmTensorField> R() const = 0;
195 
196  //- Validate the turbulence fields after construction
197  // Update derived fields as required
198  virtual void validate();
199 
200  //- Solve the turbulence equations and correct the turbulence viscosity
201  virtual void correct() = 0;
202 
203 
204  // Member Operators
205 
206  //- Disallow default bitwise assignment
207  void operator=(const turbulenceModel&) = delete;
208 };
209 
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 } // End namespace Foam
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #endif
218 
219 // ************************************************************************* //
virtual void validate()
Validate the turbulence fields after construction.
const fvMesh & mesh() const
const surfaceScalarField & phi_
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
const volVectorField & U() const
Access function to velocity field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
void operator=(const turbulenceModel &)=delete
Disallow default bitwise assignment.
virtual ~turbulenceModel()
Destructor.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > muEff() const =0
Return the effective dynamic viscosity.
virtual tmp< volSymmTensorField > R() const =0
Return the Reynolds stress tensor.
virtual tmp< volScalarField > mut() const =0
Return the turbulence dynamic viscosity.
Generic GeometricField class.
TypeName("turbulenceModel")
Runtime type information.
const fvMesh & mesh_
word GName() const
Helper function to return the name of the turbulence G field.
Abstract base class for turbulence models (RAS, LES and laminar).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const surfaceScalarField & alphaRhoPhi_
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
const nearWallDist & y() const
Return the near wall distances.
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest...
Definition: nearWallDist.H:51
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
turbulenceModel(const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const word &propertiesName)
Construct from components.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
const Time & time() const
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
Forward declarations of fvMatrix specializations.
label patchi
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
nearWallDist y_
Near wall distance boundary field.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool read()=0
Read model coefficients if they have changed.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
Namespace for OpenFOAM.
const volVectorField & U_