ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / homogeneousDynOneEqEddy / homogeneousDynOneEqEddy.H
blob6de3992574ff8a7b5f07706448eef090ea063adc
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::homogeneousDynOneEqEddy
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     homogeneousDynOneEqEddy.C
54 \*---------------------------------------------------------------------------*/
56 #ifndef homogeneousDynOneEqEddy_H
57 #define homogeneousDynOneEqEddy_H
59 #include "GenEddyVisc.H"
60 #include "LESfilter.H"
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 namespace Foam
66 namespace incompressible
68 namespace LESModels
71 /*---------------------------------------------------------------------------*\
72                            Class homogeneousDynOneEqEddy Declaration
73 \*---------------------------------------------------------------------------*/
75 class homogeneousDynOneEqEddy
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         homogeneousDynOneEqEddy(const homogeneousDynOneEqEddy&);
98         homogeneousDynOneEqEddy& operator=(const homogeneousDynOneEqEddy&);
101 public:
103     //- Runtime type information
104     TypeName("homogeneousDynOneEqEddy");
106     // Constructors
108         //- Construct from components
109         homogeneousDynOneEqEddy
110         (
111             const volVectorField& U,
112             const surfaceScalarField& phi,
113             transportModel& transport,
114             const word& turbulenceModelName = turbulenceModel::typeName,
115             const word& modelName = typeName
116         );
119     //- Destructor
120     virtual ~homogeneousDynOneEqEddy()
121     {}
124     // Member Functions
126         //- Return SGS kinetic energy
127         virtual tmp<volScalarField> k() const
128         {
129             return k_;
130         }
132         //- Return the effective diffusivity for k
133         tmp<volScalarField> DkEff() const
134         {
135             return tmp<volScalarField>
136             (
137                 new volScalarField("DkEff", nuSgs_ + nu())
138             );
139         }
141         //- Correct Eddy-Viscosity and related properties
142         virtual void correct(const tmp<volTensorField>& gradU);
144         //- Read LESProperties dictionary
145         virtual bool read();
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 } // End namespace LESModels
152 } // End namespace incompressible
153 } // End namespace Foam
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 #endif
159 // ************************************************************************* //