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-2020 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  }
156 
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  //- Return the (estimated) specific dissipation rate
170  virtual tmp<volScalarField> omega() const
171  {
173  (
174  IOobject
175  (
176  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
177  this->runTime_.timeName(),
178  this->mesh_
179  ),
180  epsilon_/(0.09*k_)
181  );
182  }
183 
184  //- Solve the turbulence equations and correct the turbulence viscosity
185  virtual void correct();
186 };
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace RASModels
192 } // End namespace incompressible
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
Foam::incompressible::RASModels::ShihQuadraticKE::k_
volScalarField k_
Definition: ShihQuadraticKE.H:97
Foam::incompressible::RASModels::ShihQuadraticKE::sigmaEps_
dimensionedScalar sigmaEps_
Definition: ShihQuadraticKE.H:86
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::incompressible::RASModels::ShihQuadraticKE::sigmak_
dimensionedScalar sigmak_
Definition: ShihQuadraticKE.H:85
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::incompressible::RASModels::ShihQuadraticKE::Cmu1_
dimensionedScalar Cmu1_
Definition: ShihQuadraticKE.H:87
turbulentTransportModel.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
Foam::incompressible::RASModels::ShihQuadraticKE::ShihQuadraticKE
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.
Definition: ShihQuadraticKE.C:83
Foam::incompressible::RASModels::ShihQuadraticKE::~ShihQuadraticKE
virtual ~ShihQuadraticKE()=default
Destructor.
rho
rho
Definition: readInitialConditions.H:88
Foam::incompressible::RASModels::ShihQuadraticKE::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: ShihQuadraticKE.H:163
Foam::incompressible::RASModels::ShihQuadraticKE::correctNonlinearStress
virtual void correctNonlinearStress(const volTensorField &gradU)
Definition: ShihQuadraticKE.C:57
Foam::incompressible::RASModels::ShihQuadraticKE::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: ShihQuadraticKE.C:235
Foam::nonlinearEddyViscosity
Eddy viscosity turbulence model with non-linear correction base class.
Definition: nonlinearEddyViscosity.H:55
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::incompressible::RASModels::ShihQuadraticKE::Cmu2_
dimensionedScalar Cmu2_
Definition: ShihQuadraticKE.H:88
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta_
dimensionedScalar Cbeta_
Definition: ShihQuadraticKE.H:89
Foam::incompressible::RASModels::ShihQuadraticKE::TypeName
TypeName("ShihQuadraticKE")
Runtime type information.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::incompressible::RASModels::ShihQuadraticKE
Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for incompressible flows.
Definition: ShihQuadraticKE.H:72
Foam::dimensioned< scalar >
Foam::incompressible::RASModels::ShihQuadraticKE::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: ShihQuadraticKE.C:257
Foam::incompressible::RASModels::ShihQuadraticKE::Ceps2_
dimensionedScalar Ceps2_
Definition: ShihQuadraticKE.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressible::RASModels::ShihQuadraticKE::epsilon_
volScalarField epsilon_
Definition: ShihQuadraticKE.H:98
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:53
Foam::incompressible::RASModels::ShihQuadraticKE::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: ShihQuadraticKE.H:157
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta1_
dimensionedScalar Cbeta1_
Definition: ShihQuadraticKE.H:90
Foam::incompressible::RASModels::ShihQuadraticKE::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: ShihQuadraticKE.H:148
U
U
Definition: pEqn.H:72
Foam::incompressible::RASModels::ShihQuadraticKE::DkEff
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: ShihQuadraticKE.H:139
nonlinearEddyViscosity.H
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::incompressible::RASModels::ShihQuadraticKE::Ceps1_
dimensionedScalar Ceps1_
Definition: ShihQuadraticKE.H:83
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta3_
dimensionedScalar Cbeta3_
Definition: ShihQuadraticKE.H:92
Foam::eddyViscosity< incompressible::RASModel >::nut_
volScalarField nut_
Definition: eddyViscosity.H:66
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::incompressible::RASModels::ShihQuadraticKE::Cbeta2_
dimensionedScalar Cbeta2_
Definition: ShihQuadraticKE.H:91
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::incompressible::RASModels::ShihQuadraticKE::correctNut
virtual void correctNut()
Definition: ShihQuadraticKE.C:51
Foam::incompressible::RASModels::ShihQuadraticKE::omega
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
Definition: ShihQuadraticKE.H:169