ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / LES / homogeneousDynOneEqEddy / homogeneousDynOneEqEddy.H
blob8464dd8cf77d67f60769cd7735ecf4d0532bb4fb
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::compressible::LESModels::homogeneousDynOneEqEddy
27 Description
28     One Equation Eddy Viscosity Model for compressible flows.
30     Eddy viscosity SGS model using a modeled balance equation to simulate
31     the behaviour of k.
32     Thus
33     \verbatim
34         d/dt(k) + div(U*k) - div(nuSgs*grad(k))
35         =
36         -rho*B*L - ce*rho*k^3/2/delta
38     and
40         B = 2/3*k*I - 2*nuSgs*dev(D)
42     where
44         D = symm(grad(U));
45         nuSgs = ck*sqrt(k)*delta
46     \endverbatim
48 SourceFiles
49     homogeneousDynOneEqEddy.C
51 \*---------------------------------------------------------------------------*/
53 #ifndef compressibleHomogeneousDynOneEqEddy_H
54 #define compressibleHomogeneousDynOneEqEddy_H
56 #include "GenEddyVisc.H"
57 #include "LESfilter.H"
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 namespace Foam
63 namespace compressible
65 namespace LESModels
68 /*---------------------------------------------------------------------------*\
69                            Class homogeneousDynOneEqEddy Declaration
70 \*---------------------------------------------------------------------------*/
72 class homogeneousDynOneEqEddy
74     public GenEddyVisc
76     // Private data
78         volScalarField k_;
80         autoPtr<LESfilter> filterPtr_;
81         LESfilter& filter_;
84     // Private Member Functions
86         //- Update sub-grid scale fields
87         void updateSubGridScaleFields(const volSymmTensorField& D);
89         //- Calculate ck, ce by filtering the velocity field U.
90         dimensionedScalar ck_(const volSymmTensorField& D) const;
91         dimensionedScalar ce_(const volSymmTensorField& D) const;
93         // Disallow default bitwise copy construct and assignment
94         homogeneousDynOneEqEddy(const homogeneousDynOneEqEddy&);
95         homogeneousDynOneEqEddy& operator=(const homogeneousDynOneEqEddy&);
98 public:
100     //- Runtime type information
101     TypeName("homogeneousDynOneEqEddy");
104     // Constructors
106         //- Constructor from components
107         homogeneousDynOneEqEddy
108         (
109             const volScalarField& rho,
110             const volVectorField& U,
111             const surfaceScalarField& phi,
112             const basicThermo& thermoPhysicalModel,
113             const word& turbulenceModelName = turbulenceModel::typeName,
114             const word& modelName = typeName
115         );
118     //- Destructor
119     virtual ~homogeneousDynOneEqEddy()
120     {}
123     // Member Functions
125         //- Return SGS kinetic energy
126         virtual tmp<volScalarField> k() const
127         {
128             return k_;
129         }
131         //- Return the effective diffusivity for k
132         tmp<volScalarField> DkEff() const
133         {
134             return tmp<volScalarField>
135             (
136                 new volScalarField("DkEff", muSgs_ + mu())
137             );
138         }
140         //- Correct Eddy-Viscosity and related properties
141         virtual void correct(const tmp<volTensorField>& gradU);
143         //- Read LESProperties dictionary
144         virtual bool read();
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 } // End namespace LESModels
151 } // End namespace compressible
152 } // End namespace Foam
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 #endif
158 // ************************************************************************* //