ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / LES / LESModel / LESModel.H
blob73515d920c00aa027f14cd0730d74a3d828fbedf
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 Namespace
25     Foam::compressible::LESModels
27 Description
28     Namespace for compressible LES models.
31 Class
32     Foam::compressible::LESModel
34 Description
35     Base class for all compressible flow LES SGS models.
37     This class defines the basic interface for a compressible flow SGS
38     model, and encapsulates data of value to all possible models.
39     In particular this includes references to all the dependent fields
40     (rho, U, phi), the physical viscosity mu, and the LESProperties
41     dictionary, which contains the model selection and model coefficients.
43 SourceFiles
44     LESModel.C
46 \*---------------------------------------------------------------------------*/
48 #ifndef compressibleLESModel_H
49 #define compressibleLESModel_H
51 #include "compressible/turbulenceModel/turbulenceModel.H"
52 #include "LESdelta.H"
53 #include "fvm.H"
54 #include "fvc.H"
55 #include "fvMatrices.H"
56 #include "basicThermo.H"
57 #include "bound.H"
58 #include "autoPtr.H"
59 #include "runTimeSelectionTables.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 namespace Foam
65 namespace compressible
68 /*---------------------------------------------------------------------------*\
69                            Class LESModel Declaration
70 \*---------------------------------------------------------------------------*/
72 class LESModel
74     public turbulenceModel,
75     public IOdictionary
78 protected:
80     // Protected data
82         Switch printCoeffs_;
83         dictionary coeffDict_;
85         dimensionedScalar kMin_;
87         autoPtr<LESdelta> delta_;
90     // Protected Member Functions
92         //- Print model coefficients
93         virtual void printCoeffs();
96 private:
98     // Private Member Functions
100         //- Disallow default bitwise copy construct
101         LESModel(const LESModel&);
103         //- Disallow default bitwise assignment
104         LESModel& operator=(const LESModel&);
107 public:
109     //- Runtime type information
110     TypeName("LESModel");
113     // Declare run-time constructor selection table
115         declareRunTimeSelectionTable
116         (
117             autoPtr,
118             LESModel,
119             dictionary,
120             (
121                 const volScalarField& rho,
122                 const volVectorField& U,
123                 const surfaceScalarField& phi,
124                 const basicThermo& thermoPhysicalModel,
125                 const word& turbulenceModelName
126             ),
127             (rho, U, phi, thermoPhysicalModel, turbulenceModelName)
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             const word& turbulenceModelName = turbulenceModel::typeName
142         );
145     // Selectors
147         //- Return a reference to the selected LES model
148         static autoPtr<LESModel> New
149         (
150             const volScalarField& rho,
151             const volVectorField& U,
152             const surfaceScalarField& phi,
153             const basicThermo& thermoPhysicalModel,
154             const word& turbulenceModelName = turbulenceModel::typeName
155         );
158     //- Destructor
159     virtual ~LESModel()
160     {}
163     // Member Functions
165         // Access
167             //- Const access to the coefficients dictionary,
168             //  which provides info. about choice of models,
169             //  and all related data (particularly model coefficients).
170             inline const dictionary& coeffDict() const
171             {
172                 return coeffDict_;
173             }
175             //- Return the lower allowable limit for k (default: SMALL)
176             const dimensionedScalar& kMin() const
177             {
178                 return kMin_;
179             }
181             //- Allow kMin to be changed
182             dimensionedScalar& kMin()
183             {
184                 return kMin_;
185             }
187             //- Access function to filter width
188             inline const volScalarField& delta() const
189             {
190                 return delta_();
191             }
194         //- Return the SGS turbulent viscosity
195         virtual tmp<volScalarField> muSgs() const = 0;
197         //- Return the effective viscosity
198         virtual tmp<volScalarField> muEff() const
199         {
200             return tmp<volScalarField>
201             (
202                 new volScalarField("muEff", muSgs() + mu())
203             );
204         }
206         //- Return the SGS turbulent thermal diffusivity
207         virtual tmp<volScalarField> alphaSgs() const = 0;
209         //- Return the effective thermal diffusivity
210         virtual tmp<volScalarField> alphaEff() const
211         {
212             return tmp<volScalarField>
213             (
214                 new volScalarField("alphaEff", alphaSgs() + alpha())
215             );
216         }
218         //- Return the effective turbulence thermal diffusivity for a patch
219         virtual tmp<scalarField> alphaEff(const label patchI) const
220         {
221             return
222                 alphaSgs()().boundaryField()[patchI]
223               + alpha().boundaryField()[patchI];
224         }
226         //- Return the sub-grid stress tensor.
227         virtual tmp<volSymmTensorField> B() const = 0;
229         //- Return the deviatoric part of the effective sub-grid
230         //  turbulence stress tensor including the laminar stress
231         virtual tmp<volSymmTensorField> devRhoBeff() const = 0;
233         //- Returns div(rho*dev(B)).
234         // This is the additional term due to the filtering of the NSE.
235         virtual tmp<fvVectorMatrix> divDevRhoBeff(volVectorField& U) const = 0;
238         // RAS compatibility functions for the turbulenceModel base class
240             //- Return the turbulence viscosity
241             virtual tmp<volScalarField> mut() const
242             {
243                 return muSgs();
244             }
246             //- Return the turbulence thermal diffusivity
247             virtual tmp<volScalarField> alphat() const
248             {
249                 return alphaSgs();
250             }
252             //- Return the Reynolds stress tensor
253             virtual tmp<volSymmTensorField> R() const
254             {
255                 return B();
256             }
258             //- Return the effective stress tensor including the laminar stress
259             virtual tmp<volSymmTensorField> devRhoReff() const
260             {
261                 return devRhoBeff();
262             }
264             //- Return the source term for the momentum equation
265             virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const
266             {
267                 return divDevRhoBeff(U);
268             }
271         //- Correct Eddy-Viscosity and related properties.
272         //  This calls correct(const tmp<volTensorField>& gradU) by supplying
273         //  gradU calculated locally.
274         virtual void correct();
276         //- Correct Eddy-Viscosity and related properties
277         virtual void correct(const tmp<volTensorField>& gradU);
279         //- Read LESProperties dictionary
280         virtual bool read();
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 } // End namespace compressible
287 } // End namespace Foam
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 #endif
293 // ************************************************************************* //