BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / kOmega / kOmega.H
blob10729f334e714713cfb9f710530415bb6a2e2554
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Class
25     Foam::incompressible::RASModels::kOmega
27 Description
28     Standard high Reynolds-number k-omega turbulence model for
29     incompressible flows.
31     References:
32     \verbatim
33         "Turbulence Modeling for CFD"
34         D. C. Wilcox,
35         DCW Industries, Inc., La Canada,
36         California, 1988.
38         See also:
39         http://www.cfd-online.com/Wiki/Wilcox's_k-omega_model
40     \endverbatim
42     The default model coefficients correspond to the following:
43     \verbatim
44         kOmegaCoeffs
45         {
46             Cmu         0.09;  // Equivalent to betaStar
47             alpha       0.52;
48             beta        0.072;
49             alphak      0.5;
50             alphaOmega  0.5;
51         }
52     \endverbatim
54 SourceFiles
55     kOmega.C
57 \*---------------------------------------------------------------------------*/
59 #ifndef kOmega_H
60 #define kOmega_H
62 #include "RASModel.H"
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 namespace Foam
68 namespace incompressible
70 namespace RASModels
73 /*---------------------------------------------------------------------------*\
74                            Class kOmega Declaration
75 \*---------------------------------------------------------------------------*/
77 class kOmega
79     public RASModel
82 protected:
84     // Protected data
86         // Model coefficients
88             dimensionedScalar Cmu_;
89             dimensionedScalar beta_;
90             dimensionedScalar alpha_;
91             dimensionedScalar alphaK_;
92             dimensionedScalar alphaOmega_;
95         // Fields
97             volScalarField k_;
98             volScalarField omega_;
99             volScalarField nut_;
102 public:
104     //- Runtime type information
105     TypeName("kOmega");
107     // Constructors
109         //- Construct from components
110         kOmega
111         (
112             const volVectorField& U,
113             const surfaceScalarField& phi,
114             transportModel& transport,
115             const word& turbulenceModelName = turbulenceModel::typeName,
116             const word& modelName = typeName
117         );
120     //- Destructor
121     virtual ~kOmega()
122     {}
125     // Member Functions
127         //- Return the turbulence viscosity
128         virtual tmp<volScalarField> nut() const
129         {
130             return nut_;
131         }
133         //- Return the effective diffusivity for k
134         tmp<volScalarField> DkEff() const
135         {
136             return tmp<volScalarField>
137             (
138                 new volScalarField("DkEff", alphaK_*nut_ + nu())
139             );
140         }
142         //- Return the effective diffusivity for omega
143         tmp<volScalarField> DomegaEff() const
144         {
145             return tmp<volScalarField>
146             (
147                 new volScalarField("DomegaEff", alphaOmega_*nut_ + nu())
148             );
149         }
151         //- Return the turbulence kinetic energy
152         virtual tmp<volScalarField> k() const
153         {
154             return k_;
155         }
157         //- Return the turbulence specific dissipation rate
158         virtual tmp<volScalarField> omega() const
159         {
160             return omega_;
161         }
163         //- Return the turbulence kinetic energy dissipation rate
164         virtual tmp<volScalarField> epsilon() const
165         {
166             return tmp<volScalarField>
167             (
168                 new volScalarField
169                 (
170                     IOobject
171                     (
172                         "epsilon",
173                         mesh_.time().timeName(),
174                         mesh_
175                     ),
176                     Cmu_*k_*omega_,
177                     omega_.boundaryField().types()
178                 )
179             );
180         }
182         //- Return the Reynolds stress tensor
183         virtual tmp<volSymmTensorField> R() const;
185         //- Return the effective stress tensor including the laminar stress
186         virtual tmp<volSymmTensorField> devReff() const;
188         //- Return the source term for the momentum equation
189         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
191         //- Solve the turbulence equations and correct the turbulence viscosity
192         virtual void correct();
194         //- Read RASProperties dictionary
195         virtual bool read();
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 } // End namespace RASModels
202 } // End namespace incompressible
203 } // End namespace Foam
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 #endif
209 // ************************************************************************* //