ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / incompressible / LES / LESModel / LESModel.H
blob45bb40ddedc4268ab69e5334005b6f90d98c88a9
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 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 k0_;
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& lamTransportModel
124             ),
125             (U, phi, lamTransportModel)
126         );
129     // Constructors
131         //- Construct from components
132         LESModel
133         (
134             const word& type,
135             const volVectorField& U,
136             const surfaceScalarField& phi,
137             transportModel& lamTransportModel
138         );
141     // Selectors
143         //- Return a reference to the selected LES model
144         static autoPtr<LESModel> New
145         (
146             const volVectorField& U,
147             const surfaceScalarField& phi,
148             transportModel& lamTransportModel
149         );
152     //- Destructor
153     virtual ~LESModel()
154     {}
157     // Member Functions
159         //- Const access to the coefficients dictionary,
160         //  which provides info. about choice of models,
161         //  and all related data (particularly model coefficients).
162         inline const dictionary& coeffDict() const
163         {
164             return coeffDict_;
165         }
167         //- Access function to filter width
168         virtual const volScalarField& delta() const
169         {
170             return delta_();
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         }
186         //- Return the SGS turbulent kinetic energy.
187         virtual tmp<volScalarField> k() const = 0;
189         //- Return the SGS turbulent dissipation.
190         virtual tmp<volScalarField> epsilon() const = 0;
192         //- Return the SGS viscosity.
193         virtual tmp<volScalarField> nuSgs() const = 0;
195         //- Return the effective viscosity
196         virtual tmp<volScalarField> nuEff() const
197         {
198             return tmp<volScalarField>
199             (
200                 new volScalarField("nuEff", nuSgs() + nu())
201             );
202         }
204         //- Return the sub-grid stress tensor.
205         virtual tmp<volSymmTensorField> B() const = 0;
207         //- Return the deviatoric part of the effective sub-grid
208         //  turbulence stress tensor including the laminar stress
209         virtual tmp<volSymmTensorField> devBeff() const = 0;
211         //- Returns div(dev(Beff)).
212         //  This is the additional term due to the filtering of the NSE.
213         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const = 0;
216         // RAS compatibility functions for the turbulenceModel base class
218             //- Return the turbulence viscosity
219             virtual tmp<volScalarField> nut() const
220             {
221                 return nuSgs();
222             }
224             //- Return the Reynolds stress tensor
225             virtual tmp<volSymmTensorField> R() const
226             {
227                 return B();
228             }
230             //- Return the effective stress tensor including the laminar stress
231             virtual tmp<volSymmTensorField> devReff() const
232             {
233                 return devBeff();
234             }
236             //- Return the source term for the momentum equation
237             virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const
238             {
239                 return divDevBeff(U);
240             }
243         //- Correct Eddy-Viscosity and related properties.
244         //  This calls correct(const tmp<volTensorField>& gradU) by supplying
245         //  gradU calculated locally.
246         void correct();
248         //- Correct Eddy-Viscosity and related properties
249         virtual void correct(const tmp<volTensorField>& gradU);
251         //- Read LESProperties dictionary
252         virtual bool read() = 0;
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 } // End namespace incompressible
259 } // End namespace Foam
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 #endif
265 // ************************************************************************* //