ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / compressible / LES / GenSGSStress / GenSGSStress.H
blob1d5423e663dca157d5b8e8902569ae46ea44327f
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::compressible::LESModels::GenSGSStress
27 Description
28     General base class for all compressible 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 compressibleGenSGSStress_H
41 #define compressibleGenSGSStress_H
43 #include "LESModel.H"
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
49 namespace compressible
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     // Model coefficients
73         dimensionedScalar ce_;
74         dimensionedScalar Prt_;
76     // Fields
78         volSymmTensorField B_;
79         volScalarField muSgs_;
80         volScalarField alphaSgs_;
83 public:
85     // Constructors
87         //- Constructor from components
88         GenSGSStress
89         (
90             const volScalarField& rho,
91             const volVectorField& U,
92             const surfaceScalarField& phi,
93             const basicThermo& thermoPhysicalModel
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             volScalarField K = k();
114             return ce_*K*sqrt(K)/delta();
115         }
117         //- Return the SGS viscosity
118         virtual tmp<volScalarField> muSgs() const
119         {
120             return muSgs_;
121         }
123         //- Return the SGS thermal diffusivity
124         virtual tmp<volScalarField> alphaSgs() const
125         {
126             return alphaSgs_;
127         }
129         //- Return thermal conductivity
130         virtual tmp<volScalarField> alphaEff() const
131         {
132             return tmp<volScalarField>
133             (
134                 new volScalarField("alphaEff", alphaSgs_ + alpha())
135             );
136         }
138         //- Return the sub-grid stress tensor
139         virtual tmp<volSymmTensorField> B() const
140         {
141             return B_;
142         }
144         //- Return the deviatoric part of the effective sub-grid
145         //  turbulence stress tensor including the laminar stress
146         virtual tmp<volSymmTensorField> devRhoBeff() const;
148         //- Returns divergence of B : i.e. the additional term in the
149         //  filtered NSE
150         virtual tmp<fvVectorMatrix> divDevRhoBeff(volVectorField& U) const;
152         //- Correct Eddy-Viscosity and related properties
153         virtual void correct(const tmp<volTensorField>& gradU);
155         //- Read LESProperties dictionary
156         virtual bool read();
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 } // End namespace LESModels
163 } // End namespace compressible
164 } // End namespace Foam
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 #endif
170 // ************************************************************************* //