initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / LES / LESModel / LESModel.H
blobcf019ea9a8aca9b93a709e94507a74c08a9ef6c1
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 Namespace
26     Foam::compressible::LESModels
28 Description
29     Namespace for compressible LES models.
32 Class
33     Foam::compressible::LESModel
35 Description
36     Base class for all compressible flow LES SGS models.
38     This class defines the basic interface for a compressible flow SGS
39     model, and encapsulates data of value to all possible models.
40     In particular this includes references to all the dependent fields
41     (rho, U, phi), the physical viscosity mu, and the LESProperties
42     dictionary, which contains the model selection and model coefficients.
44 SourceFiles
45     LESModel.C
47 \*---------------------------------------------------------------------------*/
49 #ifndef compressibleLESModel_H
50 #define compressibleLESModel_H
52 #include "compressible/turbulenceModel/turbulenceModel.H"
53 #include "LESdelta.H"
54 #include "fvm.H"
55 #include "fvc.H"
56 #include "fvMatrices.H"
57 #include "basicThermo.H"
58 #include "bound.H"
59 #include "autoPtr.H"
60 #include "runTimeSelectionTables.H"
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 namespace Foam
66 namespace compressible
69 /*---------------------------------------------------------------------------*\
70                            Class LESModel Declaration
71 \*---------------------------------------------------------------------------*/
73 class LESModel
75     public turbulenceModel,
76     public IOdictionary
79 protected:
81     // Protected data
83         Switch printCoeffs_;
84         dictionary coeffDict_;
86         dimensionedScalar k0_;
88         autoPtr<LESdelta> delta_;
91     // Protected Member Functions
93         //- Print model coefficients
94         virtual void printCoeffs();
97 private:
99     // Private Member Functions
101         //- Disallow default bitwise copy construct
102         LESModel(const LESModel&);
104         //- Disallow default bitwise assignment
105         LESModel& operator=(const LESModel&);
108 public:
110     //- Runtime type information
111     TypeName("LESModel");
114     // Declare run-time constructor selection table
116         declareRunTimeSelectionTable
117         (
118             autoPtr,
119             LESModel,
120             dictionary,
121             (
122                 const volScalarField& rho,
123                 const volVectorField& U,
124                 const surfaceScalarField& phi,
125                 const basicThermo& thermoPhysicalModel
126             ),
127             (rho, U, phi, thermoPhysicalModel)
128         );
131     // Constructors
133         //- Construct from components
134         LESModel
135         (
136             const word& type,
137             const volScalarField& rho,
138             const volVectorField& U,
139             const surfaceScalarField& phi,
140             const basicThermo& thermoPhysicalModel
141         );
144     // Selectors
146         //- Return a reference to the selected LES model
147         static autoPtr<LESModel> New
148         (
149             const volScalarField& rho,
150             const volVectorField& U,
151             const surfaceScalarField& phi,
152             const basicThermo& thermoPhysicalModel
153         );
156     //- Destructor
157     virtual ~LESModel()
158     {}
161     // Member Functions
163         // Access
165             //- Const access to the coefficients dictionary,
166             //  which provides info. about choice of models,
167             //  and all related data (particularly model coefficients).
168             inline const dictionary& coeffDict() const
169             {
170                 return coeffDict_;
171             }
173             //- Return the value of k0 which k is not allowed to be less than
174             const dimensionedScalar& k0() const
175             {
176                 return k0_;
177             }
179             //- Allow k0 to be changed
180             dimensionedScalar& k0()
181             {
182                 return k0_;
183             }
185             //- Access function to filter width
186             inline const volScalarField& delta() const
187             {
188                 return delta_();
189             }
192         //- Return the SGS turbulent kinetic energy.
193         virtual tmp<volScalarField> k() const = 0;
195         //- Return the SGS turbulent dissipation.
196         virtual tmp<volScalarField> epsilon() const = 0;
198         //- Return the SGS turbulent viscosity
199         virtual tmp<volScalarField> muSgs() const = 0;
201         //- Return the effective viscosity
202         virtual tmp<volScalarField> muEff() const
203         {
204             return tmp<volScalarField>
205             (
206                 new volScalarField("muEff", muSgs() + mu())
207             );
208         }
210         //- Return the SGS turbulent thermal diffusivity
211         virtual tmp<volScalarField> alphaSgs() const = 0;
213         //- Return the SGS thermal conductivity.
214         virtual tmp<volScalarField> alphaEff() const = 0;
216         //- Return the sub-grid stress tensor.
217         virtual tmp<volSymmTensorField> B() const = 0;
219         //- Return the deviatoric part of the effective sub-grid
220         //  turbulence stress tensor including the laminar stress
221         virtual tmp<volSymmTensorField> devRhoBeff() const = 0;
223         //- Returns div(rho*dev(B)).
224         // This is the additional term due to the filtering of the NSE.
225         virtual tmp<fvVectorMatrix> divDevRhoBeff(volVectorField& U) const = 0;
228         // RAS compatibility functions for the turbulenceModel base class
230             //- Return the turbulence viscosity
231             virtual tmp<volScalarField> mut() const
232             {
233                 return muSgs();
234             }
236             //- Return the turbulence thermal diffusivity
237             virtual tmp<volScalarField> alphat() const
238             {
239                 return alphaSgs();
240             }
242             //- Return the Reynolds stress tensor
243             virtual tmp<volSymmTensorField> R() const
244             {
245                 return B();
246             }
248             //- Return the effective stress tensor including the laminar stress
249             virtual tmp<volSymmTensorField> devRhoReff() const
250             {
251                 return devRhoBeff();
252             }
254             //- Return the source term for the momentum equation
255             virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const
256             {
257                 return divDevRhoBeff(U);
258             }
261         //- Correct Eddy-Viscosity and related properties.
262         //  This calls correct(const tmp<volTensorField>& gradU) by supplying
263         //  gradU calculated locally.
264         void correct();
266         //- Correct Eddy-Viscosity and related properties
267         virtual void correct(const tmp<volTensorField>& gradU);
269         //- Read LESProperties dictionary
270         virtual bool read() = 0;
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 } // End namespace compressible
277 } // End namespace Foam
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 #endif
283 // ************************************************************************* //