Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / RAS / kOmegaSST / kOmegaSST.H
blob1bd2e4f630d49a861b7aa25661d21bbfa4acd061
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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     with updated coefficients from
40     @verbatim
41         Menter, F. R., Kuntz, M., and Langtry, R.,
42         "Ten Years of Industrial Experience with the SST Turbulence Model",
43         Turbulence, Heat and Mass Transfer 4, 2003,
44         pp. 625 - 632.
45     @endverbatim
47     but with the consistent production terms from the 2001 paper as form in the
48     2003 paper is a typo, see
49     @verbatim
50         http://turbmodels.larc.nasa.gov/sst.html
51     @endverbatim
53     and the addition of the optional F3 term for rough walls from
54     \verbatim
55         Hellsten, A.
56         "Some Improvements in Menter’s k-omega-SST turbulence model"
57         29th AIAA Fluid Dynamics Conference,
58         AIAA-98-2554,
59         June 1998.
60     \endverbatim
62     Note that this implementation is written in terms of alpha diffusion
63     coefficients rather than the more traditional sigma (alpha = 1/sigma) so
64     that the blending can be applied to all coefficuients in a consistent
65     manner.  The paper suggests that sigma is blended but this would not be
66     consistent with the blending of the k-epsilon and k-omega models.
68     Also note that the error in the last term of equation (2) relating to
69     sigma has been corrected.
71     The default model coefficients correspond to the following:
72     @verbatim
73         kOmegaSSTCoeffs
74         {
75             alphaK1     0.85;
76             alphaK2     1.0;
77             alphaOmega1 0.5;
78             alphaOmega2 0.856;
79             beta1       0.075;
80             beta2       0.0828;
81             betaStar    0.09;
82             gamma1      5/9;
83             gamma2      0.44;
84             a1          0.31;
85             b1          1.0;
86             c1          10.0;
87             F3          no;
88         }
89     @endverbatim
91 SourceFiles
92     kOmegaSST.C
94 \*---------------------------------------------------------------------------*/
96 #ifndef kOmegaSST_H
97 #define kOmegaSST_H
99 #include "RASModel.H"
100 #include "wallDist.H"
102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 namespace Foam
106 namespace incompressible
108 namespace RASModels
111 /*---------------------------------------------------------------------------*\
112                            Class kOmega Declaration
113 \*---------------------------------------------------------------------------*/
115 class kOmegaSST
117     public RASModel
119     // Private data
121         // Model coefficients
122             dimensionedScalar alphaK1_;
123             dimensionedScalar alphaK2_;
125             dimensionedScalar alphaOmega1_;
126             dimensionedScalar alphaOmega2_;
128             dimensionedScalar gamma1_;
129             dimensionedScalar gamma2_;
131             dimensionedScalar beta1_;
132             dimensionedScalar beta2_;
134             dimensionedScalar betaStar_;
136             dimensionedScalar a1_;
137             dimensionedScalar b1_;
138             dimensionedScalar c1_;
140             Switch F3_;
143         //- Wall distance field
144         //  Note: different to wall distance in parent RASModel
145         wallDist y_;
148         // Fields
150             volScalarField k_;
151             volScalarField omega_;
152             volScalarField nut_;
155     // Private member functions
157         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
158         tmp<volScalarField> F2() const;
159         tmp<volScalarField> F3() const;
160         tmp<volScalarField> F23() const;
162         tmp<volScalarField> blend
163         (
164             const volScalarField& F1,
165             const dimensionedScalar& psi1,
166             const dimensionedScalar& psi2
167         ) const
168         {
169             return F1*(psi1 - psi2) + psi2;
170         }
172         tmp<volScalarField> alphaK
173         (
174             const volScalarField& F1
175         ) const
176         {
177             return blend(F1, alphaK1_, alphaK2_);
178         }
180         tmp<volScalarField> alphaOmega
181         (
182             const volScalarField& F1
183         ) const
184         {
185             return blend(F1, alphaOmega1_, alphaOmega2_);
186         }
188         tmp<volScalarField> beta
189         (
190             const volScalarField& F1
191         ) const
192         {
193             return blend(F1, beta1_, beta2_);
194         }
196         tmp<volScalarField> gamma
197         (
198             const volScalarField& F1
199         ) const
200         {
201             return blend(F1, gamma1_, gamma2_);
202         }
205 public:
207     //- Runtime type information
208     TypeName("kOmegaSST");
211     // Constructors
213         //- Construct from components
214         kOmegaSST
215         (
216             const volVectorField& U,
217             const surfaceScalarField& phi,
218             transportModel& transport
219         );
222     //- Destructor
223     virtual ~kOmegaSST()
224     {}
227     // Member Functions
229         //- Return the turbulence viscosity
230         virtual tmp<volScalarField> nut() const
231         {
232             return nut_;
233         }
235         //- Return the effective diffusivity for k
236         tmp<volScalarField> DkEff(const volScalarField& F1) const
237         {
238             return tmp<volScalarField>
239             (
240                 new volScalarField("DkEff", alphaK(F1)*nut_ + nu())
241             );
242         }
244         //- Return the effective diffusivity for omega
245         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
246         {
247             return tmp<volScalarField>
248             (
249                 new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu())
250             );
251         }
253         //- Return the turbulence kinetic energy
254         virtual tmp<volScalarField> k() const
255         {
256             return k_;
257         }
259         //- Return the turbulence specific dissipation rate
260         virtual tmp<volScalarField> omega() const
261         {
262             return omega_;
263         }
265         //- Return the turbulence kinetic energy dissipation rate
266         virtual tmp<volScalarField> epsilon() const
267         {
268             return tmp<volScalarField>
269             (
270                 new volScalarField
271                 (
272                     IOobject
273                     (
274                         "epsilon",
275                         mesh_.time().timeName(),
276                         mesh_
277                     ),
278                     betaStar_*k_*omega_,
279                     omega_.boundaryField().types()
280                 )
281             );
282         }
284         //- Return the Reynolds stress tensor
285         virtual tmp<volSymmTensorField> R() const;
287         //- Return the effective stress tensor including the laminar stress
288         virtual tmp<volSymmTensorField> devReff() const;
290         //- Return the source term for the momentum equation
291         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
293         //- Solve the turbulence equations and correct the turbulence viscosity
294         virtual void correct();
296         //- Read RASProperties dictionary
297         virtual bool read();
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 } // End namespace RASModels
304 } // namespace incompressible
305 } // End namespace Foam
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 #endif
311 // ************************************************************************* //