ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / oneEqEddy / oneEqEddy.H
blob0f040d829a3cdc7164a7ec8600fa016cc79892ee
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::oneEqEddy
27 Description
28     One Equation Eddy Viscosity Model for incompressible flows
30     Eddy viscosity SGS model using a modeled balance equation to simulate the
31     behaviour of k, hence,
32     \verbatim
33         d/dt(k) + div(U*k) - div(nuEff*grad(k))
34         =
35         -B*L - ce*k^3/2/delta
37     and
39         B = 2/3*k*I - 2*nuSgs*dev(D)
40         Beff = 2/3*k*I - 2*nuEff*dev(D)
42     where
44         D = symm(grad(U));
45         nuSgs = ck*sqrt(k)*delta
46         nuEff = nuSgs + nu
47     \endverbatim
49 SourceFiles
50     oneEqEddy.C
52 \*---------------------------------------------------------------------------*/
54 #ifndef oneEqEddy_H
55 #define oneEqEddy_H
57 #include "GenEddyVisc.H"
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 namespace Foam
63 namespace incompressible
65 namespace LESModels
68 /*---------------------------------------------------------------------------*\
69                            Class oneEqEddy Declaration
70 \*---------------------------------------------------------------------------*/
72 class oneEqEddy
74     public GenEddyVisc
76     // Private data
78         volScalarField k_;
80         dimensionedScalar ck_;
83     // Private Member Functions
85         //- Update sub-grid scale fields
86         void updateSubGridScaleFields();
88         // Disallow default bitwise copy construct and assignment
89         oneEqEddy(const oneEqEddy&);
90         oneEqEddy& operator=(const oneEqEddy&);
93 public:
95     //- Runtime type information
96     TypeName("oneEqEddy");
98     // Constructors
100         //- Construct from components
101         oneEqEddy
102         (
103             const volVectorField& U,
104             const surfaceScalarField& phi,
105             transportModel& transport,
106             const word& turbulenceModelName = turbulenceModel::typeName,
107             const word& modelName = typeName
108         );
111     //- Destructor
112     virtual ~oneEqEddy()
113     {}
116     // Member Functions
118         //- Return SGS kinetic energy
119         virtual tmp<volScalarField> k() const
120         {
121             return k_;
122         }
124         //- Return the effective diffusivity for k
125         tmp<volScalarField> DkEff() const
126         {
127             return tmp<volScalarField>
128             (
129                 new volScalarField("DkEff", nuSgs_ + nu())
130             );
131         }
133         //- Correct Eddy-Viscosity and related properties
134         virtual void correct(const tmp<volTensorField>& gradU);
136         //- Read LESProperties dictionary
137         virtual bool read();
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 } // End namespace LESModels
144 } // End namespace incompressible
145 } // End namespace Foam
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 #endif
151 // ************************************************************************* //