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 / mixedSmagorinsky / mixedSmagorinsky.H
blobe14d74271d4c920cb613e18db9a12ec1ca1bf1dd
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::mixedSmagorinsky
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     (Smagorinsky) with a scale similarity model.  Hence
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*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.
54 SourceFiles
55     mixedSmagorinsky.C
57 \*---------------------------------------------------------------------------*/
59 #ifndef mixedSmagorinsky_H
60 #define mixedSmagorinsky_H
62 #include "scaleSimilarity.H"
63 #include "Smagorinsky.H"
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 namespace Foam
69 namespace incompressible
71 namespace LESModels
74 /*---------------------------------------------------------------------------*\
75                            Class mixedSmagorinsky Declaration
76 \*---------------------------------------------------------------------------*/
78 class mixedSmagorinsky
80     public scaleSimilarity,
81     public Smagorinsky
83     // Private Member Functions
85         // Disallow default bitwise copy construct and assignment
86         mixedSmagorinsky(const mixedSmagorinsky&);
87         mixedSmagorinsky& operator=(const mixedSmagorinsky&);
89 public:
91     //- Runtime type information
92     TypeName("mixedSmagorinsky");
95     // Constructors
97         //- Construct from components
98         mixedSmagorinsky
99         (
100             const volVectorField& U,
101             const surfaceScalarField& phi,
102             transportModel& transport
103         );
106     //- Destructor
107     virtual ~mixedSmagorinsky()
108     {}
111     // Member Functions
113         //- Return the SGS turbulent kinetic energy.
114         virtual tmp<volScalarField> k() const;
116         //- Return the SGS turbulent disipation rate.
117         virtual tmp<volScalarField> epsilon() const;
119         //- Return the SGS viscosity.
120         virtual tmp<volScalarField> nuSgs() const
121         {
122             return nuSgs_;
123         }
125         //- Return the sub-grid stress tensor.
126         virtual tmp<volSymmTensorField> B() const;
128         //- Return the effective sub-grid turbulence stress tensor
129         //  including the laminar stress
130         virtual tmp<volSymmTensorField> devBeff() const;
132         //- Implementation of div(B). This is necessary to override
133         // (and include) the div(B) terms from both the parent classes.
134         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
136         //- Correct Eddy-Viscosity and related properties
137         virtual void correct(const tmp<volTensorField>& gradU);
139         //- Read LESProperties dictionary
140         virtual bool read();
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 } // End namespace LESModels
147 } // End namespace incompressible
148 } // End namespace Foam
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 #endif
154 // ************************************************************************* //