Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / dynMixedSmagorinsky / dynMixedSmagorinsky.H
blob7367c917dd73a0ea94f113bb7aac3bcd23685dd3
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::dynMixedSmagorinsky
27 Description
28     The Mixed Isochoric Smagorinsky Model for incompressible flows.
30     The mixed model is a linear combination of an eddy viscosity model
31     with a scale similarity model.
32     @verbatim
33           B = (L + C) + R = (F(v*v) - F(v)*F(v)) + R
34     @endverbatim
36     The algebraic eddy viscosity SGS model is founded on the assumption
37     that local equilibrium prevails, hence
38     @verbatim
39         R = 2/3*rho*k*I - 2*nuEff*dev(D)
40     where
41         k = cI*delta^2*||D||^2
42         nuEff = ck*sqrt(k)*delta + nu
43     @endverbatim
45     The Leonard and cross contributions are incorporated
46     by adding,
47     @verbatim
48         + div(((filter(U*U) - filter(U)*filter(U)) -
49           0.333*I*tr(filter(U*U) - filter(U)*filter(U))))
50         + div((filter(U*epsilon) - filter(U)*filter(epsilon)))
51     @endverbatim
52     to the rhs. of the equations.  This version implements filtering to
53     evaluate the coefficients in the model.
55 SourceFiles
56     dynMixedSmagorinsky.C
58 \*---------------------------------------------------------------------------*/
60 #ifndef dynMixedSmagorinsky_H
61 #define dynMixedSmagorinsky_H
63 #include "dynSmagorinsky.H"
64 #include "scaleSimilarity.H"
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 namespace Foam
70 namespace incompressible
72 namespace LESModels
75 /*---------------------------------------------------------------------------*\
76                            Class dynMixedSmagorinsky Declaration
77 \*---------------------------------------------------------------------------*/
79 class dynMixedSmagorinsky
81     public scaleSimilarity,
82     public dynSmagorinsky
84     // Private Member Functions
86         // Disallow default bitwise copy construct and assignment
87         dynMixedSmagorinsky(const dynMixedSmagorinsky&);
88         dynMixedSmagorinsky& operator=(const dynMixedSmagorinsky&);
90 public:
92     //- Runtime type information
93     TypeName("dynMixedSmagorinsky");
95     // Constructors
97         //- Constructors from components
98         dynMixedSmagorinsky
99         (
100             const volVectorField& U,
101             const surfaceScalarField& phi,
102             transportModel& transport
103         );
106     // Destructor
108         ~dynMixedSmagorinsky()
109         {}
112     // Member Functions
114         //- Return SGS kinetic energy
115         tmp<volScalarField> k() const;
117         //- Return sub-grid disipation rate
118         tmp<volScalarField> epsilon() const;
120         //- Return the sub-grid stress tensor.
121         tmp<volSymmTensorField> B() const;
123         //- Return the effective sub-grid turbulence stress tensor
124         //  including the laminar stress
125         tmp<volSymmTensorField> devBeff() const;
127         //- Returns div(B).
128         // This is the additional term due to the filtering of the NSE.
129         tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
131         //- Correct Eddy-Viscosity and related properties
132         void correct(const tmp<volTensorField>& gradU);
134         //- Read LESProperties dictionary
135         bool read();
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 } // End namespace LESModels
142 } // End namespace incompressible
143 } // End namespace Foam
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 #endif
149 // ************************************************************************* //