LienLeschziner.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::incompressible::RASModels::LienLeschziner
26 
27 Description
28  Lien and Leschziner low-Reynolds number k-epsilon turbulence model for
29  incompressible flows.
30 
31  This turbulence model is described in:
32  \verbatim
33  Lien, F. S., & Leschziner, M. A. (1993).
34  A pressure-velocity solution strategy for compressible flow
35  and its application to shock/boundary-layer interaction
36  using second-moment turbulence closure.
37  Journal of fluids engineering, 115(4), 717-725.
38  \endverbatim
39 
40  Implemented according to the specification in:
41  <a href=
42  "https://personalpages.manchester.ac.uk/staff/david.d.apsley/turbmod.pdf"
43  >Apsley: Turbulence Models 2002</a>
44 
45  In addition to the low-Reynolds number damping functions support for
46  wall-functions is also included to allow for low- and high-Reynolds number
47  operation.
48 
49 SourceFiles
50  LienLeschziner.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef LienLeschziner_H
55 #define LienLeschziner_H
56 
58 #include "eddyViscosity.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace incompressible
65 {
66 namespace RASModels
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class LienLeschziner Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class LienLeschziner
74 :
75  public eddyViscosity<incompressible::RASModel>
76 {
77 
78 protected:
79 
80  // Protected data
81 
82  // Model coefficients
83 
88 
91 
95 
96 
97  // Fields
98 
101 
102 
103  // Protected Member Functions
104 
107  tmp<volScalarField> E(const volScalarField& f2) const;
108 
109  //- Bound epsilon and return Cmu*sqr(k) for nut
111 
112  //- Correct the eddy-viscosity nut
113  virtual void correctNut();
114 
115 
116 public:
117 
118  TypeName("LienLeschziner");
119 
120  // Constructors
121 
122  //- Construct from components
124  (
125  const geometricOneField& alpha,
126  const geometricOneField& rho,
127  const volVectorField& U,
128  const surfaceScalarField& alphaRhoPhi,
129  const surfaceScalarField& phi,
130  const viscosity& viscosity,
131  const word& type = typeName
132  );
133 
134 
135  //- Destructor
136  virtual ~LienLeschziner()
137  {}
138 
139 
140  // Member Functions
141 
142  //- Read RASProperties dictionary
143  virtual bool read();
144 
145  //- Return the effective diffusivity for k
146  tmp<volScalarField> DkEff() const
147  {
148  return volScalarField::New
149  (
150  "DkEff",
151  nut_/sigmak_ + nu()
152  );
153  }
154 
155  //- Return the effective diffusivity for epsilon
157  {
158  return volScalarField::New
159  (
160  "DepsilonEff",
161  nut_/sigmaEps_ + nu()
162  );
163  }
164 
165  //- Return the turbulence kinetic energy
166  virtual tmp<volScalarField> k() const
167  {
168  return k_;
169  }
170 
171  //- Return the turbulence kinetic energy dissipation rate
172  virtual tmp<volScalarField> epsilon() const
173  {
174  return epsilon_;
175  }
176 
177  //- Return the turbulence specific dissipation rate
178  virtual tmp<volScalarField> omega() const
179  {
180  return volScalarField::New
181  (
182  "omega",
183  epsilon_/(Cmu_*k_)
184  );
185  }
186 
187  //- Solve the turbulence equations and correct the turbulence viscosity
188  virtual void correct();
189 };
190 
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 } // End namespace RASModels
195 } // End namespace incompressible
196 } // End namespace Foam
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 #endif
201 
202 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
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 correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
tmp< volScalarField > E(const volScalarField &f2) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
LienLeschziner(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual bool read()
Read RASProperties dictionary.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488