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

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