Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / Smagorinsky / Smagorinsky.H
blob3a7c44d7109281861a5d1dfe9781bc198aeeb25d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::incompressible::LESModels::Smagorinsky
27 Description
28     The Isochoric Smagorinsky Model for incompressible flows.
30     Algebraic eddy viscosity SGS model founded on the assumption that
31     local equilibrium prevails.
32     Thus,
33     @verbatim
34         B = 2/3*k*I - 2*nuSgs*dev(D)
35         Beff = 2/3*k*I - 2*nuEff*dev(D)
37     where
39         D = symm(grad(U));
40         k = (2*ck/ce)*delta^2*||D||^2
41         nuSgs = ck*sqrt(k)*delta
42         nuEff = nuSgs + nu
43     @endverbatim
45 SourceFiles
46     Smagorinsky.C
48 \*---------------------------------------------------------------------------*/
50 #ifndef Smagorinsky_H
51 #define Smagorinsky_H
53 #include "GenEddyVisc.H"
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 namespace Foam
59 namespace incompressible
61 namespace LESModels
64 /*---------------------------------------------------------------------------*\
65                            Class Smagorinsky Declaration
66 \*---------------------------------------------------------------------------*/
68 class Smagorinsky
70     public GenEddyVisc
72     // Private data
74         dimensionedScalar ck_;
77     // Private Member Functions
79         //- Update sub-grid scale fields
80         void updateSubGridScaleFields(const volTensorField& gradU);
82         // Disallow default bitwise copy construct and assignment
83         Smagorinsky(const Smagorinsky&);
84         Smagorinsky& operator=(const Smagorinsky&);
87 public:
89     //- Runtime type information
90     TypeName("Smagorinsky");
92     // Constructors
94         //- Construct from components
95         Smagorinsky
96         (
97             const volVectorField& U,
98             const surfaceScalarField& phi,
99             transportModel& transport
100         );
103     //- Destructor
104     virtual ~Smagorinsky()
105     {}
108     // Member Functions
110         //- Return SGS kinetic energy
111         //  calculated from the given velocity gradient
112         tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
113         {
114             return (2.0*ck_/ce_)*sqr(delta())*magSqr(dev(symm(gradU)));
115         }
117         //- Return SGS kinetic energy
118         virtual tmp<volScalarField> k() const
119         {
120             return k(fvc::grad(U()));
121         }
124         //- Correct Eddy-Viscosity and related properties
125         virtual void correct(const tmp<volTensorField>& gradU);
127         //- Read LESProperties dictionary
128         virtual bool read();
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 } // End namespace LESModels
135 } // End namespace incompressible
136 } // End namespace Foam
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 #endif
142 // ************************************************************************* //