All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2021 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  betaStar 0.09;
43  gamma 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 BasicMomentumTransportModel>
73 class kOmega
74 :
75  public eddyViscosity<RASModel<BasicMomentumTransportModel>>
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 BasicMomentumTransportModel::alphaField alphaField;
107  typedef typename BasicMomentumTransportModel::rhoField rhoField;
108 
109 
110  //- Runtime type information
111  TypeName("kOmega");
112 
113 
114  // Constructors
115 
116  //- Construct from components
117  kOmega
118  (
119  const alphaField& alpha,
120  const rhoField& rho,
121  const volVectorField& U,
122  const surfaceScalarField& alphaRhoPhi,
123  const surfaceScalarField& phi,
124  const viscosity& viscosity,
125  const word& type = typeName
126  );
127 
128 
129  //- Destructor
130  virtual ~kOmega()
131  {}
132 
133 
134  // Member Functions
135 
136  //- Read RASProperties dictionary
137  virtual bool read();
138 
139  //- Return the effective diffusivity for k
140  tmp<volScalarField> DkEff() const
141  {
142  return volScalarField::New
143  (
144  "DkEff",
145  alphaK_*this->nut_ + this->nu()
146  );
147  }
148 
149  //- Return the effective diffusivity for omega
151  {
152  return volScalarField::New
153  (
154  "DomegaEff",
155  alphaOmega_*this->nut_ + this->nu()
156  );
157  }
158 
159  //- Return the turbulence kinetic energy
160  virtual tmp<volScalarField> k() const
161  {
162  return k_;
163  }
164 
165  //- Return the turbulence specific dissipation rate
166  virtual tmp<volScalarField> omega() const
167  {
168  return omega_;
169  }
170 
171  //- Return the turbulence kinetic energy dissipation rate
172  virtual tmp<volScalarField> epsilon() const
173  {
174  return volScalarField::New
175  (
176  "epsilon",
177  betaStar_*k_*omega_,
178  omega_.boundaryField().types()
179  );
180  }
181 
182  //- Solve the turbulence equations and correct the turbulence viscosity
183  virtual void correct();
184 };
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace RASModels
190 } // End namespace Foam
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 #ifdef NoRepository
194  #include "kOmega.C"
195 #endif
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 #endif
200 
201 // ************************************************************************* //
U
Definition: pEqn.H:72
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kOmega.H:159
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:171
dimensionedScalar alphaK_
Definition: kOmega.H:86
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
TypeName("kOmega")
Runtime type information.
wordList types() const
Return a list of the patch field types.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmega.H:106
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kOmega.H:139
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:149
dimensionedScalar alphaOmega_
Definition: kOmega.H:87
phi
Definition: correctPhi.H:3
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmega.H:105
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:207
volScalarField k_
Definition: kOmega.H:92
dimensionedScalar beta_
Definition: kOmega.H:84
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
virtual void correctNut()
Definition: kOmega.C:41
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:187
virtual ~kOmega()
Destructor.
Definition: kOmega.H:129
dimensionedScalar betaStar_
Definition: kOmega.H:83
kOmega(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: kOmega.C:82
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:50
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
Definition: kOmega.H:165
A class for managing temporary objects.
Definition: PtrList.H:53
volScalarField omega_
Definition: kOmega.H:93
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: RASModel.H:208
Standard high Reynolds-number k-omega turbulence model for incompressible and compressible flows...
Definition: kOmega.H:72
dimensionedScalar gamma_
Definition: kOmega.H:85
Namespace for OpenFOAM.
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmega.C:65