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

lowReOneEqEddy.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::lowReOneEqEddy
00026 
00027 Description
00028     One Equation Eddy Viscosity Model for compressible flow
00029 
00030     @verbatim
00031         d/dt(rho*k) + div(rho*U*k) - div(muEff*grad(k))
00032         =
00033         -rho*B*L - ce*rho*k^3/2/delta
00034 
00035     and
00036 
00037         B = 2/3*k*I - 2*nuSgs*dev(D)
00038 
00039     where
00040 
00041         nuSgsHiRe = ck*sqrt(k)*delta
00042         nuSgs = (nu/beta)*(1 - exp(-beta*nuSgsHiRe/nu));
00043     @endverbatim
00044 
00045 SourceFiles
00046     lowReOneEqEddy.C
00047 
00048 \*---------------------------------------------------------------------------*/
00049 
00050 #ifndef compressibleLowReOneEqEddy_H
00051 #define compressibleLowReOneEqEddy_H
00052 
00053 #include <compressibleLESModels/GenEddyVisc.H>
00054 
00055 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00056 
00057 namespace Foam
00058 {
00059 namespace compressible
00060 {
00061 namespace LESModels
00062 {
00063 
00064 /*---------------------------------------------------------------------------*\
00065                            Class lowReOneEqEddy Declaration
00066 \*---------------------------------------------------------------------------*/
00067 
00068 class lowReOneEqEddy
00069 :
00070     public GenEddyVisc
00071 {
00072     // Private data
00073 
00074         dimensionedScalar ck_;
00075         dimensionedScalar beta_;
00076 
00077     // Private Member Functions
00078 
00079         //- Update sub-grid scale fields
00080         void updateSubGridScaleFields();
00081 
00082         // Disallow default bitwise copy construct and assignment
00083         lowReOneEqEddy(const lowReOneEqEddy&);
00084         lowReOneEqEddy& operator=(const lowReOneEqEddy&);
00085 
00086 
00087 public:
00088 
00089     //- Runtime type information
00090     TypeName("lowReOneEqEddy");
00091 
00092 
00093     // Constructors
00094 
00095         //- Constructor from components
00096         lowReOneEqEddy
00097         (
00098             const volScalarField& rho,
00099             const volVectorField& U,
00100             const surfaceScalarField& phi,
00101             const basicThermo& thermoPhysicalModel
00102         );
00103 
00104 
00105     //- Destructor
00106     virtual ~lowReOneEqEddy()
00107     {}
00108 
00109 
00110     // Member Functions
00111 
00112         //- Return the effective diffusivity for k
00113         tmp<volScalarField> DkEff() const
00114         {
00115             return tmp<volScalarField>
00116             (
00117                 new volScalarField("DkEff", muSgs_ + mu())
00118             );
00119         }
00120 
00121         //- Correct Eddy-Viscosity and related properties
00122         virtual void correct(const tmp<volTensorField>& gradU);
00123 
00124         //- Read LESProperties dictionary
00125         virtual bool read();
00126 };
00127 
00128 
00129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00130 
00131 } // End namespace LESModels
00132 } // End namespace compressible
00133 } // End namespace Foam
00134 
00135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00136 
00137 #endif
00138 
00139 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines