BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / dynLagrangian / dynLagrangian.H
blob2d83299d4ac401d2d12df74640071d9f4ca52509
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::dynLagrangian
27 Description
28     Dynamic eddy-viscosity model with Lagrangian averaging for incompressible
29     flow
31     \verbatim
32         B = 2/3*k*I - 2*nuSgs*dev(D)
33         Beff = 2/3*k*I - 2*nuEff*dev(D)
35     where
37         D = symm(grad(U))
38         nuSgs = (flm/fmm)*delta^2*sqrt(2)*|D|
39         nuEff = nuSgs + nu
41         Two relaxation equations are used to evaluate flm and fmm:
43         d/dt(flm) + div(U*flm)
44         =
45         (1/T)*(L && M - flm)
47         d/dt(fmm) + div(U*fmm)
48         =
49         (1/T)*(M && M - flm)
51     where
53         L = F(U.U) - F(U).F(U)
54         M = 2.0 delta^2 (F(|D|.D) - 4 F(|D|).F(D))
55         T = 1.5*delta*(flm.fmm)^(-1.0/8.0)
57     \endverbatim
59     Reference:
60     \verbatim
61         "A Lagrangian dynamic subgrid-scale model of turbulence"
62         Charles Meneveau,
63         Thomas Lund,
64         William Cabot,
65         J. Fluid Mech (1996), vol 319, pp. 353-385
66     \endverbatim
68 SourceFiles
69     dynLagrangian.C
71 \*---------------------------------------------------------------------------*/
73 #ifndef dynLagrangian_H
74 #define dynLagrangian_H
76 #include "GenEddyVisc.H"
77 #include "simpleFilter.H"
78 #include "LESfilter.H"
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 namespace Foam
84 namespace incompressible
86 namespace LESModels
89 /*---------------------------------------------------------------------------*\
90                        Class dynLagrangian Declaration
91 \*---------------------------------------------------------------------------*/
93 class dynLagrangian
95     public GenEddyVisc
97     // Private data
99         volScalarField flm_;
100         volScalarField fmm_;
102         dimensionedScalar theta_;
104         simpleFilter simpleFilter_;
105         autoPtr<LESfilter> filterPtr_;
106         LESfilter& filter_;
108         dimensionedScalar flm0_;
109         dimensionedScalar fmm0_;
112     // Private Member Functions
114         //- Update sub-grid scale fields
115         void updateSubGridScaleFields
116         (
117             const tmp<volTensorField>& gradU
118         );
120         // Disallow default bitwise copy construct and assignment
121         dynLagrangian(const dynLagrangian&);
122         dynLagrangian& operator=(const dynLagrangian&);
125 public:
127     //- Runtime type information
128     TypeName("dynLagrangian");
130     // Constructors
132         //- Construct from components
133         dynLagrangian
134         (
135             const volVectorField& U,
136             const surfaceScalarField& phi,
137             transportModel& transport,
138             const word& turbulenceModelName = turbulenceModel::typeName,
139             const word& modelName = typeName
140         );
143     //- Destructor
144     virtual ~dynLagrangian()
145     {}
148     // Member Functions
150         //- Return SGS kinetic energy
151         tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
152         {
153             return 2.0*sqr(delta())*magSqr(dev(symm(gradU)));
154         }
156         //- Return SGS kinetic energy
157         virtual tmp<volScalarField> k() const
158         {
159             return k(fvc::grad(U()));
160         }
162         //- Return the effective diffusivity for k
163         tmp<volScalarField> DkEff() const
164         {
165             return tmp<volScalarField>
166             (
167                 new volScalarField("DkEff", nuSgs_ + nu())
168             );
169         }
171         //- Correct Eddy-Viscosity and related properties
172         virtual void correct(const tmp<volTensorField>& gradU);
175         //- Read LESProperties dictionary
176         virtual bool read();
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 } // End namespace LESModels
183 } // End namespace incompressible
184 } // End namespace Foam
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 #endif
190 // ************************************************************************* //