OpenFOAM logo
Open Source CFD Toolkit

LienCubicKE.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     LienCubicKE
00027 
00028 Description
00029     Lien cubic non-linear k-epsilon turbulence model.
00030 
00031 SourceFiles
00032     LienCubicKE.C
00033     LienCubicKECorrect.C
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #ifndef LienCubicKE_H
00038 #define LienCubicKE_H
00039 
00040 #include "turbulenceModel.H"
00041 
00042 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00043 
00044 namespace Foam
00045 {
00046 namespace turbulenceModels
00047 {
00048 
00049 /*---------------------------------------------------------------------------*\
00050                            Class LienCubicKE Declaration
00051 \*---------------------------------------------------------------------------*/
00052 
00053 class LienCubicKE
00054 :
00055     public turbulenceModel
00056 {
00057     // Private data
00058 
00059         dimensionedScalar C1;
00060         dimensionedScalar C2;
00061         dimensionedScalar alphak;
00062         dimensionedScalar alphaEps;
00063         dimensionedScalar A1;
00064         dimensionedScalar A2;
00065         dimensionedScalar Ctau1;
00066         dimensionedScalar Ctau2;
00067         dimensionedScalar Ctau3;
00068         dimensionedScalar alphaKsi;
00069 
00070         volScalarField k_;
00071         volScalarField epsilon_;
00072 
00073         volTensorField gradU;
00074         volScalarField eta;
00075         volScalarField ksi;
00076         volScalarField Cmu;
00077         volScalarField fEta;
00078         volScalarField C5viscosity;
00079 
00080         volScalarField nut_;
00081 
00082         volTensorField nonlinearStress;
00083 
00084 
00085 public:
00086 
00087     //- Runtime type information
00088     TypeName("LienCubicKE");
00089 
00090     // Constructors
00091 
00092         //- from components
00093         LienCubicKE
00094         (
00095             const volVectorField& U,
00096             const surfaceScalarField& phi,
00097             transportModel& transport
00098         );
00099 
00100 
00101     // Destructor
00102 
00103         ~LienCubicKE(){}
00104 
00105 
00106     // Member Functions
00107 
00108         tmp<volScalarField> nut() const
00109         {
00110             return nut_;
00111         }
00112 
00113         //- Return the effective diffusivity for k
00114         tmp<volScalarField> DkEff() const
00115         {
00116             return tmp<volScalarField> 
00117             (
00118                 new volScalarField("DkEff", alphak*nut_ + nu())
00119             );
00120         }
00121 
00122         //- Return the effective diffusivity for epsilon
00123         tmp<volScalarField> DepsilonEff() const
00124         {
00125             return tmp<volScalarField> 
00126             (
00127                 new volScalarField("DepsilonEff", alphaEps*nut_ + nu())
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 
00142         tmp<volTensorField> R() const;
00143 
00144         tmp<fvVectorMatrix> divR(volVectorField& U) const;
00145 
00146         void correct();
00147 
00148         //- Read turbulenceProperties dictionary
00149         bool read();
00150 };
00151 
00152 
00153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00154 
00155 } // End namespace turbulenceModels
00156 } // End namespace Foam
00157 
00158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00159 
00160 #endif
00161 
00162 // ************************************************************************* //
For further information go to www.openfoam.org