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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::incompressible::RASModels::ShihQuadraticKE
29 
30 Group
31  grpIcoRASTurbulence
32 
33 Description
34  Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for
35  incompressible flows
36 
37  This turbulence model is described in:
38  \verbatim
39  Shih, T. H., Zhu, J., & Lumley, J. L. (1993).
40  A realizable Reynolds stress algebraic equation model.
41  NASA technical memorandum 105993.
42  \endverbatim
43 
44  Implemented according to the specification in:
45  <a href=
46  "http://personalpages.manchester.ac.uk/staff/david.d.apsley/specturb.pdf"
47  >Apsley: Turbulence Models 2002</a>
48 
49 SourceFiles
50  ShihQuadraticKE.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef ShihQuadraticKE_H
55 #define ShihQuadraticKE_H
56 
58 #include "nonlinearEddyViscosity.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace incompressible
65 {
66 namespace RASModels
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class ShihQuadraticKE Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class ShihQuadraticKE
74 :
75  public nonlinearEddyViscosity<incompressible::RASModel>
76 {
77 
78 protected:
79 
80  // Protected data
81 
82  // Model coefficients
83 
94 
95 
96  // Fields
97 
100 
101 
102  // Protected Member Functions
103 
104  virtual void correctNut();
105  virtual void correctNonlinearStress(const volTensorField& gradU);
106 
107 
108 public:
109 
110  //- Runtime type information
111  TypeName("ShihQuadraticKE");
112 
113 
114  // Constructors
115 
116  //- Construct from components
118  (
119  const geometricOneField& alpha,
120  const geometricOneField& rho,
121  const volVectorField& U,
122  const surfaceScalarField& alphaRhoPhi,
123  const surfaceScalarField& phi,
124  const transportModel& transport,
125  const word& propertiesName = turbulenceModel::propertiesName,
126  const word& type = typeName
127  );
128 
129 
130  //- Destructor
131  virtual ~ShihQuadraticKE() = default;
132 
133 
134  // Member Functions
135 
136  //- Re-read model coefficients if they have changed
137  virtual bool read();
138 
139  //- Return the effective diffusivity for k
140  tmp<volScalarField> DkEff() const
141  {
142  return tmp<volScalarField>
143  (
144  new volScalarField("DkEff", nut_/sigmak_ + nu())
145  );
146  }
147 
148  //- Return the effective diffusivity for epsilon
150  {
151  return tmp<volScalarField>
152  (
153  new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu())
154  );
155  }
157  //- Return the turbulence kinetic energy
158  virtual tmp<volScalarField> k() const
159  {
160  return k_;
161  }
162 
163  //- Return the turbulence kinetic energy dissipation rate
164  virtual tmp<volScalarField> epsilon() const
165  {
166  return epsilon_;
167  }
168 
169  //- Solve the turbulence equations and correct the turbulence viscosity
170  virtual void correct();
171 };
172 
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 } // End namespace RASModels
177 } // End namespace incompressible
178 } // End namespace Foam
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 #endif
183 
184 // ************************************************************************* //
TypeName("ShihQuadraticKE")
Runtime type information.
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...
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.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:74
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
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 Foam::string.
Definition: word.H:63
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
U
Definition: pEqn.H:72
virtual ~ShihQuadraticKE()=default
Destructor.
virtual void correctNonlinearStress(const volTensorField &gradU)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.