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