Merge remote-tracking branch 'origin/nr/multiSolverFix' into nextRelease
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / mixedSmagorinsky / mixedSmagorinsky.H
blob94052cc7fb5ba7db0008f2858c840a6f02772ca6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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::mixedSmagorinsky
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     (Smagorinsky) with a scale similarity model.  Hence
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*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.
55 SourceFiles
56     mixedSmagorinsky.C
58 \*---------------------------------------------------------------------------*/
60 #ifndef mixedSmagorinsky_H
61 #define mixedSmagorinsky_H
63 #include "scaleSimilarity.H"
64 #include "Smagorinsky.H"
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 namespace Foam
70 namespace incompressible
72 namespace LESModels
75 /*---------------------------------------------------------------------------*\
76                            Class mixedSmagorinsky Declaration
77 \*---------------------------------------------------------------------------*/
79 class mixedSmagorinsky
81     public scaleSimilarity,
82     public Smagorinsky
84     // Private Member Functions
86         // Disallow default bitwise copy construct and assignment
87         mixedSmagorinsky(const mixedSmagorinsky&);
88         mixedSmagorinsky& operator=(const mixedSmagorinsky&);
90 public:
92     //- Runtime type information
93     TypeName("mixedSmagorinsky");
96     // Constructors
98         //- Construct from components
99         mixedSmagorinsky
100         (
101             const volVectorField& U,
102             const surfaceScalarField& phi,
103             transportModel& transport
104         );
107     //- Destructor
108     virtual ~mixedSmagorinsky()
109     {}
112     // Member Functions
114         //- Return the SGS turbulent kinetic energy.
115         virtual tmp<volScalarField> k() const;
117         //- Return the SGS turbulent disipation rate.
118         virtual tmp<volScalarField> epsilon() const;
120         //- Return the SGS viscosity.
121         virtual tmp<volScalarField> nuSgs() const
122         {
123             return nuSgs_;
124         }
126         //- Return the sub-grid stress tensor.
127         virtual tmp<volSymmTensorField> B() const;
129         //- Return the effective sub-grid turbulence stress tensor
130         //  including the laminar stress
131         virtual tmp<volSymmTensorField> devBeff() const;
133         //- Implementation of div(B). This is necessary to override
134         // (and include) the div(B) terms from both the parent classes.
135         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
137         //- Correct Eddy-Viscosity and related properties
138         virtual void correct(const tmp<volTensorField>& gradU);
140         //- Read LESProperties dictionary
141         virtual bool read();
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 } // End namespace LESModels
148 } // End namespace incompressible
149 } // End namespace Foam
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 #endif
155 // ************************************************************************* //