![]() |
|
00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2005 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 the 00013 Free Software Foundation; either version 2 of the License, or (at your 00014 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, write to the Free Software Foundation, 00023 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA 00024 00025 Class 00026 realizableKE 00027 00028 Description 00029 Realizable k-epsilon turbulence model from the paper: 00030 00031 "A New k-epsilon Eddy Viscosity Model for High Reynolds Number 00032 Turbulent Flows" 00033 00034 Tsan-Hsing Shih, William W. Liou, Aamir Shabbir, Zhigang Tang and 00035 Jiang Zhu 00036 00037 Computers and Fluids Vol. 24, No. 3, pp. 227-238, 1995 00038 00039 SourceFiles 00040 realizableKE.C 00041 00042 \*---------------------------------------------------------------------------*/ 00043 00044 #ifndef realizableKE_H 00045 #define realizableKE_H 00046 00047 #include "turbulenceModel.H" 00048 00049 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00050 00051 namespace Foam 00052 { 00053 namespace compressible 00054 { 00055 namespace turbulenceModels 00056 { 00057 00058 /*---------------------------------------------------------------------------*\ 00059 Class realizableKE Declaration 00060 \*---------------------------------------------------------------------------*/ 00061 00062 class realizableKE 00063 : 00064 public turbulenceModel 00065 { 00066 // Private data 00067 00068 dimensionedScalar Cmu; 00069 dimensionedScalar A0; 00070 dimensionedScalar C2; 00071 dimensionedScalar alphak; 00072 dimensionedScalar alphaEps; 00073 dimensionedScalar alphah; 00074 00075 volScalarField k_; 00076 volScalarField epsilon_; 00077 volScalarField mut_; 00078 00079 tmp<volScalarField> rCmu 00080 ( 00081 const volTensorField& gradU, 00082 const volScalarField& S2, 00083 const volScalarField& magS 00084 ); 00085 00086 tmp<volScalarField> rCmu 00087 ( 00088 const volTensorField& gradU 00089 ); 00090 00091 public: 00092 00093 //- Runtime type information 00094 TypeName("realizableKE"); 00095 00096 // Constructors 00097 00098 //- from components 00099 realizableKE 00100 ( 00101 const volScalarField& rho, 00102 const volVectorField& U, 00103 const surfaceScalarField& phi, 00104 basicThermo& thermophysicalModel 00105 ); 00106 00107 00108 // Destructor 00109 00110 ~realizableKE(){} 00111 00112 00113 // Member Functions 00114 00115 tmp<volScalarField> mut() const 00116 { 00117 return mut_; 00118 } 00119 00120 //- Return the effective diffusivity for k 00121 tmp<volScalarField> DkEff() const 00122 { 00123 return tmp<volScalarField> 00124 ( 00125 new volScalarField("DkEff", alphak*mut_ + 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", alphaEps*mut_ + mu()) 00135 ); 00136 } 00137 00138 //- Return the effective turbulent thermal diffusivity 00139 tmp<volScalarField> alphaEff() const 00140 { 00141 return tmp<volScalarField> 00142 ( 00143 new volScalarField("alphaEff", alphah*mut_ + alpha()) 00144 ); 00145 } 00146 00147 tmp<volScalarField> k() const 00148 { 00149 return k_; 00150 } 00151 00152 tmp<volScalarField> epsilon() const 00153 { 00154 return epsilon_; 00155 } 00156 00157 tmp<volTensorField> R() const; 00158 00159 tmp<fvVectorMatrix> divRhoR(volVectorField& U) const; 00160 00161 tmp<fvScalarMatrix> divRhoUh(volScalarField& h) const; 00162 00163 void correct(); 00164 00165 //- Read turbulenceProperties dictionary 00166 bool read(); 00167 }; 00168 00169 00170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00171 00172 } // End namespace turbulenceModels 00173 } // End namespace compressible 00174 } // End namespace Foam 00175 00176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00177 00178 #endif 00179 00180 // ************************************************************************* //