ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / RAS / kOmegaSST / kOmegaSST.H
blobf4798d852d2e0e7054772508187617a6d4d6a92b
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::compressible::RASModels::kOmegaSST
27 Description
28     Implementation of the k-omega-SST turbulence model for compressible flows.
30     Turbulence model described in:
31     \verbatim
32         Menter, F., Esch, T.
33         "Elements of Industrial Heat Transfer Prediction"
34         16th Brazilian Congress of Mechanical Engineering (COBEM),
35         Nov. 2001
36     \endverbatim
38     Note that this implementation is written in terms of alpha diffusion
39     coefficients rather than the more traditional sigma (alpha = 1/sigma) so
40     that the blending can be applied to all coefficuients in a consistent
41     manner.  The paper suggests that sigma is blended but this would not be
42     consistent with the blending of the k-epsilon and k-omega models.
44     Also note that the error in the last term of equation (2) relating to
45     sigma has been corrected.
47     Wall-functions are applied in this implementation by using equations (14)
48     to specify the near-wall omega as appropriate.
50     The blending functions (15) and (16) are not currently used because of the
51     uncertainty in their origin, range of applicability and that is y+ becomes
52     sufficiently small blending u_tau in this manner clearly becomes nonsense.
54     The default model coefficients correspond to the following:
55     \verbatim
56         kOmegaSSTCoeffs
57         {
58             alphaK1     0.85034;
59             alphaK2     1.0;
60             alphaOmega1 0.5;
61             alphaOmega2 0.85616;
62             Prt         1.0;    // only for compressible
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 compressiblekOmegaSST_H
79 #define compressiblekOmegaSST_H
81 #include "RASModel.H"
82 #include "wallDist.H"
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 namespace Foam
88 namespace compressible
90 namespace RASModels
93 /*---------------------------------------------------------------------------*\
94                           Class kOmegaSST Declaration
95 \*---------------------------------------------------------------------------*/
97 class kOmegaSST
99     public RASModel
102 protected:
104     // Protected data
106         // Model coefficients
108             dimensionedScalar alphaK1_;
109             dimensionedScalar alphaK2_;
111             dimensionedScalar alphaOmega1_;
112             dimensionedScalar alphaOmega2_;
114             dimensionedScalar Prt_;
116             dimensionedScalar gamma1_;
117             dimensionedScalar gamma2_;
119             dimensionedScalar beta1_;
120             dimensionedScalar beta2_;
122             dimensionedScalar betaStar_;
124             dimensionedScalar a1_;
125             dimensionedScalar c1_;
128         //- Wall distance
129         //  Note: different to wall distance in parent RASModel
130         wallDist y_;
132         // Fields
134             volScalarField k_;
135             volScalarField omega_;
136             volScalarField mut_;
137             volScalarField alphat_;
140     // Private Member Functions
142         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
143         tmp<volScalarField> F2() const;
145         tmp<volScalarField> blend
146         (
147             const volScalarField& F1,
148             const dimensionedScalar& psi1,
149             const dimensionedScalar& psi2
150         ) const
151         {
152             return F1*(psi1 - psi2) + psi2;
153         }
155         tmp<volScalarField> alphaK(const volScalarField& F1) const
156         {
157             return blend(F1, alphaK1_, alphaK2_);
158         }
160         tmp<volScalarField> alphaOmega(const volScalarField& F1) const
161         {
162             return blend(F1, alphaOmega1_, alphaOmega2_);
163         }
165         tmp<volScalarField> beta(const volScalarField& F1) const
166         {
167             return blend(F1, beta1_, beta2_);
168         }
170         tmp<volScalarField> gamma(const volScalarField& F1) const
171         {
172             return blend(F1, gamma1_, gamma2_);
173         }
176 public:
178     //- Runtime type information
179     TypeName("kOmegaSST");
182     // Constructors
184         //- Construct from components
185         kOmegaSST
186         (
187             const volScalarField& rho,
188             const volVectorField& U,
189             const surfaceScalarField& phi,
190             const basicThermo& thermophysicalModel,
191             const word& turbulenceModelName = turbulenceModel::typeName,
192             const word& modelName = typeName
193         );
196     //- Destructor
197     virtual ~kOmegaSST()
198     {}
201     // Member Functions
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)*mut_ + mu())
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)*mut_ + mu())
218             );
219         }
221         virtual tmp<volScalarField> mut() const
222         {
223             return mut_;
224         }
226         //- Return the turbulence thermal diffusivity
227         virtual tmp<volScalarField> alphat() const
228         {
229             return alphat_;
230         }
232         //- Return the turbulence kinetic energy
233         virtual tmp<volScalarField> k() const
234         {
235             return k_;
236         }
238         virtual tmp<volScalarField> omega() const
239         {
240             return omega_;
241         }
243         //- Return the turbulence kinetic energy dissipation rate
244         virtual tmp<volScalarField> epsilon() const
245         {
246             return tmp<volScalarField>
247             (
248                 new volScalarField
249                 (
250                     IOobject
251                     (
252                         "epsilon",
253                         mesh_.time().timeName(),
254                         mesh_
255                     ),
256                     betaStar_*k_*omega_,
257                     omega_.boundaryField().types()
258                 )
259             );
260         }
262         //- Return the Reynolds stress tensor
263         virtual tmp<volSymmTensorField> R() const;
265         //- Return the effective stress tensor including the laminar stress
266         virtual tmp<volSymmTensorField> devRhoReff() const;
268         //- Return the source term for the momentum equation
269         virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const;
271         //- Solve the turbulence equations and correct the turbulence viscosity
272         virtual void correct();
274         //- Read RASProperties dictionary
275         virtual bool read();
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 } // End namespace RASModels
282 } // End namespace compressible
283 } // End namespace Foam
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 #endif
289 // ************************************************************************* //