ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / GenSGSStress / GenSGSStress.H
blob36d84f9c771a1d6f3aa23c4a107201acbf8b0e5f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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::GenSGSStress
27 Description
28     General base class for all incompressible models that directly
29     solve for the SGS stress tensor B.
31     Contains tensor fields B (the SGS stress tensor) as well as scalar
32     fields for k (SGS turbulent energy) gamma (SGS viscosity) and epsilon
33     (SGS dissipation).
35 SourceFiles
36     GenSGSStress.C
38 \*---------------------------------------------------------------------------*/
40 #ifndef GenSGSStress_H
41 #define GenSGSStress_H
43 #include "LESModel.H"
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
49 namespace incompressible
51 namespace LESModels
54 /*---------------------------------------------------------------------------*\
55                            Class GenSGSStress Declaration
56 \*---------------------------------------------------------------------------*/
58 class GenSGSStress
60     virtual public LESModel
62     // Private Member Functions
64         // Disallow default bitwise copy construct and assignment
65         GenSGSStress(const GenSGSStress&);
66         GenSGSStress& operator=(const GenSGSStress&);
69 protected:
71         dimensionedScalar ce_;
73         dimensionedScalar couplingFactor_;
75         volSymmTensorField B_;
76         volScalarField nuSgs_;
79 public:
81     //- Partial Runtime type information
82     static const word typeName;
84     // Constructors
86         //- Construct from components
87         GenSGSStress
88         (
89             const volVectorField& U,
90             const surfaceScalarField& phi,
91             transportModel& transport,
92             const word& turbulenceModelName = turbulenceModel::typeName,
93             const word& modelName = typeName
94         );
97     //- Destructor
98     virtual ~GenSGSStress()
99     {}
102     // Member Functions
104         //- Return the SGS turbulent kinetic energy.
105         virtual tmp<volScalarField> k() const
106         {
107             return 0.5*tr(B_);
108         }
110         //- Return the SGS turbulent dissipation.
111         virtual tmp<volScalarField> epsilon() const
112         {
113             const volScalarField K(k());
114             return ce_*K*sqrt(K)/delta();
115         }
117         //- Return the SGS viscosity.
118         virtual tmp<volScalarField> nuSgs() const
119         {
120             return nuSgs_;
121         }
123         //- Return the sub-grid stress tensor.
124         virtual tmp<volSymmTensorField> B() const
125         {
126             return B_;
127         }
129         //- Return the effective sub-grid turbulence stress tensor
130         //  including the laminar stress
131         virtual tmp<volSymmTensorField> devBeff() const;
133         //- Returns div(B).
134         // This is the additional term due to the filtering of the NSE.
135         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
137         //- Read LESProperties dictionary
138         virtual bool read();
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 } // End namespace LESModels
145 } // End namespace incompressible
146 } // End namespace Foam
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 #endif
152 // ************************************************************************* //