BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / homogeneousDynSmagorinsky / homogeneousDynSmagorinsky.H
blob05ac587a5b6d94efde10ae55cadb035356105484
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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::homogeneousDynSmagorinsky
27 Description
28     The Isochoric homogeneous dynamic Smagorinsky Model for
29     incompressible flows.
31     Algebraic eddy viscosity SGS model founded on the assumption that
32     local equilibrium prevails.
33     Thus,
34     \verbatim
35         B = 2/3*k*I - 2*nuSgs*dev(D)
36         Beff = 2/3*k*I - 2*nuEff*dev(D)
38     where
40         k = cI*delta^2*||D||^2
41         nuSgs = cD*delta^2*||D||
42         nuEff = nuSgs + nu
44     In the dynamic version of the choric  Smagorinsky model
45     the coefficients cI and cD are calculated during the simulation,
47         cI=<K*m>/<m*m>
49     and
51         cD=1/2*<L.M>/<M.M>,
53     where
55         K = 0.5*(F(U.U) - F(U).F(U))
56         m = delta^2*(4*||F(D)||^2 - F(||D||^2))
57         L = dev(F(U*U) - F(U)*F(U))
58         M = delta^2*(F(||D||*dev(D)) - 4*||F(D)||*F(dev(D)))
60     The averaging <...> is over the whole domain, i.e. homogeneous turbulence
61     is assumed
62     \endverbatim
64 SourceFiles
65     homogeneousDynSmagorinsky.C
67 \*---------------------------------------------------------------------------*/
69 #ifndef homogeneousDynSmagorinsky_H
70 #define homogeneousDynSmagorinsky_H
72 #include "Smagorinsky.H"
73 #include "LESfilter.H"
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 namespace Foam
79 namespace incompressible
81 namespace LESModels
84 /*---------------------------------------------------------------------------*\
85                            Class homogeneousDynSmagorinsky Declaration
86 \*---------------------------------------------------------------------------*/
88 class homogeneousDynSmagorinsky
90     public GenEddyVisc
92     // Private data
94         volScalarField k_;
96         autoPtr<LESfilter> filterPtr_;
97         LESfilter& filter_;
100     // Private Member Functions
102         //- Update sub-grid scale fields
103         void updateSubGridScaleFields(const volSymmTensorField& D);
105         //- Calculate coefficients cD, cI from filtering velocity field
106         dimensionedScalar cD(const volSymmTensorField& D) const;
107         dimensionedScalar cI(const volSymmTensorField& D) const;
109         // Disallow default bitwise copy construct and assignment
110         homogeneousDynSmagorinsky(const homogeneousDynSmagorinsky&);
111         homogeneousDynSmagorinsky& operator=(const homogeneousDynSmagorinsky&);
114 public:
116     //- Runtime type information
117     TypeName("homogeneousDynSmagorinsky");
119     // Constructors
121         //- Construct from components
122         homogeneousDynSmagorinsky
123         (
124             const volVectorField& U,
125             const surfaceScalarField& phi,
126             transportModel& transport,
127             const word& turbulenceModelName = turbulenceModel::typeName,
128             const word& modelName = typeName
129         );
132     //- Destructor
133     virtual ~homogeneousDynSmagorinsky()
134     {}
137     // Member Functions
139         //- Return SGS kinetic energy
140         virtual tmp<volScalarField> k() const
141         {
142             return k_;
143         }
145         //- Correct Eddy-Viscosity and related properties
146         virtual void correct(const tmp<volTensorField>& gradU);
148         //- Read LESProperties dictionary
149         virtual bool read();
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 } // End namespace LESModels
156 } // End namespace incompressible
157 } // End namespace Foam
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 #endif
163 // ************************************************************************* //