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-2018 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 
63 protected:
64 
65  // Protected data
66 
67  const Time& runTime_;
68  const fvMesh& mesh_;
69 
72  const surfaceScalarField& phi_;
73 
74  //- Near wall distance boundary field
76 
77 
78 private:
79 
80  // Private Member Functions
81 
82  //- Disallow default bitwise copy construct
84 
85  //- Disallow default bitwise assignment
86  void operator=(const turbulenceModel&);
87 
88 
89 public:
90 
91  //- Runtime type information
92  TypeName("turbulenceModel");
93 
94  //- Default name of the turbulence properties dictionary
95  static const word propertiesName;
96 
97 
98  // Constructors
99 
100  //- Construct from components
102  (
103  const volVectorField& U,
105  const surfaceScalarField& phi,
106  const word& propertiesName
107  );
108 
109 
110  //- Destructor
111  virtual ~turbulenceModel()
112  {}
113 
114 
115  // Member Functions
116 
117  //- Read model coefficients if they have changed
118  virtual bool read() = 0;
120  const Time& time() const
121  {
122  return runTime_;
123  }
125  const fvMesh& mesh() const
126  {
127  return mesh_;
128  }
129 
130  //- Const access to the coefficients dictionary
131  virtual const dictionary& coeffDict() const = 0;
132 
133  //- Helper function to return the name of the turbulence G field
134  inline word GName() const
135  {
136  return word(type() + ":G");
137  }
138 
139  //- Access function to velocity field
140  inline const volVectorField& U() const
141  {
142  return U_;
143  }
144 
145  //- Access function to phase flux field
146  inline const surfaceScalarField& alphaRhoPhi() const
147  {
148  return alphaRhoPhi_;
149  }
150 
151  //- Return the volumetric flux field
152  virtual tmp<surfaceScalarField> phi() const;
153 
154  //- Return the near wall distances
155  const nearWallDist& y() const
156  {
157  return y_;
158  }
159 
160  //- Return the laminar viscosity
161  virtual tmp<volScalarField> nu() const = 0;
162 
163  //- Return the laminar viscosity on patch
164  virtual tmp<scalarField> nu(const label patchi) const = 0;
165 
166  //- Return the turbulence viscosity
167  virtual tmp<volScalarField> nut() const = 0;
168 
169  //- Return the turbulence viscosity on patch
170  virtual tmp<scalarField> nut(const label patchi) const = 0;
171 
172  //- Return the effective viscosity
173  virtual tmp<volScalarField> nuEff() const = 0;
174 
175  //- Return the effective viscosity on patch
176  virtual tmp<scalarField> nuEff(const label patchi) const = 0;
177 
178  //- Return the laminar dynamic viscosity
179  virtual tmp<volScalarField> mu() const = 0;
180 
181  //- Return the laminar dynamic viscosity on patch
182  virtual tmp<scalarField> mu(const label patchi) const = 0;
183 
184  //- Return the turbulence dynamic viscosity
185  virtual tmp<volScalarField> mut() const = 0;
186 
187  //- Return the turbulence dynamic viscosity on patch
188  virtual tmp<scalarField> mut(const label patchi) const = 0;
189 
190  //- Return the effective dynamic viscosity
191  virtual tmp<volScalarField> muEff() const = 0;
192 
193  //- Return the effective dynamic viscosity on patch
194  virtual tmp<scalarField> muEff(const label patchi) const = 0;
195 
196  //- Return the turbulence kinetic energy
197  virtual tmp<volScalarField> k() const = 0;
198 
199  //- Return the turbulence kinetic energy dissipation rate
200  virtual tmp<volScalarField> epsilon() const = 0;
201 
202  //- Return the Reynolds stress tensor
203  virtual tmp<volSymmTensorField> R() const = 0;
204 
205  //- Validate the turbulence fields after construction
206  // Update derived fields as required
207  virtual void validate();
208 
209  //- Solve the turbulence equations and correct the turbulence viscosity
210  virtual void correct() = 0;
211 };
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
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:137
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.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
const Time & time() const
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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
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_