ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / incompressible / LES / DeardorffDiffStress / DeardorffDiffStress.H
blobf97fe17d57ca503809dc7035e2480056d62e8d00
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::incompressible::LESModels::DeardorffDiffStress
27 Description
28     Differential SGS Stress Equation Model for incompressible 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(B) + div(U*B) - div(nuSgs*grad(B))
35         =
36         P - c1*epsilon/k*B - 0.667*(1 - c1)*epsilon*I - c2*(P - 0.333*trP*I)
38     where
40         k = 0.5*tr(B),
41         epsilon = ce*k^3/2/delta,
42         epsilon/k = ce*k^1/2/delta
43         P = -(B'L + L'B)
44         nuSgs = ck*sqrt(k)*delta
45         nuEff = nuSgs + nu
46     @endverbatim
48 SourceFiles
49     DeardorffDiffStress.C
51 \*---------------------------------------------------------------------------*/
53 #ifndef DeardorffDiffStress_H
54 #define DeardorffDiffStress_H
56 #include "GenSGSStress.H"
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 namespace Foam
62 namespace incompressible
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         //- Construct from components
99         DeardorffDiffStress
100         (
101             const volVectorField& U,
102             const surfaceScalarField& phi,
103             transportModel& transport
104         );
107     //- Destructor
108     virtual ~DeardorffDiffStress()
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 // ************************************************************************* //