LienCubicKE.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-2022 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::LienCubicKE
26 
27 Description
28  Lien cubic non-linear low-Reynolds k-epsilon turbulence models for
29  incompressible flows.
30 
31  This turbulence model is described in:
32  \verbatim
33  Lien, F.S., Chen, W.L. & Leschziner, M.A. (1996).
34  Low-Reynolds-number eddy-viscosity modeling based on non-linear
35  stress-strain/vorticity relations.
36  Engineering Turbulence Modelling and Experiments 3, 91-100.
37  \endverbatim
38 
39  Implemented according to the specification in:
40  <a href=
41  "https://personalpages.manchester.ac.uk/staff/david.d.apsley/turbmod.pdf"
42  >Apsley: Turbulence Models 2002</a>
43 
44  In addition to the low-Reynolds number damping functions support for
45  wall-functions is also included to allow for low- and high-Reynolds number
46  operation.
47 
48 See also
49  Foam::incompressible::RASModels::ShihQuadraticKE
50 
51 SourceFiles
52  LienCubicKE.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef LienCubicKE_H
57 #define LienCubicKE_H
58 
60 #include "nonlinearEddyViscosity.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 namespace incompressible
67 {
68 namespace RASModels
69 {
70 
71 /*---------------------------------------------------------------------------*\
72  Class LienCubicKE Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class LienCubicKE
76 :
77  public nonlinearEddyViscosity<incompressible::RASModel>
78 {
79 
80 protected:
81 
82  // Protected data
83 
84  // Model coefficients
85 
99 
105 
106 
107  // Fields
111 
112  //- Wall distance
113  // Note: different to wall distance in parent RASModel
114  // which is for near-wall cells only
115  const volScalarField& y_;
116 
117 
118  // Protected Member Functions
119 
120  tmp<volScalarField> fMu() const;
121  tmp<volScalarField> f2() const;
122  tmp<volScalarField> E(const volScalarField& f2) const;
123 
124  virtual void correctNut();
125  virtual void correctNonlinearStress(const volTensorField& gradU);
126 
127 
128 public:
129 
130  //- Runtime type information
131  TypeName("LienCubicKE");
132 
133  // Constructors
134 
135  //- Construct from components
137  (
138  const geometricOneField& alpha,
139  const geometricOneField& rho,
140  const volVectorField& U,
141  const surfaceScalarField& alphaRhoPhi,
142  const surfaceScalarField& phi,
143  const viscosity& viscosity,
144  const word& type = typeName
145  );
146 
147 
148  //- Destructor
149  virtual ~LienCubicKE()
150  {}
151 
152 
153  // Member Functions
154 
155  //- Read RASProperties dictionary
156  virtual bool read();
157 
158  //- Return the effective diffusivity for k
159  tmp<volScalarField> DkEff() const
160  {
161  return volScalarField::New
162  (
163  "DkEff",
164  nut_/sigmak_ + nu()
165  );
166  }
167 
168  //- Return the effective diffusivity for epsilon
170  {
171  return volScalarField::New
172  (
173  "DepsilonEff",
174  nut_/sigmaEps_ + nu()
175  );
176  }
177 
178  //- Return the turbulence kinetic energy
179  virtual tmp<volScalarField> k() const
180  {
181  return k_;
182  }
183 
184  //- Return the turbulence kinetic energy dissipation rate
185  virtual tmp<volScalarField> epsilon() const
186  {
187  return epsilon_;
188  }
189 
190  //- Return the turbulence specific dissipation rate
191  virtual tmp<volScalarField> omega() const
192  {
193  return volScalarField::New
194  (
195  "omega",
196  epsilon_/(Cmu_*k_),
197  epsilon_.boundaryField().types()
198  );
199  }
200 
201  //- Solve the turbulence equations and correct the turbulence viscosity
202  virtual void correct();
203 };
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace RASModels
209 } // Edn namespace incompressible
210 } // End namespace Foam
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 #endif
215 
216 // ************************************************************************* //
const volScalarField & y_
Wall distance.
Definition: LienCubicKE.H:114
TypeName("LienCubicKE")
Runtime type information.
U
Definition: pEqn.H:72
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< volScalarField > fMu() const
Definition: LienCubicKE.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
wordList types() const
Return a list of the patch field types.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Eddy viscosity turbulence model with non-linear correction base class.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: LienCubicKE.H:178
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual void correctNonlinearStress(const volTensorField &gradU)
Definition: LienCubicKE.C:85
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read RASProperties dictionary.
Definition: LienCubicKE.C:338
phi
Definition: correctPhi.H:3
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
LienCubicKE(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.
Definition: LienCubicKE.C:123
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienCubicKE.H:168
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LienCubicKE.C:369
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
tmp< volScalarField > f2() const
Definition: LienCubicKE.C:57
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: LienCubicKE.H:190
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienCubicKE.H:158
tmp< volScalarField > E(const volScalarField &f2) const
Definition: LienCubicKE.C:65
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: LienCubicKE.H:184
Namespace for OpenFOAM.
Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows...
Definition: LienCubicKE.H:74