initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / dynMixedSmagorinsky / dynMixedSmagorinsky.H
bloba69f14e6ba428beea8b1250efc42f3b51802cacb
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM 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 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::incompressible::LESModels::dynMixedSmagorinsky
28 Description
29     The Mixed Isochoric Smagorinsky Model for incompressible flows.
31     The mixed model is a linear combination of an eddy viscosity model
32     with a scale similarity model.
33     @verbatim
34           B = (L + C) + R = (F(v*v) - F(v)*F(v)) + R
35     @endverbatim
37     The algebraic eddy viscosity SGS model is founded on the assumption
38     that local equilibrium prevails, hence
39     @verbatim
40         R = 2/3*rho*k*I - 2*nuEff*dev(D)
41     where
42         k = cI*delta^2*||D||^2
43         nuEff = ck*sqrt(k)*delta + nu
44     @endverbatim
46     The Leonard and cross contributions are incorporated
47     by adding,
48     @verbatim
49         + div(((filter(U*U) - filter(U)*filter(U)) -
50           0.333*I*tr(filter(U*U) - filter(U)*filter(U))))
51         + div((filter(U*epsilon) - filter(U)*filter(epsilon)))
52     @endverbatim
53     to the rhs. of the equations.  This version implements filtering to
54     evaluate the coefficients in the model.
56 SourceFiles
57     dynMixedSmagorinsky.C
59 \*---------------------------------------------------------------------------*/
61 #ifndef dynMixedSmagorinsky_H
62 #define dynMixedSmagorinsky_H
64 #include "dynSmagorinsky.H"
65 #include "scaleSimilarity.H"
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 namespace Foam
71 namespace incompressible
73 namespace LESModels
76 /*---------------------------------------------------------------------------*\
77                            Class dynMixedSmagorinsky Declaration
78 \*---------------------------------------------------------------------------*/
80 class dynMixedSmagorinsky
82     public scaleSimilarity,
83     public dynSmagorinsky
85     // Private Member Functions
87         // Disallow default bitwise copy construct and assignment
88         dynMixedSmagorinsky(const dynMixedSmagorinsky&);
89         dynMixedSmagorinsky& operator=(const dynMixedSmagorinsky&);
91 public:
93     //- Runtime type information
94     TypeName("dynMixedSmagorinsky");
96     // Constructors
98         //- Constructors from components
99         dynMixedSmagorinsky
100         (
101             const volVectorField& U,
102             const surfaceScalarField& phi,
103             transportModel& transport
104         );
107     // Destructor
109         ~dynMixedSmagorinsky()
110         {}
113     // Member Functions
115         //- Return SGS kinetic energy
116         tmp<volScalarField> k() const;
118         //- Return sub-grid disipation rate
119         tmp<volScalarField> epsilon() const;
121         //- Return the sub-grid stress tensor.
122         tmp<volSymmTensorField> B() const;
124         //- Return the effective sub-grid turbulence stress tensor
125         //  including the laminar stress
126         tmp<volSymmTensorField> devBeff() const;
128         //- Returns div(B).
129         // This is the additional term due to the filtering of the NSE.
130         tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
132         //- Correct Eddy-Viscosity and related properties
133         void correct(const tmp<volTensorField>& gradU);
135         //- Read LESProperties dictionary
136         bool read();
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 } // End namespace LESModels
143 } // End namespace incompressible
144 } // End namespace Foam
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 #endif
150 // ************************************************************************* //