Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / mixedSmagorinsky / mixedSmagorinsky.H
blobef37bbaa2b616e147809b90097759675ff44f9ac
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2011 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, 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             const word& turbulenceModelName = turbulenceModel::typeName,
104             const word& modelName = typeName
105         );
108     //- Destructor
109     virtual ~mixedSmagorinsky()
110     {}
113     // Member Functions
115         //- Return the SGS turbulent kinetic energy.
116         virtual tmp<volScalarField> k() const;
118         //- Return the SGS turbulent disipation rate.
119         virtual tmp<volScalarField> epsilon() const;
121         //- Return the SGS viscosity.
122         virtual tmp<volScalarField> nuSgs() const
123         {
124             return nuSgs_;
125         }
127         //- Return the sub-grid stress tensor.
128         virtual tmp<volSymmTensorField> B() const;
130         //- Return the effective sub-grid turbulence stress tensor
131         //  including the laminar stress
132         virtual tmp<volSymmTensorField> devBeff() const;
134         //- Implementation of div(B). This is necessary to override
135         // (and include) the div(B) terms from both the parent classes.
136         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
138         //- Correct Eddy-Viscosity and related properties
139         virtual void correct(const tmp<volTensorField>& gradU);
141         //- Read LESProperties dictionary
142         virtual bool read();
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 } // End namespace LESModels
149 } // End namespace incompressible
150 } // End namespace Foam
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 #endif
156 // ************************************************************************* //