Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / dynOneEqEddy / dynOneEqEddy.H
blobdf8b3cc3050cf0182fe6f191b63656dcd2eaed68
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::dynOneEqEddy
27 Description
28     One Equation Eddy Viscosity Model for incompressible flows.
30     Eddy viscosity SGS model using a modeled balance equation to simulate
31     the behaviour of k.
33     Thus
34     @verbatim
35         d/dt(k) + div(U*k) - div(nuSgs*grad(k))
36         =
37         -B*L - ce*k^3/2/delta
39     and
41         B = 2/3*k*I - 2*nuSgs*dev(D)
42         Beff = 2/3*k*I - 2*nuEff*dev(D)
44     where
46         D = symm(grad(U));
47         nuSgs = ck*sqrt(k)*delta
48         nuEff = nuSgs + nu
49     @endverbatim
51 SourceFiles
52     dynOneEqEddy.C
54 \*---------------------------------------------------------------------------*/
56 #ifndef dynOneEqEddy_H
57 #define dynOneEqEddy_H
59 #include "GenEddyVisc.H"
60 #include "LESfilter.H"
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 namespace Foam
66 namespace incompressible
68 namespace LESModels
71 /*---------------------------------------------------------------------------*\
72                            Class dynOneEqEddy Declaration
73 \*---------------------------------------------------------------------------*/
75 class dynOneEqEddy
77     public GenEddyVisc
79     // Private data
81         volScalarField k_;
83         autoPtr<LESfilter> filterPtr_;
84         LESfilter& filter_;
87     // Private Member Functions
89         //- Update sub-grid scale fields
90         void updateSubGridScaleFields(const volSymmTensorField& D);
92         //- Calculate ck, ce by filtering the velocity field U.
93         dimensionedScalar ck(const volSymmTensorField& D) const;
94         dimensionedScalar ce(const volSymmTensorField& D) const;
96         // Disallow default bitwise copy construct and assignment
97         dynOneEqEddy(const dynOneEqEddy&);
98         dynOneEqEddy& operator=(const dynOneEqEddy&);
101 public:
103     //- Runtime type information
104     TypeName("dynOneEqEddy");
106     // Constructors
108         //- Construct from components
109         dynOneEqEddy
110         (
111             const volVectorField& U,
112             const surfaceScalarField& phi,
113             transportModel& transport
114         );
117     //- Destructor
118     virtual ~dynOneEqEddy()
119     {}
122     // Member Functions
124         //- Return SGS kinetic energy
125         virtual tmp<volScalarField> k() const
126         {
127             return k_;
128         }
130         //- Return the effective diffusivity for k
131         tmp<volScalarField> DkEff() const
132         {
133             return tmp<volScalarField>
134             (
135                 new volScalarField("DkEff", nuSgs_ + nu())
136             );
137         }
139         //- Correct Eddy-Viscosity and related properties
140         virtual void correct(const tmp<volTensorField>& gradU);
142         //- Read LESProperties dictionary
143         virtual bool read();
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 } // End namespace LESModels
150 } // End namespace incompressible
151 } // End namespace Foam
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 #endif
157 // ************************************************************************* //