OpenFOAM logo
Open Source CFD Toolkit

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-2005 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 the
00013     Free Software Foundation; either version 2 of the License, or (at your
00014     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, write to the Free Software Foundation,
00023     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
00024 
00025 Class
00026     Smagorinsky
00027 
00028 Description
00029 <pre>
00030     The Isochoric Smagorinsky Model
00031     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00032     Algebraic eddy viscosity SGS model founded on the assumption that
00033     local equilibrium prevails, hence
00034 
00035         B = 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 </pre>
00044 
00045 SourceFiles
00046     Smagorinsky.C
00047 
00048 \*---------------------------------------------------------------------------*/
00049 
00050 #ifndef Smagorinsky_H
00051 #define Smagorinsky_H
00052 
00053 #include "GenEddyVisc.H"
00054 
00055 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00056 
00057 namespace Foam
00058 {
00059 namespace LESmodels
00060 {
00061 
00062 /*---------------------------------------------------------------------------*\
00063                            Class Smagorinsky Declaration
00064 \*---------------------------------------------------------------------------*/
00065 
00066 class Smagorinsky
00067 :
00068     public GenEddyVisc
00069 {
00070     // Private data
00071 
00072         dimensionedScalar ck_;
00073 
00074 
00075     // Private Member Functions
00076 
00077         // Disallow default bitwise copy construct and assignment
00078         Smagorinsky(const Smagorinsky&);
00079         Smagorinsky& operator=(const Smagorinsky&);
00080 
00081 
00082 public:
00083 
00084     //- Runtime type information
00085     TypeName("Smagorinsky");
00086 
00087     // Constructors
00088 
00089         //- Construct from components
00090         Smagorinsky
00091         (
00092             const volVectorField& U,
00093             const surfaceScalarField& phi,
00094             transportModel& transport
00095         );
00096 
00097 
00098     // Destructor
00099 
00100         ~Smagorinsky()
00101         {}
00102 
00103 
00104     // Member Functions
00105 
00106         //- Return SGS kinetic energy
00107         //  calculated from the given velocity gradient
00108         tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
00109         {
00110             return (2.0*ck_/ce_)*sqr(delta())*magSqr(dev(symm(gradU)));
00111         }
00112 
00113         //- Return SGS kinetic energy
00114         tmp<volScalarField> k() const
00115         {
00116             return k(fvc::grad(U()));
00117         }
00118 
00119 
00120         //- Correct Eddy-Viscosity and related properties
00121         virtual void correct(const tmp<volTensorField>& gradU);
00122 
00123         //- Read turbulenceProperties dictionary
00124         virtual bool read();
00125 };
00126 
00127 
00128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00129 
00130 } // End namespace LESmodels
00131 } // End namespace Foam
00132 
00133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00134 
00135 #endif
00136 
00137 // ************************************************************************* //
For further information go to www.openfoam.org