turbulenceModel.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2015 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::turbulenceModel
29
30Group
31 grpTurbulence
32
33Description
34 Abstract base class for turbulence models (RAS, LES and laminar).
35
36SourceFiles
37 turbulenceModel.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef turbulenceModel_H
42#define turbulenceModel_H
43
44#include "IOdictionary.H"
45#include "primitiveFieldsFwd.H"
46#include "volFieldsFwd.H"
47#include "surfaceFieldsFwd.H"
48#include "fvMatricesFwd.H"
49#include "nearWallDist.H"
50#include "geometricOneField.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57// Forward declarations
58class fvMesh;
59
60/*---------------------------------------------------------------------------*\
61 Class turbulenceModel Declaration
62\*---------------------------------------------------------------------------*/
65:
66 public IOdictionary
67{
68
69protected:
70
71 // Protected data
73 const Time& runTime_;
74 const fvMesh& mesh_;
79
80 //- Near wall distance boundary field
82
83
84private:
85
86 // Private Member Functions
87
88 //- No copy construct
89 turbulenceModel(const turbulenceModel&) = delete;
90
91 //- No copy assignment
92 void operator=(const turbulenceModel&) = delete;
93
94
95public:
96
97 //- Runtime type information
98 TypeName("turbulenceModel");
99
100 //- Default name of the turbulence properties dictionary
101 static const word propertiesName;
102
103
104 // Constructors
105
106 //- Construct from components
108 (
109 const volVectorField& U,
111 const surfaceScalarField& phi,
112 const word& propertiesName
113 );
114
115
116 //- Destructor
117 virtual ~turbulenceModel() = default;
118
119
120 // Member Functions
121
122 //- Read model coefficients if they have changed
123 virtual bool read() = 0;
125 const Time& time() const
126 {
127 return runTime_;
128 }
130 const fvMesh& mesh() const
131 {
132 return mesh_;
133 }
134
135 //- Const access to the coefficients dictionary
136 virtual const dictionary& coeffDict() const = 0;
137
138 //- Helper function to return the name of the turbulence G field
139 inline word GName() const
140 {
141 return word(type() + ":G");
142 }
143
144 //- Access function to velocity field
145 inline const volVectorField& U() const
146 {
147 return U_;
148 }
149
150 //- Access function to phase flux field
151 inline const surfaceScalarField& alphaRhoPhi() const
152 {
153 return alphaRhoPhi_;
154 }
155
156 //- Return the volumetric flux field
157 virtual tmp<surfaceScalarField> phi() const;
158
159 //- Return the near wall distances
160 const nearWallDist& y() const
161 {
162 return y_;
163 }
164
165 //- Return the laminar viscosity
166 virtual tmp<volScalarField> nu() const = 0;
167
168 //- Return the laminar viscosity on patch
169 virtual tmp<scalarField> nu(const label patchi) const = 0;
170
171 //- Return the turbulence viscosity
172 virtual tmp<volScalarField> nut() const = 0;
173
174 //- Return the turbulence viscosity on patch
175 virtual tmp<scalarField> nut(const label patchi) const = 0;
176
177 //- Return the effective viscosity
178 virtual tmp<volScalarField> nuEff() const = 0;
179
180 //- Return the effective viscosity on patch
181 virtual tmp<scalarField> nuEff(const label patchi) const = 0;
182
183 //- Return the laminar dynamic viscosity
184 virtual tmp<volScalarField> mu() const = 0;
185
186 //- Return the laminar dynamic viscosity on patch
187 virtual tmp<scalarField> mu(const label patchi) const = 0;
188
189 //- Return the turbulence dynamic viscosity
190 virtual tmp<volScalarField> mut() const = 0;
191
192 //- Return the turbulence dynamic viscosity on patch
193 virtual tmp<scalarField> mut(const label patchi) const = 0;
194
195 //- Return the effective dynamic viscosity
196 virtual tmp<volScalarField> muEff() const = 0;
197
198 //- Return the effective dynamic viscosity on patch
199 virtual tmp<scalarField> muEff(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 specific dissipation rate
208 virtual tmp<volScalarField> omega() const = 0;
209
210 //- Return the Reynolds stress tensor
211 virtual tmp<volSymmTensorField> R() const = 0;
212
213 //- Validate the turbulence fields after construction
214 // Update derived fields as required
215 virtual void validate();
216
217 //- Solve the turbulence equations and correct the turbulence viscosity
218 virtual void correct() = 0;
219};
220
221
222// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223
224} // End namespace Foam
225
226// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227
228#endif
229
230// ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:56
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
virtual tmp< scalarField > nut(const label patchi) const =0
Return the turbulence viscosity on patch.
virtual ~turbulenceModel()=default
Destructor.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > mut() const =0
Return the turbulence dynamic viscosity.
const fvMesh & mesh_
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
const Time & time() const
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< scalarField > mut(const label patchi) const =0
Return the turbulence dynamic viscosity on patch.
nearWallDist y_
Near wall distance boundary field.
virtual void validate()
Validate the turbulence fields after construction.
word GName() const
Helper function to return the name of the turbulence G field.
const surfaceScalarField & phi_
virtual tmp< scalarField > nu(const label patchi) const =0
Return the laminar viscosity on patch.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
virtual bool read()=0
Read model coefficients if they have changed.
virtual tmp< scalarField > muEff(const label patchi) const =0
Return the effective dynamic viscosity on patch.
const nearWallDist & y() const
Return the near wall distances.
const volVectorField & U_
virtual tmp< volScalarField > omega() const =0
Return the specific dissipation rate.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
virtual tmp< scalarField > nuEff(const label patchi) const =0
Return the effective viscosity on patch.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
const fvMesh & mesh() const
const surfaceScalarField & alphaRhoPhi_
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
virtual tmp< volSymmTensorField > R() const =0
Return the Reynolds stress tensor.
virtual tmp< volScalarField > muEff() const =0
Return the effective dynamic viscosity.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
TypeName("turbulenceModel")
Runtime type information.
virtual tmp< scalarField > mu(const label patchi) const =0
Return the laminar dynamic viscosity on patch.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
A class for handling words, derived from Foam::string.
Definition: word.H:68
Forward declarations of fvMatrix specializations.
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73