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