RASModel.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-2017 OpenFOAM Foundation
9 Copyright (C) 2021 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::RASModel
29
30Description
31 Templated abstract base class for RAS turbulence models
32
33SourceFiles
34 RASModel.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef RASModel_H
39#define RASModel_H
40
41#include "TurbulenceModel.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48/*---------------------------------------------------------------------------*\
49 Class RASModel Declaration
50\*---------------------------------------------------------------------------*/
51
52template<class BasicTurbulenceModel>
53class RASModel
54:
55 public BasicTurbulenceModel
56{
57
58protected:
59
60 // Protected data
61
62 //- RAS coefficients dictionary
64
65 //- Turbulence on/off flag
67
68 //- Flag to print the model coeffs at run-time
70
71 //- Model coefficients dictionary
73
74 //- Lower limit of k
76
77 //- Lower limit of epsilon
79
80 //- Lower limit for omega
82
83
84 // Protected Member Functions
85
86 //- Print model coefficients
87 virtual void printCoeffs(const word& type);
88
89 //- No copy construct
90 RASModel(const RASModel&) = delete;
91
92 //- No copy assignment
93 void operator=(const RASModel&) = delete;
94
95
96public:
98 typedef typename BasicTurbulenceModel::alphaField alphaField;
99 typedef typename BasicTurbulenceModel::rhoField rhoField;
100 typedef typename BasicTurbulenceModel::transportModel transportModel;
101
102
103 //- Runtime type information
104 TypeName("RAS");
105
106
107 // Declare run-time constructor selection table
110 (
111 autoPtr,
112 RASModel,
114 (
115 const alphaField& alpha,
116 const rhoField& rho,
117 const volVectorField& U,
118 const surfaceScalarField& alphaRhoPhi,
119 const surfaceScalarField& phi,
120 const transportModel& transport,
121 const word& propertiesName
122 ),
123 (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
124 );
125
126
127 // Constructors
128
129 //- Construct from components
131 (
132 const word& type,
133 const alphaField& alpha,
134 const rhoField& rho,
135 const volVectorField& U,
136 const surfaceScalarField& alphaRhoPhi,
137 const surfaceScalarField& phi,
138 const transportModel& transport,
139 const word& propertiesName
140 );
141
142
143 // Selectors
144
145 //- Return a reference to the selected RAS model
147 (
148 const alphaField& alpha,
149 const rhoField& rho,
150 const volVectorField& U,
151 const surfaceScalarField& alphaRhoPhi,
152 const surfaceScalarField& phi,
153 const transportModel& transport,
154 const word& propertiesName = turbulenceModel::propertiesName
155 );
156
157
158 //- Destructor
159 virtual ~RASModel() = default;
160
161
162 // Member Functions
163
164 //- Read model coefficients if they have changed
165 virtual bool read();
166
167
168 // Access
169
170 //- Return the lower allowable limit for k (default: SMALL)
171 const dimensionedScalar& kMin() const
172 {
173 return kMin_;
174 }
175
176 //- Return the lower allowable limit for epsilon (default: SMALL)
177 const dimensionedScalar& epsilonMin() const
178 {
179 return epsilonMin_;
180 }
181
182 //- Return the lower allowable limit for omega (default: SMALL)
183 const dimensionedScalar& omegaMin() const
184 {
185 return omegaMin_;
186 }
187
188 //- Allow kMin to be changed
190 {
191 return kMin_;
192 }
193
194 //- Allow epsilonMin to be changed
196 {
197 return epsilonMin_;
198 }
199
200 //- Allow omegaMin to be changed
202 {
203 return omegaMin_;
204 }
205
206 //- Const access to the coefficients dictionary
207 virtual const dictionary& coeffDict() const
208 {
209 return coeffDict_;
210 }
211
212
213 //- Return the effective viscosity
214 virtual tmp<volScalarField> nuEff() const
215 {
217 (
219 (
220 IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
221 this->nut() + this->nu()
222 )
223 );
224 }
225
226 //- Return the effective viscosity on patch
227 virtual tmp<scalarField> nuEff(const label patchi) const
228 {
229 return this->nut(patchi) + this->nu(patchi);
230 }
231
232 //- Return the turbulence kinetic energy dissipation rate
233 virtual tmp<volScalarField> epsilon() const;
234
235 //- Return the specific dissipation rate
236 virtual tmp<volScalarField> omega() const;
237
238 //- Solve the turbulence equations and correct the turbulence viscosity
239 virtual void correct();
240};
241
242
243// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244
245} // End namespace Foam
246
247// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248
249#ifdef NoRepository
250 #include "RASModel.C"
251#endif
252
253// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254
255#endif
256
257// ************************************************************************* //
surfaceScalarField & phi
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:55
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: RASModel.C:192
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:97
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:213
dictionary coeffDict_
Model coefficients dictionary.
Definition: RASModel.H:71
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:98
RASModel(const RASModel &)=delete
No copy construct.
virtual tmp< scalarField > nuEff(const label patchi) const
Return the effective viscosity on patch.
Definition: RASModel.H:226
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:230
dictionary RASDict_
RAS coefficients dictionary.
Definition: RASModel.H:62
const dimensionedScalar & epsilonMin() const
Return the lower allowable limit for epsilon (default: SMALL)
Definition: RASModel.H:176
virtual ~RASModel()=default
Destructor.
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: RASModel.H:80
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:206
dimensionedScalar & omegaMin()
Allow omegaMin to be changed.
Definition: RASModel.H:200
dimensionedScalar & epsilonMin()
Allow epsilonMin to be changed.
Definition: RASModel.H:194
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:77
static autoPtr< RASModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected RAS model.
Definition: RASModel.C:115
declareRunTimeSelectionTable(autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: SMALL)
Definition: RASModel.H:170
void operator=(const RASModel &)=delete
No copy assignment.
Switch turbulence_
Turbulence on/off flag.
Definition: RASModel.H:65
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:74
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:99
const dimensionedScalar & omegaMin() const
Return the lower allowable limit for omega (default: SMALL)
Definition: RASModel.H:182
TypeName("RAS")
Runtime type information.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: RASModel.C:211
dimensionedScalar & kMin()
Allow kMin to be changed.
Definition: RASModel.H:188
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: RASModel.H:68
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:171
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
scalar nut
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
volScalarField & nu
volScalarField & alpha
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73