initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / DeardorffDiffStress / DeardorffDiffStress.H
blob9305473130877ee5cb4d9a7eda7576302a27b090
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::incompressible::LESModels::DeardorffDiffStress
28 Description
29     Differential SGS Stress Equation Model for incompressible flows
31     The DSEM uses a model version of the full balance equation for the SGS
32     stress tensor to simulate the behaviour of B.
33     Thus,
34     @verbatim
35         d/dt(B) + div(U*B) - div(nuSgs*grad(B))
36         =
37         P - c1*epsilon/k*B - 0.667*(1 - c1)*epsilon*I - c2*(P - 0.333*trP*I)
39     where
41         k = 0.5*tr(B),
42         epsilon = ce*k^3/2/delta,
43         epsilon/k = ce*k^1/2/delta
44         P = -(B'L + L'B)
45         nuSgs = ck*sqrt(k)*delta
46         nuEff = nuSgs + nu
47     @endverbatim
49 SourceFiles
50     DeardorffDiffStress.C
52 \*---------------------------------------------------------------------------*/
54 #ifndef DeardorffDiffStress_H
55 #define DeardorffDiffStress_H
57 #include "GenSGSStress.H"
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 namespace Foam
63 namespace incompressible
65 namespace LESModels
68 /*---------------------------------------------------------------------------*\
69                            Class DeardorffDiffStress Declaration
70 \*---------------------------------------------------------------------------*/
72 class DeardorffDiffStress
74     public GenSGSStress
76     // Private data
78         dimensionedScalar ck_;
79         dimensionedScalar cm_;
82     // Private Member Functions
84         //- Update sub-grid scale fields
85         void updateSubGridScaleFields(const volScalarField& K);
87         // Disallow default bitwise copy construct and assignment
88         DeardorffDiffStress(const DeardorffDiffStress&);
89         DeardorffDiffStress& operator=(const DeardorffDiffStress&);
92 public:
94     //- Runtime type information
95     TypeName("DeardorffDiffStress");
97     // Constructors
99         //- Construct from components
100         DeardorffDiffStress
101         (
102             const volVectorField& U,
103             const surfaceScalarField& phi,
104             transportModel& transport
105         );
108     //- Destructor
109     virtual ~DeardorffDiffStress()
110     {}
113     // Member Functions
115         //- Return the effective diffusivity for B
116         tmp<volScalarField> DBEff() const
117         {
118             return tmp<volScalarField>
119             (
120                 new volScalarField("DBEff", nuSgs_ + nu())
121             );
122         }
124         //- Correct Eddy-Viscosity and related properties
125         virtual void correct(const tmp<volTensorField>& gradU);
127         //- Read LESProperties dictionary
128         virtual bool read();
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 } // End namespace LESModels
135 } // End namespace incompressible
136 } // End namespace Foam
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 #endif
142 // ************************************************************************* //