initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / LES / lowReOneEqEddy / lowReOneEqEddy.H
blob6b892b3cc9769165da53c228cd97af610ddeeb11
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::compressible::LESModels::lowReOneEqEddy
28 Description
29     One Equation Eddy Viscosity Model for compressible flow
31     @verbatim
32         d/dt(rho*k) + div(rho*U*k) - div(muEff*grad(k))
33         =
34         -rho*B*L - ce*rho*k^3/2/delta
36     and
38         B = 2/3*k*I - 2*nuSgs*dev(D)
40     where
42         nuSgsHiRe = ck*sqrt(k)*delta
43         nuSgs = (nu/beta)*(1 - exp(-beta*nuSgsHiRe/nu));
44     @endverbatim
46 SourceFiles
47     lowReOneEqEddy.C
49 \*---------------------------------------------------------------------------*/
51 #ifndef compressibleLowReOneEqEddy_H
52 #define compressibleLowReOneEqEddy_H
54 #include "GenEddyVisc.H"
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 namespace Foam
60 namespace compressible
62 namespace LESModels
65 /*---------------------------------------------------------------------------*\
66                            Class lowReOneEqEddy Declaration
67 \*---------------------------------------------------------------------------*/
69 class lowReOneEqEddy
71     public GenEddyVisc
73     // Private data
75         dimensionedScalar ck_;
76         dimensionedScalar beta_;
78     // Private Member Functions
80         //- Update sub-grid scale fields
81         void updateSubGridScaleFields();
83         // Disallow default bitwise copy construct and assignment
84         lowReOneEqEddy(const lowReOneEqEddy&);
85         lowReOneEqEddy& operator=(const lowReOneEqEddy&);
88 public:
90     //- Runtime type information
91     TypeName("lowReOneEqEddy");
94     // Constructors
96         //- Constructor from components
97         lowReOneEqEddy
98         (
99             const volScalarField& rho,
100             const volVectorField& U,
101             const surfaceScalarField& phi,
102             const basicThermo& thermoPhysicalModel
103         );
106     //- Destructor
107     virtual ~lowReOneEqEddy()
108     {}
111     // Member Functions
113         //- Return the effective diffusivity for k
114         tmp<volScalarField> DkEff() const
115         {
116             return tmp<volScalarField>
117             (
118                 new volScalarField("DkEff", muSgs_ + mu())
119             );
120         }
122         //- Correct Eddy-Viscosity and related properties
123         virtual void correct(const tmp<volTensorField>& gradU);
125         //- Read LESProperties dictionary
126         virtual bool read();
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 } // End namespace LESModels
133 } // End namespace compressible
134 } // End namespace Foam
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 #endif
140 // ************************************************************************* //