BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / RASModel / RASModel.H
blob60fc04c2862465d5138f78f1db7592d204859909
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::incompressible::RASModels
27 Description
28     Namespace for incompressible RAS turbulence models.
30 Class
31     Foam::incompressible::RASModel
33 Description
34     Abstract base class for incompressible turbulence models.
36 SourceFiles
37     RASModel.C
39 \*---------------------------------------------------------------------------*/
41 #ifndef RASModel_H
42 #define RASModel_H
44 #include "incompressible/turbulenceModel/turbulenceModel.H"
45 #include "volFields.H"
46 #include "surfaceFields.H"
47 #include "nearWallDist.H"
48 #include "fvm.H"
49 #include "fvc.H"
50 #include "fvMatrices.H"
51 #include "incompressible/transportModel/transportModel.H"
52 #include "IOdictionary.H"
53 #include "Switch.H"
54 #include "bound.H"
55 #include "autoPtr.H"
56 #include "runTimeSelectionTables.H"
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 namespace Foam
62 namespace incompressible
65 /*---------------------------------------------------------------------------*\
66                            Class RASModel Declaration
67 \*---------------------------------------------------------------------------*/
69 class RASModel
71     public turbulenceModel,
72     public IOdictionary
75 protected:
77     // Protected data
79         //- Turbulence on/off flag
80         Switch turbulence_;
82         //- Flag to print the model coeffs at run-time
83         Switch printCoeffs_;
85         //- Model coefficients dictionary
86         dictionary coeffDict_;
88         //- Lower limit of k
89         dimensionedScalar kMin_;
91         //- Lower limit of epsilon
92         dimensionedScalar epsilonMin_;
94         //- Lower limit for omega
95         dimensionedScalar omegaMin_;
97         //- Near wall distance boundary field
98         nearWallDist y_;
101     // Protected Member Functions
103         //- Print model coefficients
104         virtual void printCoeffs();
107 private:
109     // Private Member Functions
111         //- Disallow default bitwise copy construct
112         RASModel(const RASModel&);
114         //- Disallow default bitwise assignment
115         void operator=(const RASModel&);
118 public:
120     //- Runtime type information
121     TypeName("RASModel");
124     // Declare run-time constructor selection table
126         declareRunTimeSelectionTable
127         (
128             autoPtr,
129             RASModel,
130             dictionary,
131             (
132                 const volVectorField& U,
133                 const surfaceScalarField& phi,
134                 transportModel& transport,
135                 const word& turbulenceModelName
136             ),
137             (U, phi, transport, turbulenceModelName)
138         );
141     // Constructors
143         //- Construct from components
144         RASModel
145         (
146             const word& type,
147             const volVectorField& U,
148             const surfaceScalarField& phi,
149             transportModel& transport,
150             const word& turbulenceModelName = turbulenceModel::typeName
151         );
154     // Selectors
156         //- Return a reference to the selected RAS model
157         static autoPtr<RASModel> New
158         (
159             const volVectorField& U,
160             const surfaceScalarField& phi,
161             transportModel& transport,
162             const word& turbulenceModelName = turbulenceModel::typeName
163         );
166     //- Destructor
167     virtual ~RASModel()
168     {}
171     // Member Functions
173         // Access
175             //- Return the lower allowable limit for k (default: SMALL)
176             const dimensionedScalar& kMin() const
177             {
178                 return kMin_;
179             }
181             //- Return the lower allowable limit for epsilon (default: SMALL)
182             const dimensionedScalar& epsilonMin() const
183             {
184                 return epsilonMin_;
185             }
187             //- Return the lower allowable limit for omega (default: SMALL)
188             const dimensionedScalar& omegaMin() const
189             {
190                 return omegaMin_;
191             }
193             //- Allow kMin to be changed
194             dimensionedScalar& kMin()
195             {
196                 return kMin_;
197             }
199             //- Allow epsilonMin to be changed
200             dimensionedScalar& epsilonMin()
201             {
202                 return epsilonMin_;
203             }
205             //- Allow omegaMin to be changed
206             dimensionedScalar& omegaMin()
207             {
208                 return omegaMin_;
209             }
211             //- Return the near wall distances
212             const nearWallDist& y() const
213             {
214                 return y_;
215             }
217             //- Calculate y+ at the edge of the laminar sublayer
218             scalar yPlusLam(const scalar kappa, const scalar E) const;
220             //- Const access to the coefficients dictionary
221             const dictionary& coeffDict() const
222             {
223                 return coeffDict_;
224             }
227         //- Return the effective viscosity
228         virtual tmp<volScalarField> nuEff() const
229         {
230             return tmp<volScalarField>
231             (
232                 new volScalarField("nuEff", nut() + nu())
233             );
234         }
236         //- Return yPlus for the given patch
237         virtual tmp<scalarField> yPlus
238         (
239             const label patchI,
240             const scalar Cmu
241         ) const;
243         //- Solve the turbulence equations and correct the turbulence viscosity
244         virtual void correct();
246         //- Read RASProperties dictionary
247         virtual bool read();
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 } // End namespace incompressible
254 } // End namespace Foam
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 #endif
260 // ************************************************************************* //