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::incompressible::RASModels::NonlinearKEShih 00026 00027 Description 00028 Shih's quadratic non-linear k-epsilon turbulence model for 00029 incompressible flows 00030 00031 SourceFiles 00032 NonlinearKEShih.C 00033 00034 \*---------------------------------------------------------------------------*/ 00035 00036 #ifndef NonlinearKEShih_H 00037 #define NonlinearKEShih_H 00038 00039 #include <incompressibleRASModels/RASModel.H> 00040 00041 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00042 00043 namespace Foam 00044 { 00045 namespace incompressible 00046 { 00047 namespace RASModels 00048 { 00049 00050 /*---------------------------------------------------------------------------*\ 00051 Class NonlinearKEShih Declaration 00052 \*---------------------------------------------------------------------------*/ 00053 00054 class NonlinearKEShih 00055 : 00056 public RASModel 00057 { 00058 // Private data 00059 00060 // Model coefficients 00061 00062 dimensionedScalar C1_; 00063 dimensionedScalar C2_; 00064 dimensionedScalar sigmak_; 00065 dimensionedScalar sigmaEps_; 00066 dimensionedScalar A1_; 00067 dimensionedScalar A2_; 00068 dimensionedScalar Ctau1_; 00069 dimensionedScalar Ctau2_; 00070 dimensionedScalar Ctau3_; 00071 dimensionedScalar alphaKsi_; 00072 00073 dimensionedScalar kappa_; 00074 dimensionedScalar E_; 00075 00076 00077 // Fields 00078 00079 volScalarField k_; 00080 volScalarField epsilon_; 00081 00082 volTensorField gradU_; 00083 volScalarField eta_; 00084 volScalarField ksi_; 00085 volScalarField Cmu_; 00086 volScalarField fEta_; 00087 00088 volScalarField nut_; 00089 00090 volSymmTensorField nonlinearStress_; 00091 00092 00093 public: 00094 00095 //- Runtime type information 00096 TypeName("NonlinearKEShih"); 00097 00098 // Constructors 00099 00100 //- Construct from components 00101 NonlinearKEShih 00102 ( 00103 const volVectorField& U, 00104 const surfaceScalarField& phi, 00105 transportModel& transport 00106 ); 00107 00108 00109 //- Destructor 00110 virtual ~NonlinearKEShih() 00111 {} 00112 00113 00114 // Member Functions 00115 00116 //- Return the turbulence viscosity 00117 virtual tmp<volScalarField> nut() const 00118 { 00119 return nut_; 00120 } 00121 00122 //- Return the effective diffusivity for k 00123 tmp<volScalarField> DkEff() const 00124 { 00125 return tmp<volScalarField> 00126 ( 00127 new volScalarField("DkEff", nut_/sigmak_ + nu()) 00128 ); 00129 } 00130 00131 //- Return the effective diffusivity for epsilon 00132 tmp<volScalarField> DepsilonEff() const 00133 { 00134 return tmp<volScalarField> 00135 ( 00136 new volScalarField("DepsilonEff", nut_/sigmaEps_ + nu()) 00137 ); 00138 } 00139 00140 //- Return the turbulence kinetic energy 00141 virtual tmp<volScalarField> k() const 00142 { 00143 return k_; 00144 } 00145 00146 //- Return the turbulence kinetic energy dissipation rate 00147 virtual tmp<volScalarField> epsilon() const 00148 { 00149 return epsilon_; 00150 } 00151 00152 //- Return the Reynolds stress tensor 00153 virtual tmp<volSymmTensorField> R() const; 00154 00155 //- Return the effective stress tensor including the laminar stress 00156 virtual tmp<volSymmTensorField> devReff() const; 00157 00158 //- Return the source term for the momentum equation 00159 virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const; 00160 00161 //- Solve the turbulence equations and correct the turbulence viscosity 00162 virtual void correct(); 00163 00164 //- Read RASProperties dictionary 00165 virtual bool read(); 00166 }; 00167 00168 00169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00170 00171 } // End namespace RASModels 00172 } // End namespace incompressible 00173 } // End namespace Foam 00174 00175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00176 00177 #endif 00178 00179 // ************************ vim: set sw=4 sts=4 et: ************************ //