LienCubicKE.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 | Copyright (C) 2011-2016 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 Group
28  grpIcoRASTurbulence
29 
30 Description
31  Lien cubic non-linear low-Reynolds k-epsilon turbulence models for
32  incompressible flows.
33 
34  This turbulence model is described in:
35  \verbatim
36  Lien, F.S., Chen, W.L. & Leschziner, M.A. (1996).
37  Low-Reynolds-number eddy-viscosity modeling based on non-linear
38  stress-strain/vorticity relations.
39  Engineering Turbulence Modelling and Experiments 3, 91-100.
40  \endverbatim
41 
42  Implemented according to the specification in:
43  <a href=
44  "http://personalpages.manchester.ac.uk/staff/david.d.apsley/specturb.pdf"
45  >Apsley: Turbulence Models 2002</a>
46 
47  In addition to the low-Reynolds number damping functions support for
48  wall-functions is also included to allow for low- and high-Reynolds number
49  operation.
50 
51 See also
52  Foam::incompressible::RASModels::ShihQuadraticKE
53 
54 SourceFiles
55  LienCubicKE.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef LienCubicKE_H
60 #define LienCubicKE_H
61 
63 #include "nonlinearEddyViscosity.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 namespace incompressible
70 {
71 namespace RASModels
72 {
73 
74 /*---------------------------------------------------------------------------*\
75  Class LienCubicKE Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class LienCubicKE
79 :
80  public nonlinearEddyViscosity<incompressible::RASModel>
81 {
82 
83 protected:
84 
85  // Protected data
86 
87  // Model coefficients
88 
108 
109 
110  // Fields
114 
115  //- Wall distance
116  // Note: different to wall distance in parent RASModel
117  // which is for near-wall cells only
118  const volScalarField& y_;
119 
120 
121  // Protected Member Functions
122 
123  tmp<volScalarField> fMu() const;
124  tmp<volScalarField> f2() const;
125  tmp<volScalarField> E(const volScalarField& f2) const;
126 
127  virtual void correctNut();
128  virtual void correctNonlinearStress(const volTensorField& gradU);
129 
130 
131 public:
132 
133  //- Runtime type information
134  TypeName("LienCubicKE");
135 
136  // Constructors
137 
138  //- Construct from components
140  (
141  const geometricOneField& alpha,
142  const geometricOneField& rho,
143  const volVectorField& U,
144  const surfaceScalarField& alphaRhoPhi,
145  const surfaceScalarField& phi,
146  const transportModel& transport,
147  const word& propertiesName = turbulenceModel::propertiesName,
148  const word& type = typeName
149  );
150 
151 
152  //- Destructor
153  virtual ~LienCubicKE()
154  {}
155 
156 
157  // Member Functions
158 
159  //- Read RASProperties dictionary
160  virtual bool read();
161 
162  //- Return the effective diffusivity for k
163  tmp<volScalarField> DkEff() const
164  {
165  return tmp<volScalarField>
166  (
167  new volScalarField("DkEff", nut_/sigmak_ + nu())
168  );
169  }
170 
171  //- Return the effective diffusivity for epsilon
173  {
174  return tmp<volScalarField>
175  (
176  new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
177  );
178  }
179 
180  //- Return the turbulence kinetic energy
181  virtual tmp<volScalarField> k() const
182  {
183  return k_;
184  }
185 
186  //- Return the turbulence kinetic energy dissipation rate
187  virtual tmp<volScalarField> epsilon() const
188  {
189  return epsilon_;
190  }
191 
192  //- Solve the turbulence equations and correct the turbulence viscosity
193  virtual void correct();
194 };
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace RASModels
200 } // Edn namespace incompressible
201 } // End namespace Foam
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 #endif
206 
207 // ************************************************************************* //
const volScalarField & y_
Wall distance.
Definition: LienCubicKE.H:117
surfaceScalarField & phi
U
Definition: pEqn.H:83
TypeName("LienCubicKE")
Runtime type information.
tmp< volScalarField > fMu() const
Definition: LienCubicKE.C:47
Eddy viscosity turbulence model with non-linear correction base class.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienCubicKE.H:162
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: LienCubicKE.H:186
incompressible::RASModel::transportModel transportModel
Definition: eddyViscosity.H:75
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual void correctNonlinearStress(const volTensorField &gradU)
Definition: LienCubicKE.C:85
LienCubicKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: LienCubicKE.C:123
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read RASProperties dictionary.
Definition: LienCubicKE.C:340
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienCubicKE.H:171
tmp< volScalarField > f2() const
Definition: LienCubicKE.C:57
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: LienCubicKE.H:180
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LienCubicKE.C:371
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
tmp< volScalarField > E(const volScalarField &f2) const
Definition: LienCubicKE.C:65
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.
Lien cubic non-linear low-Reynolds k-epsilon turbulence models for incompressible flows...
Definition: LienCubicKE.H:77