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

homogeneousDynSmagorinsky.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-2011 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::incompressible::LESModels::homogeneousDynSmagorinsky
00026 
00027 Description
00028     The Isochoric homogeneous dynamic Smagorinsky Model for
00029     incompressible flows.
00030 
00031     Algebraic eddy viscosity SGS model founded on the assumption that
00032     local equilibrium prevails.
00033     Thus,
00034     @verbatim
00035         B = 2/3*k*I - 2*nuSgs*dev(D)
00036         Beff = 2/3*k*I - 2*nuEff*dev(D)
00037 
00038     where
00039 
00040         k = cI*delta^2*||D||^2
00041         nuSgs = cD*delta^2*||D||
00042         nuEff = nuSgs + nu
00043 
00044     In the dynamic version of the choric  Smagorinsky model
00045     the coefficients cI and cD are calculated during the simulation,
00046 
00047         cI=<K*m>/<m*m>
00048 
00049     and
00050 
00051         cD=<L.M>/<M.M>,
00052 
00053     where
00054 
00055         K = 0.5*(F(U.U) - F(U).F(U))
00056         m = delta^2*(4*||F(D)||^2 - F(||D||^2))
00057         L = dev(F(U*U) - F(U)*F(U))
00058         M = delta^2*(F(||D||*dev(D)) - 4*||F(D)||*F(dev(D)))
00059 
00060     The averaging <...> is over the whole domain, i.e. homogeneous turbulence
00061     is assumed.
00062     @endverbatim
00063 
00064 SourceFiles
00065     homogeneousDynSmagorinsky.C
00066 
00067 \*---------------------------------------------------------------------------*/
00068 
00069 #ifndef homogeneousDynSmagorinsky_H
00070 #define homogeneousDynSmagorinsky_H
00071 
00072 #include <incompressibleLESModels/Smagorinsky.H>
00073 #include <LESfilters/LESfilter.H>
00074 
00075 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00076 
00077 namespace Foam
00078 {
00079 namespace incompressible
00080 {
00081 namespace LESModels
00082 {
00083 
00084 /*---------------------------------------------------------------------------*\
00085                            Class homogeneousDynSmagorinsky Declaration
00086 \*---------------------------------------------------------------------------*/
00087 
00088 class homogeneousDynSmagorinsky
00089 :
00090     public GenEddyVisc
00091 {
00092     // Private data
00093 
00094         volScalarField k_;
00095 
00096         autoPtr<LESfilter> filterPtr_;
00097         LESfilter& filter_;
00098 
00099 
00100     // Private Member Functions
00101 
00102         //- Update sub-grid scale fields
00103         void updateSubGridScaleFields(const volSymmTensorField& D);
00104 
00105         //- Calculate coefficients cD, cI from filtering velocity field
00106         dimensionedScalar cD(const volSymmTensorField& D) const;
00107         dimensionedScalar cI(const volSymmTensorField& D) const;
00108 
00109         // Disallow default bitwise copy construct and assignment
00110         homogeneousDynSmagorinsky(const homogeneousDynSmagorinsky&);
00111         homogeneousDynSmagorinsky& operator=(const homogeneousDynSmagorinsky&);
00112 
00113 
00114 public:
00115 
00116     //- Runtime type information
00117     TypeName("homogeneousDynSmagorinsky");
00118 
00119     // Constructors
00120 
00121         //- Construct from components
00122         homogeneousDynSmagorinsky
00123         (
00124             const volVectorField& U,
00125             const surfaceScalarField& phi,
00126             transportModel& transport
00127         );
00128 
00129 
00130     //- Destructor
00131     virtual ~homogeneousDynSmagorinsky()
00132     {}
00133 
00134 
00135     // Member Functions
00136 
00137         //- Return SGS kinetic energy
00138         virtual tmp<volScalarField> k() const
00139         {
00140             return k_;
00141         }
00142 
00143         //- Correct Eddy-Viscosity and related properties
00144         virtual void correct(const tmp<volTensorField>& gradU);
00145 
00146         //- Read LESProperties dictionary
00147         virtual bool read();
00148 };
00149 
00150 
00151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00152 
00153 } // End namespace LESModels
00154 } // End namespace incompressible
00155 } // End namespace Foam
00156 
00157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00158 
00159 #endif
00160 
00161 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines