ShihQuadraticKE.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-2019 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 Description
28  Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for
29  incompressible flows
30 
31  This turbulence model is described in:
32  \verbatim
33  Shih, T. H., Zhu, J., & Lumley, J. L. (1993).
34  A realizable Reynolds stress algebraic equation model.
35  NASA technical memorandum 105993.
36  \endverbatim
37 
38  Implemented according to the specification in:
39  <a href=
40  "https://personalpages.manchester.ac.uk/staff/david.d.apsley/turbmod.pdf"
41  >Apsley: Turbulence Models 2002</a>
42 
43 SourceFiles
44  ShihQuadraticKE.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef ShihQuadraticKE_H
49 #define ShihQuadraticKE_H
50 
52 #include "nonlinearEddyViscosity.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 namespace incompressible
59 {
60 namespace RASModels
61 {
62 
63 /*---------------------------------------------------------------------------*\
64  Class ShihQuadraticKE Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class ShihQuadraticKE
68 :
69  public nonlinearEddyViscosity<incompressible::RASModel>
70 {
71 
72 protected:
73 
74  // Protected data
75 
76  // Model coefficients
77 
88 
89 
90  // Fields
91 
94 
95 
96  // Protected Member Functions
97 
98  virtual void correctNut();
99  virtual void correctNonlinearStress(const volTensorField& gradU);
100 
101 
102 public:
103 
104  //- Runtime type information
105  TypeName("ShihQuadraticKE");
106 
107 
108  // Constructors
109 
110  //- Construct from components
112  (
113  const geometricOneField& alpha,
114  const geometricOneField& rho,
115  const volVectorField& U,
116  const surfaceScalarField& alphaRhoPhi,
117  const surfaceScalarField& phi,
118  const transportModel& transport,
119  const word& propertiesName = turbulenceModel::propertiesName,
120  const word& type = typeName
121  );
122 
123 
124  //- Destructor
125  virtual ~ShihQuadraticKE()
126  {}
127 
128 
129  // Member Functions
130 
131  //- Read RASProperties dictionary
132  virtual bool read();
133 
134  //- Return the effective diffusivity for k
135  tmp<volScalarField> DkEff() const
136  {
137  return volScalarField::New
138  (
139  "DkEff",
140  nut_/sigmak_ + nu()
141  );
142  }
143 
144  //- Return the effective diffusivity for epsilon
146  {
147  return volScalarField::New
148  (
149  "DepsilonEff",
150  nut_/sigmaEps_ + nu()
151  );
152  }
153 
154  //- Return the turbulence kinetic energy
155  virtual tmp<volScalarField> k() const
156  {
157  return k_;
158  }
159 
160  //- Return the turbulence kinetic energy dissipation rate
161  virtual tmp<volScalarField> epsilon() const
162  {
163  return epsilon_;
164  }
165 
166  //- Solve the turbulence equations and correct the turbulence viscosity
167  virtual void correct();
168 };
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace RASModels
174 } // End namespace incompressible
175 } // End namespace Foam
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 #endif
180 
181 // ************************************************************************* //
TypeName("ShihQuadraticKE")
Runtime type information.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
surfaceScalarField & phi
Shih&#39;s quadratic algebraic Reynolds stress k-epsilon turbulence model for incompressible flows...
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
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.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
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.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
U
Definition: pEqn.H:72
Base-class for all transport models used by the incompressible turbulence models. ...
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual void correctNonlinearStress(const volTensorField &gradU)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
volScalarField & nu
Namespace for OpenFOAM.