kOmega.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::RASModels::kOmega
26 
27 Description
28  Standard high Reynolds-number k-omega turbulence model for
29  incompressible and compressible flows.
30 
31  References:
32  \verbatim
33  Wilcox, D. C. (1998).
34  Turbulence modeling for CFD
35  (Vol. 2, pp. 103-217). La Canada, CA: DCW industries.
36  \endverbatim
37 
38  The default model coefficients are
39  \verbatim
40  kOmegaCoeffs
41  {
42  Cmu 0.09; // Equivalent to betaStar
43  alpha 0.52;
44  beta 0.072;
45  alphak 0.5;
46  alphaOmega 0.5;
47  }
48  \endverbatim
49 
50 SourceFiles
51  kOmega.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef kOmega_H
56 #define kOmega_H
57 
58 #include "RASModel.H"
59 #include "eddyViscosity.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 namespace RASModels
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class kOmega Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 template<class BasicTurbulenceModel>
73 class kOmega
74 :
75  public eddyViscosity<RASModel<BasicTurbulenceModel>>
76 {
77 
78 protected:
79 
80  // Protected data
81 
82  // Model coefficients
83 
89 
90 
91  // Fields
92 
95 
96 
97  // Protected Member Functions
98 
99  virtual void correctNut();
100  virtual tmp<fvScalarMatrix> kSource() const;
101  virtual tmp<fvScalarMatrix> omegaSource() const;
102 
103 
104 public:
106  typedef typename BasicTurbulenceModel::alphaField alphaField;
107  typedef typename BasicTurbulenceModel::rhoField rhoField;
108  typedef typename BasicTurbulenceModel::transportModel transportModel;
109 
110 
111  //- Runtime type information
112  TypeName("kOmega");
113 
114 
115  // Constructors
116 
117  //- Construct from components
118  kOmega
119  (
120  const alphaField& alpha,
121  const rhoField& rho,
122  const volVectorField& U,
123  const surfaceScalarField& alphaRhoPhi,
124  const surfaceScalarField& phi,
125  const transportModel& transport,
126  const word& propertiesName = turbulenceModel::propertiesName,
127  const word& type = typeName
128  );
129 
130 
131  //- Destructor
132  virtual ~kOmega()
133  {}
134 
135 
136  // Member Functions
137 
138  //- Read RASProperties dictionary
139  virtual bool read();
140 
141  //- Return the effective diffusivity for k
142  tmp<volScalarField> DkEff() const
143  {
144  return volScalarField::New
145  (
146  "DkEff",
147  alphaK_*this->nut_ + this->nu()
148  );
149  }
150 
151  //- Return the effective diffusivity for omega
153  {
154  return volScalarField::New
155  (
156  "DomegaEff",
157  alphaOmega_*this->nut_ + this->nu()
158  );
159  }
160 
161  //- Return the turbulence kinetic energy
162  virtual tmp<volScalarField> k() const
163  {
164  return k_;
165  }
166 
167  //- Return the turbulence specific dissipation rate
168  virtual tmp<volScalarField> omega() const
169  {
170  return omega_;
171  }
172 
173  //- Return the turbulence kinetic energy dissipation rate
174  virtual tmp<volScalarField> epsilon() const
175  {
176  return volScalarField::New
177  (
178  "epsilon",
179  Cmu_*k_*omega_,
180  omega_.boundaryField().types()
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 Foam
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 #ifdef NoRepository
196  #include "kOmega.C"
197 #endif
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 #endif
202 
203 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
Definition: kOmega.H:107
surfaceScalarField & phi
kOmega(const alphaField &alpha, const rhoField &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: kOmega.C:83
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kOmega.H:161
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: kOmega.H:173
dimensionedScalar alphaK_
Definition: kOmega.H:86
TypeName("kOmega")
Runtime type information.
BasicTurbulenceModel::alphaField alphaField
Definition: kOmega.H:105
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
dimensionedScalar Cmu_
Definition: kOmega.H:83
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kOmega.H:141
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.
static const word propertiesName
Default name of the turbulence properties dictionary.
BasicTurbulenceModel::rhoField rhoField
Definition: kOmega.H:106
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > DomegaEff() const
Return the effective diffusivity for omega.
Definition: kOmega.H:151
dimensionedScalar alphaOmega_
Definition: kOmega.H:87
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:210
volScalarField k_
Definition: kOmega.H:92
dimensionedScalar beta_
Definition: kOmega.H:84
virtual void correctNut()
Definition: kOmega.C:40
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:190
virtual ~kOmega()
Destructor.
Definition: kOmega.H:131
U
Definition: pEqn.H:72
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmega.C:51
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: kOmega.H:167
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 omega_
Definition: kOmega.H:93
Standard high Reynolds-number k-omega turbulence model for incompressible and compressible flows...
Definition: kOmega.H:72
volScalarField & nu
dimensionedScalar gamma_
Definition: kOmega.H:85
Namespace for OpenFOAM.
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmega.C:66