FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

PDRkEpsilon.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-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::compressible::RASModels::PDRkEpsilon
00026 
00027 Description
00028     Standard k-epsilon turbulence model with additional source terms
00029     corresponding to PDR basic drag model (\link basic.H \endlink)
00030 
00031     The default model coefficients correspond to the following:
00032     @verbatim
00033         kEpsilonCoeffs
00034         {
00035             Cmu         0.09;
00036             C1          1.44;
00037             C2          1.92;
00038             C3          -0.33;  // only for compressible
00039             sigmak      1.0;    // only for compressible
00040             sigmaEps    1.3;
00041             Prt         1.0;    // only for compressible
00042         }
00043     @endverbatim
00044 
00045     The turbulence source term \f$ G_{R} \f$ appears in the
00046     \f$ \kappa-\epsilon \f$ equation for the generation of turbulence due to
00047     interaction with unresolved obstacles.
00048 
00049     In the \f$ \epsilon  \f$ equation \f$ C_{1} G_{R} \f$ is added as a source
00050     term.
00051 
00052     In the \f$ \kappa \f$ equation \f$ G_{R} \f$ is added as a source term.
00053 
00054 SourceFiles
00055     PDRkEpsilon.C
00056     PDRkEpsilonCorrect.C
00057 
00058 \*---------------------------------------------------------------------------*/
00059 
00060 #ifndef compressiblePDRkEpsilon_H
00061 #define compressiblePDRkEpsilon_H
00062 
00063 #include <compressibleRASModels/RASModel.H>
00064 
00065 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00066 
00067 namespace Foam
00068 {
00069 namespace compressible
00070 {
00071 namespace RASModels
00072 {
00073 
00074 /*---------------------------------------------------------------------------*\
00075                            Class PDRkEpsilon Declaration
00076 \*---------------------------------------------------------------------------*/
00077 
00078 class PDRkEpsilon
00079 :
00080     public RASModel
00081 {
00082     // Private data
00083 
00084     // Model coefficients
00085 
00086         dimensionedScalar Cmu_;
00087         dimensionedScalar C1_;
00088         dimensionedScalar C2_;
00089         dimensionedScalar sigmak_;
00090         dimensionedScalar sigmaEps_;
00091         dimensionedScalar Prt_;
00092 
00093     // Fields
00094 
00095         volScalarField k_;
00096         volScalarField epsilon_;
00097         volScalarField mut_;
00098         volScalarField alphat_;
00099 
00100 
00101 public:
00102 
00103     //- Runtime type information
00104     TypeName("PDRkEpsilon");
00105 
00106 
00107     // Constructors
00108 
00109         //- Construct from components
00110         PDRkEpsilon
00111         (
00112             const volScalarField& rho,
00113             const volVectorField& U,
00114             const surfaceScalarField& phi,
00115             const basicThermo& thermophysicalModel
00116         );
00117 
00118 
00119     //- Destructor
00120     virtual ~PDRkEpsilon()
00121     {}
00122 
00123 
00124     // Member Functions
00125 
00126         tmp<volScalarField> mut() const
00127         {
00128             return mut_;
00129         }
00130 
00131         //- Return the effective diffusivity for k
00132         tmp<volScalarField> DkEff() const
00133         {
00134             return tmp<volScalarField>
00135             (
00136                 new volScalarField("DkEff", mut_/sigmak_ + mu())
00137             );
00138         }
00139 
00140         //- Return the effective diffusivity for epsilon
00141         tmp<volScalarField> DepsilonEff() const
00142         {
00143             return tmp<volScalarField>
00144             (
00145                 new volScalarField("DepsilonEff", mut_/sigmaEps_ + mu())
00146             );
00147         }
00148 
00149         //- Return the effective turbulent thermal diffusivity
00150         tmp<volScalarField> alphaEff() const
00151         {
00152             return tmp<volScalarField>
00153             (
00154                 new volScalarField("alphaEff", alphat_ + alpha())
00155             );
00156         }
00157 
00158         //- Return the turbulence kinetic energy
00159         tmp<volScalarField> k() const
00160         {
00161             return k_;
00162         }
00163 
00164         //- Return the turbulence kinetic energy dissipation rate
00165         tmp<volScalarField> epsilon() const
00166         {
00167             return epsilon_;
00168         }
00169 
00170         //- Return the Reynolds stress tensor
00171         tmp<volSymmTensorField> R() const;
00172 
00173         //- Return the effective stress tensor including the laminar stress
00174         tmp<volSymmTensorField> devRhoReff() const;
00175 
00176         //- Return the source term for the momentum equation
00177         tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const;
00178 
00179         //- Solve the turbulence equations and correct the turbulence viscosity
00180         void correct();
00181 
00182         //- Read turbulenceProperties dictionary
00183         bool read();
00184 };
00185 
00186 
00187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00188 
00189 } // End namespace RASModels
00190 } // End namespace compressible
00191 } // End namespace Foam
00192 
00193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00194 
00195 #endif
00196 
00197 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines