ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / compressible / RAS / RASModel / RASModel.H
blobc80417a7c1939c555194233fca24c052663a2019
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::compressible::RASModels
27 Description
28     Namespace for compressible RAS turbulence models.
31 Class
32     Foam::compressible::RASModel
34 Description
35     Abstract base class for turbulence models for compressible and combusting
36     flows.
38 SourceFiles
39     RASModel.C
40     newTurbulenceModel.C
42 \*---------------------------------------------------------------------------*/
44 #ifndef compressibleRASModel_H
45 #define compressibleRASModel_H
47 #include "compressible/turbulenceModel/turbulenceModel.H"
48 #include "volFields.H"
49 #include "surfaceFields.H"
50 #include "nearWallDist.H"
51 #include "fvm.H"
52 #include "fvc.H"
53 #include "fvMatrices.H"
54 #include "basicThermo.H"
55 #include "IOdictionary.H"
56 #include "Switch.H"
57 #include "bound.H"
58 #include "autoPtr.H"
59 #include "runTimeSelectionTables.H"
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 namespace Foam
65 namespace compressible
68 /*---------------------------------------------------------------------------*\
69                            Class RASModel Declaration
70 \*---------------------------------------------------------------------------*/
72 class RASModel
74     public turbulenceModel,
75     public IOdictionary
78 protected:
80     // Protected data
82         //- Turbulence on/off flag
83         Switch turbulence_;
85         //- Flag to print the model coeffs at run-time
86         Switch printCoeffs_;
88         //- Model coefficients dictionary
89         dictionary coeffDict_;
91         //- Value of y+ at the edge of the laminar sublayer
92         scalar yPlusLam_;
94         //- Lower limit of k
95         dimensionedScalar k0_;
97         //- Lower limit of epsilon
98         dimensionedScalar epsilon0_;
100         //- Small epsilon value used to avoid divide by zero
101         dimensionedScalar epsilonSmall_;
103         //- Lower limit for omega
104         dimensionedScalar omega0_;
106         //- Small omega value used to avoid divide by zero
107         dimensionedScalar omegaSmall_;
109         //- Near wall distance boundary field
110         nearWallDist y_;
113     // Protected member functions
115         //- Print model coefficients
116         virtual void printCoeffs();
119 private:
121     // Private Member Functions
123         //- Disallow default bitwise copy construct
124         RASModel(const RASModel&);
126         //- Disallow default bitwise assignment
127         void operator=(const RASModel&);
130 public:
132     //- Runtime type information
133     TypeName("RASModel");
136     // Declare run-time constructor selection table
138         declareRunTimeSelectionTable
139         (
140             autoPtr,
141             RASModel,
142             dictionary,
143             (
144                 const volScalarField& rho,
145                 const volVectorField& U,
146                 const surfaceScalarField& phi,
147                 const basicThermo& thermoPhysicalModel
148             ),
149             (rho, U, phi, thermoPhysicalModel)
150         );
153     // Constructors
155         //- Construct from components
156         RASModel
157         (
158             const word& type,
159             const volScalarField& rho,
160             const volVectorField& U,
161             const surfaceScalarField& phi,
162             const basicThermo& thermoPhysicalModel
163         );
166     // Selectors
168         //- Return a reference to the selected turbulence model
169         static autoPtr<RASModel> New
170         (
171             const volScalarField& rho,
172             const volVectorField& U,
173             const surfaceScalarField& phi,
174             const basicThermo& thermoPhysicalModel
175         );
178     //- Destructor
179     virtual ~RASModel()
180     {}
183     // Member Functions
185         // Access
187             //- Return the value of k0 which k is not allowed to be less than
188             const dimensionedScalar& k0() const
189             {
190                 return k0_;
191             }
193             //- Return the value of epsilon0 which epsilon is not allowed to be
194             //  less than
195             const dimensionedScalar& epsilon0() const
196             {
197                 return epsilon0_;
198             }
200             //- Return the value of epsilonSmall which is added to epsilon when
201             //  calculating nut
202             const dimensionedScalar& epsilonSmall() const
203             {
204                 return epsilonSmall_;
205             }
207             //- Return the value of omega0 which epsilon is not allowed to be
208             //  less than
209             const dimensionedScalar& omega0() const
210             {
211                 return omega0_;
212             }
214             //- Return the value of omegaSmall which is added to epsilon when
215             //  calculating nut
216             const dimensionedScalar& omegaSmall() const
217             {
218                 return omegaSmall_;
219             }
221             //- Allow k0 to be changed
222             dimensionedScalar& k0()
223             {
224                 return k0_;
225             }
227             //- Allow epsilon0 to be changed
228             dimensionedScalar& epsilon0()
229             {
230                 return epsilon0_;
231             }
233             //- Allow epsilonSmall to be changed
234             dimensionedScalar& epsilonSmall()
235             {
236                 return epsilonSmall_;
237             }
239             //- Allow omega0 to be changed
240             dimensionedScalar& omega0()
241             {
242                 return omega0_;
243             }
245             //- Allow omegaSmall to be changed
246             dimensionedScalar& omegaSmall()
247             {
248                 return omegaSmall_;
249             }
251             //- Return the near wall distances
252             const nearWallDist& y() const
253             {
254                 return y_;
255             }
257             //- Calculate y+ at the edge of the laminar sublayer
258             scalar yPlusLam(const scalar kappa, const scalar E) const;
260             //- Const access to the coefficients dictionary
261             const dictionary& coeffDict() const
262             {
263                 return coeffDict_;
264             }
267         //- Return the turbulence viscosity
268         virtual tmp<volScalarField> mut() const = 0;
270         //- Return the effective viscosity
271         virtual tmp<volScalarField> muEff() const
272         {
273             return tmp<volScalarField>
274             (
275                 new volScalarField("muEff", mut() + mu())
276             );
277         }
279         //- Return the effective turbulent thermal diffusivity
280         virtual tmp<volScalarField> alphaEff() const = 0;
282         //- Return the turbulence kinetic energy
283         virtual tmp<volScalarField> k() const = 0;
285         //- Return the turbulence kinetic energy dissipation rate
286         virtual tmp<volScalarField> epsilon() const = 0;
288         //- Return the Reynolds stress tensor
289         virtual tmp<volSymmTensorField> R() const = 0;
291         //- Return the effective stress tensor including the laminar stress
292         virtual tmp<volSymmTensorField> devRhoReff() const = 0;
294         //- Return the source term for the momentum equation
295         virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const = 0;
297         //- Return yPlus for the given patch
298         virtual tmp<scalarField> yPlus
299         (
300             const label patchI,
301             const scalar Cmu
302         ) const;
304         //- Solve the turbulence equations and correct the turbulence viscosity
305         virtual void correct() = 0;
307         //- Read RASProperties dictionary
308         virtual bool read() = 0;
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 } // End namespace compressible
315 } // End namespace Foam
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 #endif
321 // ************************************************************************* //