BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / kOmegaSST / kOmegaSST.H
blobbd4b3742b5e4ee48b042f12bd1a8da8f129f3f7c
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::kOmegaSST
27 Description
28     Implementation of the k-omega-SST turbulence model for incompressible
29     flows.
31     Turbulence model described in:
32     \verbatim
33         Menter, F., Esch, T.
34         "Elements of Industrial Heat Transfer Prediction"
35         16th Brazilian Congress of Mechanical Engineering (COBEM),
36         Nov. 2001
37     \endverbatim
39     Note that this implementation is written in terms of alpha diffusion
40     coefficients rather than the more traditional sigma (alpha = 1/sigma) so
41     that the blending can be applied to all coefficuients in a consistent
42     manner.  The paper suggests that sigma is blended but this would not be
43     consistent with the blending of the k-epsilon and k-omega models.
45     Also note that the error in the last term of equation (2) relating to
46     sigma has been corrected.
48     Wall-functions are applied in this implementation by using equations (14)
49     to specify the near-wall omega as appropriate.
51     The blending functions (15) and (16) are not currently used because of the
52     uncertainty in their origin, range of applicability and that is y+ becomes
53     sufficiently small blending u_tau in this manner clearly becomes nonsense.
55     The default model coefficients correspond to the following:
56     \verbatim
57         kOmegaSSTCoeffs
58         {
59             alphaK1     0.85034;
60             alphaK2     1.0;
61             alphaOmega1 0.5;
62             alphaOmega2 0.85616;
63             beta1       0.075;
64             beta2       0.0828;
65             betaStar    0.09;
66             gamma1      0.5532;
67             gamma2      0.4403;
68             a1          0.31;
69             c1          10.0;
70         }
71     \endverbatim
73 SourceFiles
74     kOmegaSST.C
76 \*---------------------------------------------------------------------------*/
78 #ifndef kOmegaSST_H
79 #define kOmegaSST_H
81 #include "RASModel.H"
82 #include "wallDist.H"
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 namespace Foam
88 namespace incompressible
90 namespace RASModels
93 /*---------------------------------------------------------------------------*\
94                           Class kOmegaSST Declaration
95 \*---------------------------------------------------------------------------*/
97 class kOmegaSST
99     public RASModel
102 protected:
104     // Protected data
106         // Model coefficients
107             dimensionedScalar alphaK1_;
108             dimensionedScalar alphaK2_;
110             dimensionedScalar alphaOmega1_;
111             dimensionedScalar alphaOmega2_;
113             dimensionedScalar gamma1_;
114             dimensionedScalar gamma2_;
116             dimensionedScalar beta1_;
117             dimensionedScalar beta2_;
119             dimensionedScalar betaStar_;
121             dimensionedScalar a1_;
122             dimensionedScalar c1_;
124         //- Wall distance field
125         //  Note: different to wall distance in parent RASModel
126         wallDist y_;
128         // Fields
130             volScalarField k_;
131             volScalarField omega_;
132             volScalarField nut_;
135     // Protected Member Functions
137         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
138         tmp<volScalarField> F2() const;
140         tmp<volScalarField> blend
141         (
142             const volScalarField& F1,
143             const dimensionedScalar& psi1,
144             const dimensionedScalar& psi2
145         ) const
146         {
147             return F1*(psi1 - psi2) + psi2;
148         }
150         tmp<volScalarField> alphaK(const volScalarField& F1) const
151         {
152             return blend(F1, alphaK1_, alphaK2_);
153         }
155         tmp<volScalarField> alphaOmega(const volScalarField& F1) const
156         {
157             return blend(F1, alphaOmega1_, alphaOmega2_);
158         }
160         tmp<volScalarField> beta(const volScalarField& F1) const
161         {
162             return blend(F1, beta1_, beta2_);
163         }
165         tmp<volScalarField> gamma(const volScalarField& F1) const
166         {
167             return blend(F1, gamma1_, gamma2_);
168         }
171 public:
173     //- Runtime type information
174     TypeName("kOmegaSST");
177     // Constructors
179         //- Construct from components
180         kOmegaSST
181         (
182             const volVectorField& U,
183             const surfaceScalarField& phi,
184             transportModel& transport,
185             const word& turbulenceModelName = turbulenceModel::typeName,
186             const word& modelName = typeName
187         );
190     //- Destructor
191     virtual ~kOmegaSST()
192     {}
195     // Member Functions
197         //- Return the turbulence viscosity
198         virtual tmp<volScalarField> nut() const
199         {
200             return nut_;
201         }
203         //- Return the effective diffusivity for k
204         tmp<volScalarField> DkEff(const volScalarField& F1) const
205         {
206             return tmp<volScalarField>
207             (
208                 new volScalarField("DkEff", alphaK(F1)*nut_ + nu())
209             );
210         }
212         //- Return the effective diffusivity for omega
213         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
214         {
215             return tmp<volScalarField>
216             (
217                 new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu())
218             );
219         }
221         //- Return the turbulence kinetic energy
222         virtual tmp<volScalarField> k() const
223         {
224             return k_;
225         }
227         //- Return the turbulence specific dissipation rate
228         virtual tmp<volScalarField> omega() const
229         {
230             return omega_;
231         }
233         //- Return the turbulence kinetic energy dissipation rate
234         virtual tmp<volScalarField> epsilon() const
235         {
236             return tmp<volScalarField>
237             (
238                 new volScalarField
239                 (
240                     IOobject
241                     (
242                         "epsilon",
243                         mesh_.time().timeName(),
244                         mesh_
245                     ),
246                     betaStar_*k_*omega_,
247                     omega_.boundaryField().types()
248                 )
249             );
250         }
252         //- Return the Reynolds stress tensor
253         virtual tmp<volSymmTensorField> R() const;
255         //- Return the effective stress tensor including the laminar stress
256         virtual tmp<volSymmTensorField> devReff() const;
258         //- Return the source term for the momentum equation
259         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
261         //- Solve the turbulence equations and correct the turbulence viscosity
262         virtual void correct();
264         //- Read RASProperties dictionary
265         virtual bool read();
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 } // End namespace RASModels
272 } // namespace incompressible
273 } // End namespace Foam
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 #endif
279 // ************************************************************************* //