1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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/>.
25 Foam::incompressible::RASModels::kOmegaSST
28 Implementation of the k-omega-SST turbulence model for incompressible
31 Turbulence model described in:
34 "Elements of Industrial Heat Transfer Prediction"
35 16th Brazilian Congress of Mechanical Engineering (COBEM),
39 with updated coefficients from
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,
47 but with the consistent production terms from the 2001 paper as form in the
48 2003 paper is a typo, see
50 http://turbmodels.larc.nasa.gov/sst.html
53 and the addition of the optional F3 term for rough walls from
56 "Some Improvements in Menter’s k-omega-SST turbulence model"
57 29th AIAA Fluid Dynamics Conference,
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:
94 \*---------------------------------------------------------------------------*/
100 #include "wallDist.H"
102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 namespace incompressible
111 /*---------------------------------------------------------------------------*\
112 Class kOmega Declaration
113 \*---------------------------------------------------------------------------*/
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_;
143 //- Wall distance field
144 // Note: different to wall distance in parent RASModel
151 volScalarField omega_;
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
164 const volScalarField& F1,
165 const dimensionedScalar& psi1,
166 const dimensionedScalar& psi2
169 return F1*(psi1 - psi2) + psi2;
172 tmp<volScalarField> alphaK
174 const volScalarField& F1
177 return blend(F1, alphaK1_, alphaK2_);
180 tmp<volScalarField> alphaOmega
182 const volScalarField& F1
185 return blend(F1, alphaOmega1_, alphaOmega2_);
188 tmp<volScalarField> beta
190 const volScalarField& F1
193 return blend(F1, beta1_, beta2_);
196 tmp<volScalarField> gamma
198 const volScalarField& F1
201 return blend(F1, gamma1_, gamma2_);
207 //- Runtime type information
208 TypeName("kOmegaSST");
213 //- Construct from components
216 const volVectorField& U,
217 const surfaceScalarField& phi,
218 transportModel& transport
229 //- Return the turbulence viscosity
230 virtual tmp<volScalarField> nut() const
235 //- Return the effective diffusivity for k
236 tmp<volScalarField> DkEff(const volScalarField& F1) const
238 return tmp<volScalarField>
240 new volScalarField("DkEff", alphaK(F1)*nut_ + nu())
244 //- Return the effective diffusivity for omega
245 tmp<volScalarField> DomegaEff(const volScalarField& F1) const
247 return tmp<volScalarField>
249 new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu())
253 //- Return the turbulence kinetic energy
254 virtual tmp<volScalarField> k() const
259 //- Return the turbulence specific dissipation rate
260 virtual tmp<volScalarField> omega() const
265 //- Return the turbulence kinetic energy dissipation rate
266 virtual tmp<volScalarField> epsilon() const
268 return tmp<volScalarField>
275 mesh_.time().timeName(),
279 omega_.boundaryField().types()
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
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 } // End namespace RASModels
304 } // namespace incompressible
305 } // End namespace Foam
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 // ************************************************************************* //