Merge commit 'd3b269b7c6ffa0cd68845adfecdfb849316dba71' into nextRelease
[foam-extend-3.2.git] / src / turbulenceModels / compressible / LES / DeardorffDiffStress / DeardorffDiffStress.H
blobf4b195ad5bdcfb3bb788d11c8a601f1a4d9771a3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     |
5     \\  /    A nd           | For copyright notice see file Copyright
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::compressible::LESModels::DeardorffDiffStress
27 Description
28     Differential SGS Stress Equation Model for compressible flows
30     The DSEM uses a model version of the full balance equation for the SGS
31     stress tensor to simulate the behaviour of B.
32     Thus,
33     @verbatim
34         d/dt(rho*B) + div(rho*U*B) - div(muSgs*grad(B))
35         =
36         P - c1*rho*epsilon/k*B - 0.667*(1 - c1)*rho*epsilon*I - c2*(P - 0.333*trP*I)
38     where
40         k = 0.5*trB,
41         epsilon = ce*k^3/2/delta,
42         epsilon/k = ce*k^1/2/delta
43         P = -rho*(B'L + L'B)
44         muSgs = ck*rho*sqrt(k)*delta
45         muEff = muSgs + mu
46     @endverbatim
48 SourceFiles
49     DeardorffDiffStress.C
51 \*---------------------------------------------------------------------------*/
53 #ifndef compressibleDeardorffDiffStress_H
54 #define compressibleDeardorffDiffStress_H
56 #include "GenSGSStress.H"
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 namespace Foam
62 namespace compressible
64 namespace LESModels
67 /*---------------------------------------------------------------------------*\
68                            Class DeardorffDiffStress Declaration
69 \*---------------------------------------------------------------------------*/
71 class DeardorffDiffStress
73     public GenSGSStress
75     // Private data
77         dimensionedScalar ck_;
78         dimensionedScalar cm_;
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         DeardorffDiffStress(const DeardorffDiffStress&);
88         DeardorffDiffStress& operator=(const DeardorffDiffStress&);
91 public:
93     //- Runtime type information
94     TypeName("DeardorffDiffStress");
96     // Constructors
98         //- Constructor from components
99         DeardorffDiffStress
100         (
101             const volScalarField& rho,
102             const volVectorField& U,
103             const surfaceScalarField& phi,
104             const basicThermo& thermoPhysicalModel
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", muSgs_ + mu())
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 compressible
136 } // End namespace Foam
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 #endif
142 // ************************************************************************* //