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

Smagorinsky.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::incompressible::LESModels::Smagorinsky
00026 
00027 Description
00028     The Isochoric Smagorinsky Model for incompressible flows.
00029 
00030     Algebraic eddy viscosity SGS model founded on the assumption that
00031     local equilibrium prevails.
00032     Thus,
00033     @verbatim
00034         B = 2/3*k*I - 2*nuSgs*dev(D)
00035         Beff = 2/3*k*I - 2*nuEff*dev(D)
00036 
00037     where
00038 
00039         D = symm(grad(U));
00040         k = (2*ck/ce)*delta^2*||D||^2
00041         nuSgs = ck*sqrt(k)*delta
00042         nuEff = nuSgs + nu
00043     @endverbatim
00044 
00045 SourceFiles
00046     Smagorinsky.C
00047 
00048 \*---------------------------------------------------------------------------*/
00049 
00050 #ifndef Smagorinsky_H
00051 #define Smagorinsky_H
00052 
00053 #include <incompressibleLESModels/GenEddyVisc.H>
00054 
00055 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00056 
00057 namespace Foam
00058 {
00059 namespace incompressible
00060 {
00061 namespace LESModels
00062 {
00063 
00064 /*---------------------------------------------------------------------------*\
00065                            Class Smagorinsky Declaration
00066 \*---------------------------------------------------------------------------*/
00067 
00068 class Smagorinsky
00069 :
00070     public GenEddyVisc
00071 {
00072     // Private data
00073 
00074         dimensionedScalar ck_;
00075 
00076 
00077     // Private Member Functions
00078 
00079         //- Update sub-grid scale fields
00080         void updateSubGridScaleFields(const volTensorField& gradU);
00081 
00082         // Disallow default bitwise copy construct and assignment
00083         Smagorinsky(const Smagorinsky&);
00084         Smagorinsky& operator=(const Smagorinsky&);
00085 
00086 
00087 public:
00088 
00089     //- Runtime type information
00090     TypeName("Smagorinsky");
00091 
00092     // Constructors
00093 
00094         //- Construct from components
00095         Smagorinsky
00096         (
00097             const volVectorField& U,
00098             const surfaceScalarField& phi,
00099             transportModel& transport
00100         );
00101 
00102 
00103     //- Destructor
00104     virtual ~Smagorinsky()
00105     {}
00106 
00107 
00108     // Member Functions
00109 
00110         //- Return SGS kinetic energy
00111         //  calculated from the given velocity gradient
00112         tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
00113         {
00114             return (2.0*ck_/ce_)*sqr(delta())*magSqr(dev(symm(gradU)));
00115         }
00116 
00117         //- Return SGS kinetic energy
00118         virtual tmp<volScalarField> k() const
00119         {
00120             return k(fvc::grad(U()));
00121         }
00122 
00123 
00124         //- Correct Eddy-Viscosity and related properties
00125         virtual void correct(const tmp<volTensorField>& gradU);
00126 
00127         //- Read LESProperties dictionary
00128         virtual bool read();
00129 };
00130 
00131 
00132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00133 
00134 } // End namespace LESModels
00135 } // End namespace incompressible
00136 } // End namespace Foam
00137 
00138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00139 
00140 #endif
00141 
00142 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines