Loading...
Searching...
No Matches
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-2022 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 Foam_RASModel_H
39#define Foam_RASModel_H
40
41#include "TurbulenceModel.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
48/*---------------------------------------------------------------------------*\
49 Class RASModelBase Declaration
50\*---------------------------------------------------------------------------*/
51
52class RASModelBase
53{
54public:
55
56 //- Runtime type information
57 ClassName("RASModelBase");
58
59
60 // Constructors
61
62 //- Constructor
63 RASModelBase() noexcept = default;
64
65
66 //- Destructor
67 virtual ~RASModelBase() = default;
68};
69
70
71/*---------------------------------------------------------------------------*\
72 Class RASModel Declaration
73\*---------------------------------------------------------------------------*/
74
75template<class BasicTurbulenceModel>
76class RASModel
78 public RASModelBase,
79 public BasicTurbulenceModel
80{
81protected:
82
83 // Protected Data
84
85 //- RAS coefficients dictionary
87
88 //- Turbulence on/off flag
90
91 //- Flag to print the model coeffs at run-time
93
94 //- Model coefficients dictionary
96
97 //- Lower limit of k
100 //- Lower limit of epsilon
102
103 //- Lower limit for omega
105
106
107 // Protected Member Functions
108
109 //- Print model coefficients
110 virtual void printCoeffs(const word& type);
111
112 //- No copy construct
113 RASModel(const RASModel&) = delete;
115 //- No copy assignment
116 void operator=(const RASModel&) = delete;
117
118
119public:
120
121 typedef typename BasicTurbulenceModel::alphaField alphaField;
122 typedef typename BasicTurbulenceModel::rhoField rhoField;
123 typedef typename BasicTurbulenceModel::transportModel transportModel;
124
125
126 //- Runtime type information
127 TypeName("RAS");
128
129
130 // Declare run-time constructor selection table
131
133 (
134 autoPtr,
135 RASModel,
138 const alphaField& alpha,
139 const rhoField& rho,
140 const volVectorField& U,
141 const surfaceScalarField& alphaRhoPhi,
143 const transportModel& transport,
144 const word& propertiesName
145 ),
146 (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
147 );
148
149
150 // Constructors
151
152 //- Construct from components
154 (
155 const word& type,
156 const alphaField& alpha,
157 const rhoField& rho,
158 const volVectorField& U,
159 const surfaceScalarField& alphaRhoPhi,
160 const surfaceScalarField& phi,
161 const transportModel& transport,
162 const word& propertiesName
163 );
164
165
166 // Selectors
167
168 //- Return a reference to the selected RAS model
170 (
171 const alphaField& alpha,
172 const rhoField& rho,
173 const volVectorField& U,
174 const surfaceScalarField& alphaRhoPhi,
175 const surfaceScalarField& phi,
176 const transportModel& transport,
177 const word& propertiesName = turbulenceModel::propertiesName
178 );
179
180
181 //- Destructor
182 virtual ~RASModel() = default;
183
184
185 // Member Functions
186
187 //- Read model coefficients if they have changed
188 virtual bool read();
189
190
191 // Access
192
193 //- Return the lower allowable limit for k (default: SMALL)
194 const dimensionedScalar& kMin() const noexcept
195 {
196 return kMin_;
197 }
198
199 //- Return the lower allowable limit for epsilon (default: SMALL)
200 const dimensionedScalar& epsilonMin() const noexcept
201 {
202 return epsilonMin_;
203 }
204
205 //- Return the lower allowable limit for omega (default: SMALL)
206 const dimensionedScalar& omegaMin() const noexcept
207 {
208 return omegaMin_;
209 }
210
211 //- Allow kMin to be changed
213 {
214 return kMin_;
215 }
216
217 //- Allow epsilonMin to be changed
218 dimensionedScalar& epsilonMin() noexcept
220 return epsilonMin_;
221 }
222
223 //- Allow omegaMin to be changed
224 dimensionedScalar& omegaMin() noexcept
225 {
226 return omegaMin_;
228
229 //- Const access to the coefficients dictionary
230 virtual const dictionary& coeffDict() const
231 {
232 return coeffDict_;
233 }
234
236 //- Return the effective viscosity
237 virtual tmp<volScalarField> nuEff() const
238 {
240 (
242 (
243 IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
244 this->nut() + this->nu()
245 )
246 );
247 }
248
249 //- Return the effective viscosity on patch
250 virtual tmp<scalarField> nuEff(const label patchi) const
252 return this->nut(patchi) + this->nu(patchi);
253 }
254
255 //- Return the turbulence kinetic energy dissipation rate
256 virtual tmp<volScalarField> epsilon() const;
257
258 //- Return the specific dissipation rate
259 virtual tmp<volScalarField> omega() const;
260
261 //- Solve the turbulence equations and correct the turbulence viscosity
262 virtual void correct();
263};
264
265
266// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268} // End namespace Foam
269
270// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271
272#ifdef NoRepository
273 #include "RASModel.C"
274#endif
276// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277
278#endif
279
280// ************************************************************************* //
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
ClassName("RASModelBase")
Runtime type information.
RASModelBase() noexcept=default
Constructor.
Templated abstract base class for RAS turbulence models.
Definition RASModel.H:81
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition RASModel.C:185
EddyDiffusivity< fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
Definition RASModel.H:142
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition RASModel.H:284
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition RASModel.C:27
dimensionedScalar & omegaMin() noexcept
Allow omegaMin to be changed.
Definition RASModel.H:267
EddyDiffusivity< fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
Definition RASModel.H:143
RASModel(const RASModel &)=delete
No copy construct.
virtual tmp< scalarField > nuEff(const label patchi) const
Return the effective viscosity on patch.
Definition RASModel.H:299
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition RASModel.C:215
virtual ~RASModel()=default
Destructor.
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition RASModel.H:275
dimensionedScalar & epsilonMin() noexcept
Allow epsilonMin to be changed.
Definition RASModel.H:259
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:108
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 & omegaMin() const noexcept
Return the lower allowable limit for omega (default: SMALL).
Definition RASModel.H:243
void operator=(const RASModel &)=delete
No copy assignment.
EddyDiffusivity< fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
Definition RASModel.H:144
const dimensionedScalar & kMin() const noexcept
Return the lower allowable limit for k (default: SMALL).
Definition RASModel.H:227
RASModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
Definition RASModel.C:40
TypeName("RAS")
Runtime type information.
const dimensionedScalar & epsilonMin() const noexcept
Return the lower allowable limit for epsilon (default: SMALL).
Definition RASModel.H:235
dimensionedScalar & kMin() noexcept
Allow kMin to be changed.
Definition RASModel.H:251
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition RASModel.C:200
virtual bool read()
Read model coefficients if they have changed.
Definition RASModel.C:164
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
U
Definition pEqn.H:72
thermo correct()
scalar epsilon
scalar nut
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const direction noexcept
Definition scalarImpl.H:265
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:68