initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / LRRDiffStress / LRRDiffStress.H
blob1436b5ac3a979998b044645542763b9561968721
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::LRRDiffStress
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, hence,
33     @verbatim
34         d/dt(B) + div(U*B) - div(nuSgs*grad(B))
35         =
36         P - c1*e/k*B - 0.667*(1 - c1)*e*I - c2*(P - 0.333*trP*I)
37     where
38         k = 0.5*trB,
39         epsilon = ce*k^3/2/delta
40         epsilon/k = ce*k^1/2/delta
41         P = -(B'L + L'B)
42         nuEff = ck*sqrt(k)*delta + nu
43     @endverbatim
45     This version from Launder, Rece & Rodi 1975
47 SourceFiles
48     LRRDiffStress.C
50 \*---------------------------------------------------------------------------*/
52 #ifndef LRRDiffStress_H
53 #define LRRDiffStress_H
55 #include "GenSGSStress.H"
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 namespace Foam
61 namespace incompressible
63 namespace LESModels
66 /*---------------------------------------------------------------------------*\
67                            Class LRRDiffStress Declaration
68 \*---------------------------------------------------------------------------*/
70 class LRRDiffStress
72     public GenSGSStress
74     // Private data
76         dimensionedScalar ck_;
77         dimensionedScalar c1_;
78         dimensionedScalar c2_;
81     // Private Member Functions
83         //- Update sub-grid scale fields
84         void updateSubGridScaleFields(const volScalarField& K);
86         // Disallow default bitwise copy construct and assignment
87         LRRDiffStress(const LRRDiffStress&);
88         LRRDiffStress& operator=(const LRRDiffStress&);
91 public:
93     //- Runtime type information
94     TypeName("LRRDiffStress");
96     // Constructors
98         //- Construct from components
99         LRRDiffStress
100         (
101             const volVectorField& U,
102             const surfaceScalarField& phi,
103             transportModel& transport
104         );
107     //- Destructor
108     virtual ~LRRDiffStress()
109     {}
112     // Member Functions
114         //- Return the effective diffusivity for B
115         tmp<volScalarField> DBEff() const
116         {
117             return tmp<volScalarField>
118             (
119                 new volScalarField("DBEff", nuSgs_ + nu())
120             );
121         }
123         //- Correct Eddy-Viscosity and related properties
124         virtual void correct(const tmp<volTensorField>& gradU);
126         //- Read LESProperties dictionary
127         virtual bool read();
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 } // End namespace LESModels
134 } // End namespace incompressible
135 } // End namespace Foam
137 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 #endif
141 // ************************************************************************* //