FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

kOmegaSST.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Class
00025     Foam::incompressible::RASModels::kOmegaSST
00026 
00027 Description
00028     Implementation of the k-omega-SST turbulence model for incompressible
00029     flows.
00030 
00031     Turbulence model described in:
00032     @verbatim
00033         Menter, F., Esch, T.
00034         "Elements of Industrial Heat Transfer Prediction"
00035         16th Brazilian Congress of Mechanical Engineering (COBEM),
00036         Nov. 2001
00037     @endverbatim
00038 
00039     Note that this implementation is written in terms of alpha diffusion
00040     coefficients rather than the more traditional sigma (alpha = 1/sigma) so
00041     that the blending can be applied to all coefficuients in a consistent
00042     manner.  The paper suggests that sigma is blended but this would not be
00043     consistent with the blending of the k-epsilon and k-omega models.
00044 
00045     Also note that the error in the last term of equation (2) relating to
00046     sigma has been corrected.
00047 
00048     Wall-functions are applied in this implementation by using equations (14)
00049     to specify the near-wall omega as appropriate.
00050 
00051     The blending functions (15) and (16) are not currently used because of the
00052     uncertainty in their origin, range of applicability and that is y+ becomes
00053     sufficiently small blending u_tau in this manner clearly becomes nonsense.
00054 
00055     The default model coefficients correspond to the following:
00056     @verbatim
00057         kOmegaSSTCoeffs
00058         {
00059             alphaK1     0.85034;
00060             alphaK2     1.0;
00061             alphaOmega1 0.5;
00062             alphaOmega2 0.85616;
00063             beta1       0.075;
00064             beta2       0.0828;
00065             betaStar    0.09;
00066             gamma1      0.5532;
00067             gamma2      0.4403;
00068             a1          0.31;
00069             c1          10.0;
00070         }
00071     @endverbatim
00072 
00073 SourceFiles
00074     kOmegaSST.C
00075 
00076 \*---------------------------------------------------------------------------*/
00077 
00078 #ifndef kOmegaSST_H
00079 #define kOmegaSST_H
00080 
00081 #include <incompressibleRASModels/RASModel.H>
00082 #include <finiteVolume/wallDist.H>
00083 
00084 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00085 
00086 namespace Foam
00087 {
00088 namespace incompressible
00089 {
00090 namespace RASModels
00091 {
00092 
00093 /*---------------------------------------------------------------------------*\
00094                            Class kOmega Declaration
00095 \*---------------------------------------------------------------------------*/
00096 
00097 class kOmegaSST
00098 :
00099     public RASModel
00100 {
00101     // Private data
00102 
00103         // Model coefficients
00104             dimensionedScalar alphaK1_;
00105             dimensionedScalar alphaK2_;
00106 
00107             dimensionedScalar alphaOmega1_;
00108             dimensionedScalar alphaOmega2_;
00109 
00110             dimensionedScalar gamma1_;
00111             dimensionedScalar gamma2_;
00112 
00113             dimensionedScalar beta1_;
00114             dimensionedScalar beta2_;
00115 
00116             dimensionedScalar betaStar_;
00117 
00118             dimensionedScalar a1_;
00119             dimensionedScalar c1_;
00120 
00121         //- Wall distance field
00122         //  Note: different to wall distance in parent RASModel
00123         wallDist y_;
00124 
00125         // Fields
00126 
00127             volScalarField k_;
00128             volScalarField omega_;
00129             volScalarField nut_;
00130 
00131 
00132     // Private member functions
00133 
00134         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
00135         tmp<volScalarField> F2() const;
00136 
00137         tmp<volScalarField> blend
00138         (
00139             const volScalarField& F1,
00140             const dimensionedScalar& psi1,
00141             const dimensionedScalar& psi2
00142         ) const
00143         {
00144             return F1*(psi1 - psi2) + psi2;
00145         }
00146 
00147         tmp<volScalarField> alphaK
00148         (
00149             const volScalarField& F1
00150         ) const
00151         {
00152             return blend(F1, alphaK1_, alphaK2_);
00153         }
00154 
00155         tmp<volScalarField> alphaOmega
00156         (
00157             const volScalarField& F1
00158         ) const
00159         {
00160             return blend(F1, alphaOmega1_, alphaOmega2_);
00161         }
00162 
00163         tmp<volScalarField> beta
00164         (
00165             const volScalarField& F1
00166         ) const
00167         {
00168             return blend(F1, beta1_, beta2_);
00169         }
00170 
00171         tmp<volScalarField> gamma
00172         (
00173             const volScalarField& F1
00174         ) const
00175         {
00176             return blend(F1, gamma1_, gamma2_);
00177         }
00178 
00179 
00180 public:
00181 
00182     //- Runtime type information
00183     TypeName("kOmegaSST");
00184 
00185 
00186     // Constructors
00187 
00188         //- Construct from components
00189         kOmegaSST
00190         (
00191             const volVectorField& U,
00192             const surfaceScalarField& phi,
00193             transportModel& transport
00194         );
00195 
00196 
00197     //- Destructor
00198     virtual ~kOmegaSST()
00199     {}
00200 
00201 
00202     // Member Functions
00203 
00204         //- Return the turbulence viscosity
00205         virtual tmp<volScalarField> nut() const
00206         {
00207             return nut_;
00208         }
00209 
00210         //- Return the effective diffusivity for k
00211         tmp<volScalarField> DkEff(const volScalarField& F1) const
00212         {
00213             return tmp<volScalarField>
00214             (
00215                 new volScalarField("DkEff", alphaK(F1)*nut_ + nu())
00216             );
00217         }
00218 
00219         //- Return the effective diffusivity for omega
00220         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
00221         {
00222             return tmp<volScalarField>
00223             (
00224                 new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu())
00225             );
00226         }
00227 
00228         //- Return the turbulence kinetic energy
00229         virtual tmp<volScalarField> k() const
00230         {
00231             return k_;
00232         }
00233 
00234         //- Return the turbulence specific dissipation rate
00235         virtual tmp<volScalarField> omega() const
00236         {
00237             return omega_;
00238         }
00239 
00240         //- Return the turbulence kinetic energy dissipation rate
00241         virtual tmp<volScalarField> epsilon() const
00242         {
00243             return tmp<volScalarField>
00244             (
00245                 new volScalarField
00246                 (
00247                     IOobject
00248                     (
00249                         "epsilon",
00250                         mesh_.time().timeName(),
00251                         mesh_
00252                     ),
00253                     betaStar_*k_*omega_,
00254                     omega_.boundaryField().types()
00255                 )
00256             );
00257         }
00258 
00259         //- Return the Reynolds stress tensor
00260         virtual tmp<volSymmTensorField> R() const;
00261 
00262         //- Return the effective stress tensor including the laminar stress
00263         virtual tmp<volSymmTensorField> devReff() const;
00264 
00265         //- Return the source term for the momentum equation
00266         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
00267 
00268         //- Solve the turbulence equations and correct the turbulence viscosity
00269         virtual void correct();
00270 
00271         //- Read RASProperties dictionary
00272         virtual bool read();
00273 };
00274 
00275 
00276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00277 
00278 } // End namespace RASModels
00279 } // namespace incompressible
00280 } // End namespace Foam
00281 
00282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00283 
00284 #endif
00285 
00286 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines