initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / LES / oneEqEddy / oneEqEddy.H
blob14a180a046ba87bffa8ed5514e7884bcb6f731c3
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::oneEqEddy
28 Description
29     One Equation Eddy Viscosity Model for incompressible flows
31     Eddy viscosity SGS model using a modeled balance equation to simulate the
32     behaviour of k, hence,
33     @verbatim
34         d/dt(rho*k) + div(rho*U*k) - div(muEff*grad(k))
35         =
36         -rho*D:B - 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         muSgs = ck*rho*sqrt(k)*delta
46     @endverbatim
49 SourceFiles
50     oneEqEddy.C
52 \*---------------------------------------------------------------------------*/
54 #ifndef compressibleOneEqEddy_H
55 #define compressibleOneEqEddy_H
57 #include "GenEddyVisc.H"
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 namespace Foam
63 namespace compressible
65 namespace LESModels
68 /*---------------------------------------------------------------------------*\
69                            Class oneEqEddy Declaration
70 \*---------------------------------------------------------------------------*/
72 class oneEqEddy
74     public GenEddyVisc
76     // Private data
78         dimensionedScalar ck_;
81     // Private Member Functions
83         //- Update sub-grid scale fields
84         void updateSubGridScaleFields();
86         // Disallow default bitwise copy construct and assignment
87         oneEqEddy(const oneEqEddy&);
88         oneEqEddy& operator=(const oneEqEddy&);
91 public:
93     //- Runtime type information
94     TypeName("oneEqEddy");
97     // Constructors
99         //- Constructor from components
100         oneEqEddy
101         (
102             const volScalarField& rho,
103             const volVectorField& U,
104             const surfaceScalarField& phi,
105             const basicThermo& thermoPhysicalModel
106         );
109     // Destructor
110     virtual ~oneEqEddy()
111     {}
114     // Member Functions
116         //- Return the effective diffusivity for k
117         tmp<volScalarField> DkEff() const
118         {
119             return tmp<volScalarField>
120             (
121                 new volScalarField("DkEff", muSgs_ + mu())
122             );
123         }
125         //- Correct Eddy-Viscosity and related properties
126         virtual void correct(const tmp<volTensorField>& gradU);
128         //- Read LESProperties dictionary
129         virtual bool read();
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 } // End namespace LESModels
136 } // End namespace compressible
137 } // End namespace Foam
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 #endif
143 // ************************************************************************* //