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::compressible::RASModels::kEpsilon 00026 00027 Description 00028 Standard k-epsilon turbulence model for compressible flows 00029 00030 The default model coefficients correspond to the following: 00031 @verbatim 00032 kEpsilonCoeffs 00033 { 00034 Cmu 0.09; 00035 C1 1.44; 00036 C2 1.92; 00037 C3 -0.33; // only for compressible 00038 sigmak 1.0; // only for compressible 00039 sigmaEps 1.3; 00040 Prt 1.0; // only for compressible 00041 } 00042 @endverbatim 00043 00044 SourceFiles 00045 kEpsilon.C 00046 kEpsilonCorrect.C 00047 00048 \*---------------------------------------------------------------------------*/ 00049 00050 #ifndef compressiblekEpsilon_H 00051 #define compressiblekEpsilon_H 00052 00053 #include <compressibleRASModels/RASModel.H> 00054 00055 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00056 00057 namespace Foam 00058 { 00059 namespace compressible 00060 { 00061 namespace RASModels 00062 { 00063 00064 /*---------------------------------------------------------------------------*\ 00065 Class kEpsilon Declaration 00066 \*---------------------------------------------------------------------------*/ 00067 00068 class kEpsilon 00069 : 00070 public RASModel 00071 { 00072 // Private data 00073 00074 // Model coefficients 00075 00076 dimensionedScalar Cmu_; 00077 dimensionedScalar C1_; 00078 dimensionedScalar C2_; 00079 dimensionedScalar C3_; 00080 dimensionedScalar sigmak_; 00081 dimensionedScalar sigmaEps_; 00082 dimensionedScalar Prt_; 00083 00084 // Fields 00085 00086 volScalarField k_; 00087 volScalarField epsilon_; 00088 volScalarField mut_; 00089 volScalarField alphat_; 00090 00091 00092 public: 00093 00094 //- Runtime type information 00095 TypeName("kEpsilon"); 00096 00097 // Constructors 00098 00099 //- Construct from components 00100 kEpsilon 00101 ( 00102 const volScalarField& rho, 00103 const volVectorField& U, 00104 const surfaceScalarField& phi, 00105 const basicThermo& thermophysicalModel 00106 ); 00107 00108 00109 //- Destructor 00110 virtual ~kEpsilon() 00111 {} 00112 00113 00114 // Member Functions 00115 00116 //- Return the effective diffusivity for k 00117 tmp<volScalarField> DkEff() const 00118 { 00119 return tmp<volScalarField> 00120 ( 00121 new volScalarField("DkEff", mut_/sigmak_ + mu()) 00122 ); 00123 } 00124 00125 //- Return the effective diffusivity for epsilon 00126 tmp<volScalarField> DepsilonEff() const 00127 { 00128 return tmp<volScalarField> 00129 ( 00130 new volScalarField("DepsilonEff", mut_/sigmaEps_ + mu()) 00131 ); 00132 } 00133 00134 //- Return the turbulence viscosity 00135 virtual tmp<volScalarField> mut() const 00136 { 00137 return mut_; 00138 } 00139 00140 //- Return the effective turbulent thermal diffusivity 00141 virtual tmp<volScalarField> alphaEff() const 00142 { 00143 return tmp<volScalarField> 00144 ( 00145 new volScalarField("alphaEff", alphat_ + alpha()) 00146 ); 00147 } 00148 00149 //- Return the turbulence kinetic energy 00150 virtual tmp<volScalarField> k() const 00151 { 00152 return k_; 00153 } 00154 00155 //- Return the turbulence kinetic energy dissipation rate 00156 virtual tmp<volScalarField> epsilon() const 00157 { 00158 return epsilon_; 00159 } 00160 00161 //- Return the Reynolds stress tensor 00162 virtual tmp<volSymmTensorField> R() const; 00163 00164 //- Return the effective stress tensor including the laminar stress 00165 virtual tmp<volSymmTensorField> devRhoReff() const; 00166 00167 //- Return the source term for the momentum equation 00168 virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const; 00169 00170 //- Solve the turbulence equations and correct the turbulence viscosity 00171 virtual void correct(); 00172 00173 //- Read RASProperties dictionary 00174 virtual bool read(); 00175 }; 00176 00177 00178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00179 00180 } // End namespace RASModels 00181 } // End namespace compressible 00182 } // End namespace Foam 00183 00184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00185 00186 #endif 00187 00188 // ************************ vim: set sw=4 sts=4 et: ************************ //