initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / GenSGSStress / GenSGSStress.H
blobbf1a0eb73891df5528ac913a175c8f7e813a2cbe
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::GenSGSStress
28 Description
29     General base class for all incompressible models that directly
30     solve for the SGS stress tensor B.
32     Contains tensor fields B (the SGS stress tensor) as well as scalar
33     fields for k (SGS turbulent energy) gamma (SGS viscosity) and epsilon
34     (SGS dissipation).
36 SourceFiles
37     GenSGSStress.C
39 \*---------------------------------------------------------------------------*/
41 #ifndef GenSGSStress_H
42 #define GenSGSStress_H
44 #include "LESModel.H"
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 namespace Foam
50 namespace incompressible
52 namespace LESModels
55 /*---------------------------------------------------------------------------*\
56                            Class GenSGSStress Declaration
57 \*---------------------------------------------------------------------------*/
59 class GenSGSStress
61     virtual public LESModel
63     // Private Member Functions
65         // Disallow default bitwise copy construct and assignment
66         GenSGSStress(const GenSGSStress&);
67         GenSGSStress& operator=(const GenSGSStress&);
70 protected:
72         dimensionedScalar ce_;
74         dimensionedScalar couplingFactor_;
76         volSymmTensorField B_;
77         volScalarField nuSgs_;
80 public:
82     // Constructors
84         //- Construct from components
85         GenSGSStress
86         (
87             const volVectorField& U,
88             const surfaceScalarField& phi,
89             transportModel& transport
90         );
93     //- Destructor
94     virtual ~GenSGSStress()
95     {}
98     // Member Functions
100         //- Return the SGS turbulent kinetic energy.
101         virtual tmp<volScalarField> k() const
102         {
103             return 0.5*tr(B_);
104         }
106         //- Return the SGS turbulent dissipation.
107         virtual tmp<volScalarField> epsilon() const
108         {
109             volScalarField K = k();
110             return ce_*K*sqrt(K)/delta();
111         }
113         //- Return the SGS viscosity.
114         virtual tmp<volScalarField> nuSgs() const
115         {
116             return nuSgs_;
117         }
119         //- Return the sub-grid stress tensor.
120         virtual tmp<volSymmTensorField> B() const
121         {
122             return B_;
123         }
125         //- Return the effective sub-grid turbulence stress tensor
126         //  including the laminar stress
127         virtual tmp<volSymmTensorField> devBeff() const;
129         //- Returns div(B).
130         // This is the additional term due to the filtering of the NSE.
131         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
133         //- Read LESProperties dictionary
134         virtual bool read();
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 } // End namespace LESModels
141 } // End namespace incompressible
142 } // End namespace Foam
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 #endif
148 // ************************************************************************* //