![]() |
|
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 RNGkEpsilon 00027 00028 Description 00029 Renormalisation group k-epsilon turbulence model. 00030 00031 SourceFiles 00032 RNGkEpsilon.C 00033 RNGkEpsilonCorrect.C 00034 00035 \*---------------------------------------------------------------------------*/ 00036 00037 #ifndef RNGkEpsilon_H 00038 #define RNGkEpsilon_H 00039 00040 #include "turbulenceModel.H" 00041 00042 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00043 00044 namespace Foam 00045 { 00046 namespace turbulenceModels 00047 { 00048 00049 /*---------------------------------------------------------------------------*\ 00050 Class RNGkEpsilon Declaration 00051 \*---------------------------------------------------------------------------*/ 00052 00053 class RNGkEpsilon 00054 : 00055 public turbulenceModel 00056 { 00057 // Private data 00058 00059 dimensionedScalar Cmu; 00060 dimensionedScalar C1; 00061 dimensionedScalar C2; 00062 dimensionedScalar alphak; 00063 dimensionedScalar alphaEps; 00064 dimensionedScalar eta0; 00065 dimensionedScalar beta; 00066 00067 volScalarField k_; 00068 volScalarField epsilon_; 00069 volScalarField nut_; 00070 00071 00072 public: 00073 00074 //- Runtime type information 00075 TypeName("RNGkEpsilon"); 00076 00077 // Constructors 00078 00079 //- from components 00080 RNGkEpsilon 00081 ( 00082 const volVectorField& U, 00083 const surfaceScalarField& phi, 00084 transportModel& transport 00085 ); 00086 00087 00088 // Destructor 00089 00090 ~RNGkEpsilon(){} 00091 00092 00093 // Member Functions 00094 00095 tmp<volScalarField> nut() const 00096 { 00097 return nut_; 00098 } 00099 00100 //- Return the effective diffusivity for k 00101 tmp<volScalarField> DkEff() const 00102 { 00103 return tmp<volScalarField> 00104 ( 00105 new volScalarField("DkEff", alphak*nut_ + nu()) 00106 ); 00107 } 00108 00109 //- Return the effective diffusivity for epsilon 00110 tmp<volScalarField> DepsilonEff() const 00111 { 00112 return tmp<volScalarField> 00113 ( 00114 new volScalarField("DepsilonEff", alphaEps*nut_ + nu()) 00115 ); 00116 } 00117 00118 tmp<volScalarField> k() const 00119 { 00120 return k_; 00121 } 00122 00123 tmp<volScalarField> epsilon() const 00124 { 00125 return epsilon_; 00126 } 00127 00128 tmp<volTensorField> R() const; 00129 00130 tmp<fvVectorMatrix> divR(volVectorField& U) const; 00131 00132 void correct(); 00133 00134 //- Read turbulenceProperties dictionary 00135 bool read(); 00136 }; 00137 00138 00139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00140 00141 } // End namespace turbulenceModels 00142 } // End namespace Foam 00143 00144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00145 00146 #endif 00147 00148 // ************************************************************************* //