ShihQuadraticKE.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-2015 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::ShihQuadraticKE
26 
27 Group
28  grpIcoRASTurbulence
29 
30 Description
31  Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for
32  incompressible flows
33 
34  This turbulence model is described in:
35  \verbatim
36  Shih, T. H., Zhu, J., & Lumley, J. L. (1993).
37  A realizable Reynolds stress algebraic equation model.
38  NASA technical memorandum 105993.
39  \endverbatim
40 
41  Implemented according to the specification in:
42  <a href=
43  "http://personalpages.manchester.ac.uk/staff/david.d.apsley/specturb.pdf"
44  >Apsley: Turbulence Models 2002</a>
45 
46 SourceFiles
47  ShihQuadraticKE.C
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #ifndef ShihQuadraticKE_H
52 #define ShihQuadraticKE_H
53 
55 #include "nonlinearEddyViscosity.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61 namespace incompressible
62 {
63 namespace RASModels
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class ShihQuadraticKE Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class ShihQuadraticKE
71 :
72  public nonlinearEddyViscosity<incompressible::RASModel>
73 {
74 
75 protected:
76 
77  // Protected data
78 
79  // Model coefficients
80 
91 
92 
93  // Fields
94 
97 
98 
99  // Protected Member Functions
100 
101  virtual void correctNut();
102  virtual void correctNonlinearStress(const volTensorField& gradU);
103 
104 
105 public:
106 
107  //- Runtime type information
108  TypeName("ShihQuadraticKE");
109 
110 
111  // Constructors
112 
113  //- Construct from components
115  (
116  const geometricOneField& alpha,
117  const geometricOneField& rho,
118  const volVectorField& U,
119  const surfaceScalarField& alphaRhoPhi,
120  const surfaceScalarField& phi,
121  const transportModel& transport,
122  const word& propertiesName = turbulenceModel::propertiesName,
123  const word& type = typeName
124  );
125 
126 
127  //- Destructor
128  virtual ~ShihQuadraticKE()
129  {}
130 
131 
132  // Member Functions
133 
134  //- Read RASProperties dictionary
135  virtual bool read();
136 
137  //- Return the effective diffusivity for k
138  tmp<volScalarField> DkEff() const
139  {
140  return tmp<volScalarField>
141  (
142  new volScalarField("DkEff", nut_/sigmak_ + nu())
143  );
144  }
145 
146  //- Return the effective diffusivity for epsilon
148  {
149  return tmp<volScalarField>
150  (
151  new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
152  );
153  }
154 
155  //- Return the turbulence kinetic energy
156  virtual tmp<volScalarField> k() const
157  {
158  return k_;
159  }
160 
161  //- Return the turbulence kinetic energy dissipation rate
162  virtual tmp<volScalarField> epsilon() const
163  {
164  return epsilon_;
165  }
166 
167  //- Solve the turbulence equations and correct the turbulence viscosity
168  virtual void correct();
169 };
170 
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 } // End namespace RASModels
175 } // End namespace incompressible
176 } // End namespace Foam
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 #endif
181 
182 // ************************************************************************* //
surfaceScalarField & phi
TypeName("ShihQuadraticKE")
Runtime type information.
U
Definition: pEqn.H:83
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Shih&#39;s quadratic algebraic Reynolds stress k-epsilon turbulence model for incompressible flows...
ShihQuadraticKE(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.
Eddy viscosity turbulence model with non-linear correction base class.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
incompressible::RASModel::transportModel transportModel
Definition: eddyViscosity.H:75
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
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.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read RASProperties dictionary.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
virtual void correctNonlinearStress(const volTensorField &gradU)
A class for managing temporary objects.
Definition: PtrList.H:54
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.