1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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/>.
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 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:
76 \*---------------------------------------------------------------------------*/
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 namespace incompressible
93 /*---------------------------------------------------------------------------*\
94 Class kOmegaSST Declaration
95 \*---------------------------------------------------------------------------*/
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
131 volScalarField omega_;
135 // Protected Member Functions
137 tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
138 tmp<volScalarField> F2() const;
140 tmp<volScalarField> blend
142 const volScalarField& F1,
143 const dimensionedScalar& psi1,
144 const dimensionedScalar& psi2
147 return F1*(psi1 - psi2) + psi2;
150 tmp<volScalarField> alphaK(const volScalarField& F1) const
152 return blend(F1, alphaK1_, alphaK2_);
155 tmp<volScalarField> alphaOmega(const volScalarField& F1) const
157 return blend(F1, alphaOmega1_, alphaOmega2_);
160 tmp<volScalarField> beta(const volScalarField& F1) const
162 return blend(F1, beta1_, beta2_);
165 tmp<volScalarField> gamma(const volScalarField& F1) const
167 return blend(F1, gamma1_, gamma2_);
173 //- Runtime type information
174 TypeName("kOmegaSST");
179 //- Construct from components
182 const volVectorField& U,
183 const surfaceScalarField& phi,
184 transportModel& transport,
185 const word& turbulenceModelName = turbulenceModel::typeName,
186 const word& modelName = typeName
197 //- Return the turbulence viscosity
198 virtual tmp<volScalarField> nut() const
203 //- Return the effective diffusivity for k
204 tmp<volScalarField> DkEff(const volScalarField& F1) const
206 return tmp<volScalarField>
208 new volScalarField("DkEff", alphaK(F1)*nut_ + nu())
212 //- Return the effective diffusivity for omega
213 tmp<volScalarField> DomegaEff(const volScalarField& F1) const
215 return tmp<volScalarField>
217 new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu())
221 //- Return the turbulence kinetic energy
222 virtual tmp<volScalarField> k() const
227 //- Return the turbulence specific dissipation rate
228 virtual tmp<volScalarField> omega() const
233 //- Return the turbulence kinetic energy dissipation rate
234 virtual tmp<volScalarField> epsilon() const
236 return tmp<volScalarField>
243 mesh_.time().timeName(),
247 omega_.boundaryField().types()
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
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 } // End namespace RASModels
272 } // namespace incompressible
273 } // End namespace Foam
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 // ************************************************************************* //