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