1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 Foam::compressible::RASModels
29 Namespace for compressible RAS turbulence models.
33 Foam::compressible::RASModel
36 Abstract base class for turbulence models for compressible and combusting
43 \*---------------------------------------------------------------------------*/
45 #ifndef compressibleRASModel_H
46 #define compressibleRASModel_H
48 #include "compressible/turbulenceModel/turbulenceModel.H"
49 #include "volFields.H"
50 #include "surfaceFields.H"
51 #include "nearWallDist.H"
54 #include "fvMatrices.H"
55 #include "basicThermo.H"
56 #include "IOdictionary.H"
60 #include "runTimeSelectionTables.H"
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 namespace compressible
69 /*---------------------------------------------------------------------------*\
70 Class RASModel Declaration
71 \*---------------------------------------------------------------------------*/
75 public turbulenceModel,
83 //- Turbulence on/off flag
86 //- Flag to print the model coeffs at run-time
89 //- Model coefficients dictionary
90 dictionary coeffDict_;
92 //- Value of y+ at the edge of the laminar sublayer
96 dimensionedScalar k0_;
98 //- Lower limit of epsilon
99 dimensionedScalar epsilon0_;
101 //- Small epsilon value used to avoid divide by zero
102 dimensionedScalar epsilonSmall_;
104 //- Lower limit for omega
105 dimensionedScalar omega0_;
107 //- Small omega value used to avoid divide by zero
108 dimensionedScalar omegaSmall_;
110 //- Near wall distance boundary field
114 // Protected member functions
116 //- Print model coefficients
117 virtual void printCoeffs();
122 // Private Member Functions
124 //- Disallow default bitwise copy construct
125 RASModel(const RASModel&);
127 //- Disallow default bitwise assignment
128 void operator=(const RASModel&);
133 //- Runtime type information
134 TypeName("RASModel");
137 // Declare run-time constructor selection table
139 declareRunTimeSelectionTable
145 const volScalarField& rho,
146 const volVectorField& U,
147 const surfaceScalarField& phi,
148 const basicThermo& thermoPhysicalModel
150 (rho, U, phi, thermoPhysicalModel)
156 //- Construct from components
160 const volScalarField& rho,
161 const volVectorField& U,
162 const surfaceScalarField& phi,
163 const basicThermo& thermoPhysicalModel
169 //- Return a reference to the selected turbulence model
170 static autoPtr<RASModel> New
172 const volScalarField& rho,
173 const volVectorField& U,
174 const surfaceScalarField& phi,
175 const basicThermo& thermoPhysicalModel
188 //- Return the value of k0 which k is not allowed to be less than
189 const dimensionedScalar& k0() const
194 //- Return the value of epsilon0 which epsilon is not allowed to be
196 const dimensionedScalar& epsilon0() const
201 //- Return the value of epsilonSmall which is added to epsilon when
203 const dimensionedScalar& epsilonSmall() const
205 return epsilonSmall_;
208 //- Return the value of omega0 which epsilon is not allowed to be
210 const dimensionedScalar& omega0() const
215 //- Return the value of omegaSmall which is added to epsilon when
217 const dimensionedScalar& omegaSmall() const
222 //- Allow k0 to be changed
223 dimensionedScalar& k0()
228 //- Allow epsilon0 to be changed
229 dimensionedScalar& epsilon0()
234 //- Allow epsilonSmall to be changed
235 dimensionedScalar& epsilonSmall()
237 return epsilonSmall_;
240 //- Allow omega0 to be changed
241 dimensionedScalar& omega0()
246 //- Allow omegaSmall to be changed
247 dimensionedScalar& omegaSmall()
252 //- Return the near wall distances
253 const nearWallDist& y() const
258 //- Calculate y+ at the edge of the laminar sublayer
259 scalar yPlusLam(const scalar kappa, const scalar E) const;
261 //- Const access to the coefficients dictionary
262 const dictionary& coeffDict() const
268 //- Return the turbulence viscosity
269 virtual tmp<volScalarField> mut() const = 0;
271 //- Return the effective viscosity
272 virtual tmp<volScalarField> muEff() const
274 return tmp<volScalarField>
276 new volScalarField("muEff", mut() + mu())
280 //- Return the effective turbulent thermal diffusivity
281 virtual tmp<volScalarField> alphaEff() const = 0;
283 //- Return the turbulence kinetic energy
284 virtual tmp<volScalarField> k() const = 0;
286 //- Return the turbulence kinetic energy dissipation rate
287 virtual tmp<volScalarField> epsilon() const = 0;
289 //- Return the Reynolds stress tensor
290 virtual tmp<volSymmTensorField> R() const = 0;
292 //- Return the effective stress tensor including the laminar stress
293 virtual tmp<volSymmTensorField> devRhoReff() const = 0;
295 //- Return the source term for the momentum equation
296 virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const = 0;
298 //- Return yPlus for the given patch
299 virtual tmp<scalarField> yPlus
305 //- Solve the turbulence equations and correct the turbulence viscosity
306 virtual void correct() = 0;
308 //- Read RASProperties dictionary
309 virtual bool read() = 0;
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 } // End namespace compressible
316 } // End namespace Foam
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 // ************************************************************************* //