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

GenSGSStress.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::LESModels::GenSGSStress
00026 
00027 Description
00028     General base class for all compressible models that directly
00029     solve for the SGS stress tensor B.
00030 
00031     Contains tensor fields B (the SGS stress tensor) as well as scalar
00032     fields for k (SGS turbulent energy) gamma (SGS viscosity) and epsilon
00033     (SGS dissipation).
00034 
00035 SourceFiles
00036     GenSGSStress.C
00037 
00038 \*---------------------------------------------------------------------------*/
00039 
00040 #ifndef compressibleGenSGSStress_H
00041 #define compressibleGenSGSStress_H
00042 
00043 #include <compressibleLESModels/LESModel.H>
00044 
00045 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00046 
00047 namespace Foam
00048 {
00049 namespace compressible
00050 {
00051 namespace LESModels
00052 {
00053 
00054 /*---------------------------------------------------------------------------*\
00055                            Class GenSGSStress Declaration
00056 \*---------------------------------------------------------------------------*/
00057 
00058 class GenSGSStress
00059 :
00060     virtual public LESModel
00061 {
00062     // Private Member Functions
00063 
00064         // Disallow default bitwise copy construct and assignment
00065         GenSGSStress(const GenSGSStress&);
00066         GenSGSStress& operator=(const GenSGSStress&);
00067 
00068 
00069 protected:
00070 
00071     // Model coefficients
00072 
00073         dimensionedScalar ce_;
00074         dimensionedScalar Prt_;
00075 
00076     // Fields
00077 
00078         volSymmTensorField B_;
00079         volScalarField muSgs_;
00080         volScalarField alphaSgs_;
00081 
00082 
00083 public:
00084 
00085     // Constructors
00086 
00087         //- Constructor from components
00088         GenSGSStress
00089         (
00090             const volScalarField& rho,
00091             const volVectorField& U,
00092             const surfaceScalarField& phi,
00093             const basicThermo& thermoPhysicalModel
00094         );
00095 
00096 
00097     //- Destructor
00098     virtual ~GenSGSStress()
00099     {}
00100 
00101 
00102     // Member Functions
00103 
00104         //- Return the SGS turbulent kinetic energy
00105         virtual tmp<volScalarField> k() const
00106         {
00107             return 0.5*tr(B_);
00108         }
00109 
00110         //- Return the SGS turbulent dissipation
00111         virtual tmp<volScalarField> epsilon() const
00112         {
00113             volScalarField K = k();
00114             return ce_*K*sqrt(K)/delta();
00115         }
00116 
00117         //- Return the SGS viscosity
00118         virtual tmp<volScalarField> muSgs() const
00119         {
00120             return muSgs_;
00121         }
00122 
00123         //- Return the SGS thermal diffusivity
00124         virtual tmp<volScalarField> alphaSgs() const
00125         {
00126             return alphaSgs_;
00127         }
00128 
00129         //- Return thermal conductivity
00130         virtual tmp<volScalarField> alphaEff() const
00131         {
00132             return tmp<volScalarField>
00133             (
00134                 new volScalarField("alphaEff", alphaSgs_ + alpha())
00135             );
00136         }
00137 
00138         //- Return the sub-grid stress tensor
00139         virtual tmp<volSymmTensorField> B() const
00140         {
00141             return B_;
00142         }
00143 
00144         //- Return the deviatoric part of the effective sub-grid
00145         //  turbulence stress tensor including the laminar stress
00146         virtual tmp<volSymmTensorField> devRhoBeff() const;
00147 
00148         //- Returns divergence of B : i.e. the additional term in the
00149         //  filtered NSE
00150         virtual tmp<fvVectorMatrix> divDevRhoBeff(volVectorField& U) const;
00151 
00152         //- Correct Eddy-Viscosity and related properties
00153         virtual void correct(const tmp<volTensorField>& gradU);
00154 
00155         //- Read LESProperties dictionary
00156         virtual bool read();
00157 };
00158 
00159 
00160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00161 
00162 } // End namespace LESModels
00163 } // End namespace compressible
00164 } // End namespace Foam
00165 
00166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00167 
00168 #endif
00169 
00170 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines