LienLeschziner.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) 2011-2015 OpenFOAM Foundation
9 Copyright (C) 2019-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::incompressible::RASModels::LienLeschziner
29
30Group
31 grpIcoRASTurbulence
32
33Description
34 Lien and Leschziner low-Reynolds number k-epsilon turbulence model for
35 incompressible flows.
36
37 This turbulence model is described in:
38 \verbatim
39 Lien, F. S., & Leschziner, M. A. (1993).
40 A pressure-velocity solution strategy for compressible flow
41 and its application to shock/boundary-layer interaction
42 using second-moment turbulence closure.
43 Journal of fluids engineering, 115(4), 717-725.
44 \endverbatim
45
46 Implemented according to the specification in:
47 <a href=
48 "http://personalpages.manchester.ac.uk/staff/david.d.apsley/specturb.pdf"
49 >Apsley: Turbulence Models 2002</a>
50
51 In addition to the low-Reynolds number damping functions support for
52 wall-functions is also included to allow for low- and high-Reynolds number
53 operation.
54
55SourceFiles
56 LienLeschziner.C
57
58\*---------------------------------------------------------------------------*/
59
60#ifndef LienLeschziner_H
61#define LienLeschziner_H
62
64#include "eddyViscosity.H"
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68namespace Foam
69{
70namespace incompressible
71{
72namespace RASModels
73{
74
75/*---------------------------------------------------------------------------*\
76 Class LienLeschziner Declaration
77\*---------------------------------------------------------------------------*/
80:
81 public eddyViscosity<incompressible::RASModel>
82{
83
84protected:
85
86 // Protected data
87
88 // Model coefficients
101
102
103 // Fields
107
108 //- Wall distance
109 // Note: different to wall distance in parent RASModel
110 // which is for near-wall cells only
111 const volScalarField& y_;
112
113
114 // Protected Member Functions
115
116 tmp<volScalarField> fMu() const;
117 tmp<volScalarField> f2() const;
119
120 virtual void correctNut();
121
122
123public:
125 TypeName("LienLeschziner");
126
127 // Constructors
128
129 //- Construct from components
131 (
133 const geometricOneField& rho,
134 const volVectorField& U,
135 const surfaceScalarField& alphaRhoPhi,
136 const surfaceScalarField& phi,
137 const transportModel& transport,
138 const word& propertiesName = turbulenceModel::propertiesName,
139 const word& type = typeName
140 );
141
142
143 //- Destructor
144 virtual ~LienLeschziner() = default;
145
146
147 // Member Functions
148
149 //- Re-read model coefficients if they have changed
150 virtual bool read();
151
152 //- Return the effective diffusivity for k
154 {
156 (
157 new volScalarField("DkEff", nut_/sigmak_ + nu())
158 );
159 }
160
161 //- Return the effective diffusivity for epsilon
163 {
165 (
166 new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
167 );
168 }
169
170 //- Return the turbulence kinetic energy
171 virtual tmp<volScalarField> k() const
172 {
173 return k_;
174 }
175
176 //- Return the turbulence kinetic energy dissipation rate
177 virtual tmp<volScalarField> epsilon() const
178 {
179 return epsilon_;
180 }
181
182 //- Solve the turbulence equations and correct the turbulence viscosity
183 virtual void correct();
184};
185
186
187// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188
189} // End namespace RASModels
190} // End namespace incompressible
191} // End namespace Foam
192
193// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194
195#endif
196
197// ************************************************************************* //
surfaceScalarField & phi
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Lien and Leschziner low-Reynolds number k-epsilon turbulence model for incompressible flows.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
const volScalarField & y_
Wall distance.
tmp< volScalarField > E(const volScalarField &f2) const
virtual ~LienLeschziner()=default
Destructor.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
BasicTurbulenceModel::transportModel transportModel
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:290
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
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 TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73