OpenFOAM logo
Open Source CFD Toolkit

kEpsilon.H

Go to the documentation of this file.
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     kEpsilon
00027 
00028 Description
00029     Standard k-epsilon turbulence model.
00030 
00031 SourceFiles
00032     kEpsilon.C
00033     kEpsilonCorrect.C
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #ifndef compressiblekEpsilon_H
00038 #define compressiblekEpsilon_H
00039 
00040 #include "turbulenceModel.H"
00041 
00042 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00043 
00044 namespace Foam
00045 {
00046 namespace compressible
00047 {
00048 namespace turbulenceModels
00049 {
00050 
00051 /*---------------------------------------------------------------------------*\
00052                            Class kEpsilon Declaration
00053 \*---------------------------------------------------------------------------*/
00054 
00055 class kEpsilon
00056 :
00057     public turbulenceModel
00058 {
00059     // Private data
00060 
00061         dimensionedScalar Cmu;
00062         dimensionedScalar C1;
00063         dimensionedScalar C2;
00064         dimensionedScalar C3;
00065         dimensionedScalar alphak;
00066         dimensionedScalar alphaEps;
00067         dimensionedScalar alphah;
00068 
00069         volScalarField k_;
00070         volScalarField epsilon_;
00071         volScalarField mut_;
00072 
00073 
00074 public:
00075 
00076     //- Runtime type information
00077     TypeName("kEpsilon");
00078 
00079     // Constructors
00080 
00081         //- from components
00082         kEpsilon
00083         (
00084             const volScalarField& rho,
00085             const volVectorField& U,
00086             const surfaceScalarField& phi,
00087             basicThermo& thermophysicalModel
00088         );
00089 
00090 
00091     // Destructor
00092 
00093         ~kEpsilon()
00094         {}
00095 
00096 
00097     // Member Functions
00098 
00099         tmp<volScalarField> mut() const
00100         {
00101             return mut_;
00102         }
00103 
00104         //- Return the effective diffusivity for k
00105         tmp<volScalarField> DkEff() const
00106         {
00107             return tmp<volScalarField> 
00108             (
00109                 new volScalarField("DkEff", alphak*mut_ + mu())
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*mut_ + mu())
00119             );
00120         }
00121 
00122         //- Return the effective turbulent thermal diffusivity
00123         tmp<volScalarField> alphaEff() const
00124         {
00125             return tmp<volScalarField> 
00126             (
00127                 new volScalarField("alphaEff", alphah*mut_ + alpha())
00128             );
00129         }
00130 
00131         tmp<volScalarField> k() const
00132         {
00133             return k_;
00134         }
00135 
00136         tmp<volScalarField> epsilon() const
00137         {
00138             return epsilon_;
00139         }
00140 
00141         tmp<volTensorField> R() const;
00142 
00143         tmp<fvVectorMatrix> divRhoR(volVectorField& U) const;
00144 
00145         void correct();
00146 
00147         //- Read turbulenceProperties dictionary
00148         bool read();
00149 };
00150 
00151 
00152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00153 
00154 } // End namespace turbulenceModels
00155 } // End namespace compressible
00156 } // End namespace Foam
00157 
00158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00159 
00160 #endif
00161 
00162 // ************************************************************************* //
For further information go to www.openfoam.org