Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / LESModel / LESModel.H
blobc6456f433a5881711690494724a12223bb8fb03f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-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 Namespace
25     Foam::incompressible::LESModels
27 Description
28     Namespace for incompressible LES models.
30 Class
31     Foam::incompressible::LESModel
33 Description
34     Base class for all incompressible flow LES SGS models.
36     This class defines the basic interface for an incompressible flow SGS
37     model, and encapsulates data of value to all possible models.
38     In particular this includes references to all the dependent fields
39     (U, phi), the physical viscosity nu, and the LESProperties
40     dictionary, which contains the model selection and model coefficients.
42 SourceFiles
43     LESModel.C
45 \*---------------------------------------------------------------------------*/
47 #ifndef LESModel_H
48 #define LESModel_H
50 #include "incompressible/turbulenceModel/turbulenceModel.H"
51 #include "LESdelta.H"
52 #include "fvm.H"
53 #include "fvc.H"
54 #include "fvMatrices.H"
55 #include "incompressible/transportModel/transportModel.H"
56 #include "wallFvPatch.H"
57 #include "bound.H"
58 #include "autoPtr.H"
59 #include "runTimeSelectionTables.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 namespace Foam
65 namespace incompressible
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 volVectorField& U,
122                 const surfaceScalarField& phi,
123                 transportModel& transport,
124                 const word& turbulenceModelName
125             ),
126             (U, phi, transport, turbulenceModelName)
127         );
130     // Constructors
132         //- Construct from components
133         LESModel
134         (
135             const word& type,
136             const volVectorField& U,
137             const surfaceScalarField& phi,
138             transportModel& transport,
139             const word& turbulenceModelName = turbulenceModel::typeName
140         );
143     // Selectors
145         //- Return a reference to the selected LES model
146         static autoPtr<LESModel> New
147         (
148             const volVectorField& U,
149             const surfaceScalarField& phi,
150             transportModel& transport,
151             const word& turbulenceModelName = turbulenceModel::typeName
152         );
155     //- Destructor
156     virtual ~LESModel()
157     {}
160     // Member Functions
162         // Access
164             //- Const access to the coefficients dictionary,
165             //  which provides info. about choice of models,
166             //  and all related data (particularly model coefficients).
167             inline const dictionary& coeffDict() const
168             {
169                 return coeffDict_;
170             }
172             //- Return the lower allowable limit for k (default: SMALL)
173             const dimensionedScalar& kMin() const
174             {
175                 return kMin_;
176             }
178             //- Allow kMin to be changed
179             dimensionedScalar& kMin()
180             {
181                 return kMin_;
182             }
184             //- Access function to filter width
185             virtual const volScalarField& delta() const
186             {
187                 return delta_();
188             }
191         //- Return the SGS viscosity.
192         virtual tmp<volScalarField> nuSgs() const = 0;
194         //- Return the effective viscosity
195         virtual tmp<volScalarField> nuEff() const
196         {
197             return tmp<volScalarField>
198             (
199                 new volScalarField("nuEff", nuSgs() + nu())
200             );
201         }
203         //- Return the sub-grid stress tensor.
204         virtual tmp<volSymmTensorField> B() const = 0;
206         //- Return the deviatoric part of the effective sub-grid
207         //  turbulence stress tensor including the laminar stress
208         virtual tmp<volSymmTensorField> devBeff() const = 0;
210         //- Returns div(dev(Beff)).
211         //  This is the additional term due to the filtering of the NSE.
212         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const = 0;
215         // RAS compatibility functions for the turbulenceModel base class
217             //- Return the turbulence viscosity
218             virtual tmp<volScalarField> nut() const
219             {
220                 return nuSgs();
221             }
223             //- Return the Reynolds stress tensor
224             virtual tmp<volSymmTensorField> R() const
225             {
226                 return B();
227             }
229             //- Return the effective stress tensor including the laminar stress
230             virtual tmp<volSymmTensorField> devReff() const
231             {
232                 return devBeff();
233             }
235             //- Return the source term for the momentum equation
236             virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const
237             {
238                 return divDevBeff(U);
239             }
242         //- Correct Eddy-Viscosity and related properties.
243         //  This calls correct(const tmp<volTensorField>& gradU) by supplying
244         //  gradU calculated locally.
245         virtual void correct();
247         //- Correct Eddy-Viscosity and related properties
248         virtual void correct(const tmp<volTensorField>& gradU);
250         //- Read LESProperties dictionary
251         virtual bool read();
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 } // End namespace incompressible
258 } // End namespace Foam
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 #endif
264 // ************************************************************************* //