OpenFOAM logo
Open Source CFD Toolkit

NonlinearKEShih.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     NonlinearKEShih
00027 
00028 Description
00029     Shih's quadratic non-linear k-epsilon turbulence model.
00030 
00031 SourceFiles
00032     NonlinearKEShih.C
00033     NonlinearKEShihCorrect.C
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #ifndef NonlinearKEShih_H
00038 #define NonlinearKEShih_H
00039 
00040 #include "turbulenceModel.H"
00041 
00042 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00043 
00044 namespace Foam
00045 {
00046 namespace turbulenceModels
00047 {
00048 
00049 /*---------------------------------------------------------------------------*\
00050                            Class NonlinearKEShih Declaration
00051 \*---------------------------------------------------------------------------*/
00052 
00053 class NonlinearKEShih
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 
00079         volScalarField nut_;
00080 
00081         volTensorField nonlinearStress;
00082 
00083 public:
00084 
00085     //- Runtime type information
00086     TypeName("NonlinearKEShih");
00087 
00088     // Constructors
00089 
00090         //- from components
00091         NonlinearKEShih
00092         (
00093             const volVectorField& U,
00094             const surfaceScalarField& phi,
00095             transportModel& transport
00096         );
00097 
00098 
00099     // Destructor
00100 
00101         ~NonlinearKEShih(){}
00102 
00103 
00104     // Member Functions
00105 
00106         tmp<volScalarField> nut() const
00107         {
00108             return nut_;
00109         }
00110 
00111         //- Return the effective diffusivity for k
00112         tmp<volScalarField> DkEff() const
00113         {
00114             return tmp<volScalarField> 
00115             (
00116                 new volScalarField("DkEff", alphak*nut_ + nu())
00117             );
00118         }
00119 
00120         //- Return the effective diffusivity for epsilon
00121         tmp<volScalarField> DepsilonEff() const
00122         {
00123             return tmp<volScalarField> 
00124             (
00125                 new volScalarField("DepsilonEff", alphaEps*nut_ + nu())
00126             );
00127         }
00128 
00129         tmp<volScalarField> k() const
00130         {
00131             return k_;
00132         }
00133 
00134         tmp<volScalarField> epsilon() const
00135         {
00136             return epsilon_;
00137         }
00138 
00139 
00140         tmp<volTensorField> R() const;
00141 
00142         tmp<fvVectorMatrix> divR(volVectorField& U) const;
00143 
00144         void correct();
00145 
00146         //- Read turbulenceProperties dictionary
00147         bool read();
00148 };
00149 
00150 
00151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00152 
00153 } // End namespace turbulenceModels
00154 } // End namespace Foam
00155 
00156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00157 
00158 #endif
00159 
00160 // ************************************************************************* //
For further information go to www.openfoam.org