RASModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::RASModel
26 
27 Description
28  Templated abstract base class for RAS turbulence models
29 
30 SourceFiles
31  RASModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef RASModel_H
36 #define RASModel_H
37 
38 #include "TurbulenceModel.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 /*---------------------------------------------------------------------------*\
46  Class RASModel Declaration
47 \*---------------------------------------------------------------------------*/
48 
49 template<class BasicTurbulenceModel>
50 class RASModel
51 :
52  public BasicTurbulenceModel
53 {
54 
55 protected:
56 
57  // Protected data
58 
59  //- RAS coefficients dictionary
61 
62  //- Turbulence on/off flag
64 
65  //- Flag to print the model coeffs at run-time
67 
68  //- Model coefficients dictionary
70 
71  //- Lower limit of k
73 
74  //- Lower limit of epsilon
76 
77  //- Lower limit for omega
79 
80 
81  // Protected Member Functions
82 
83  //- Print model coefficients
84  virtual void printCoeffs(const word& type);
85 
86 
87 public:
88 
89  typedef typename BasicTurbulenceModel::alphaField alphaField;
90  typedef typename BasicTurbulenceModel::rhoField rhoField;
91  typedef typename BasicTurbulenceModel::transportModel transportModel;
92 
93 
94  //- Runtime type information
95  TypeName("RAS");
96 
97 
98  // Declare run-time constructor selection table
99 
101  (
102  autoPtr,
103  RASModel,
104  dictionary,
105  (
106  const alphaField& alpha,
107  const rhoField& rho,
108  const volVectorField& U,
109  const surfaceScalarField& alphaRhoPhi,
110  const surfaceScalarField& phi,
111  const transportModel& transport,
112  const word& propertiesName
113  ),
114  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
115  );
116 
117 
118  // Constructors
119 
120  //- Construct from components
121  RASModel
122  (
123  const word& type,
124  const alphaField& alpha,
125  const rhoField& rho,
126  const volVectorField& U,
127  const surfaceScalarField& alphaRhoPhi,
128  const surfaceScalarField& phi,
129  const transportModel& transport,
130  const word& propertiesName
131  );
132 
133  //- Disallow default bitwise copy construction
134  RASModel(const RASModel&) = delete;
135 
136 
137  // Selectors
138 
139  //- Return a reference to the selected RAS model
140  static autoPtr<RASModel> New
141  (
142  const alphaField& alpha,
143  const rhoField& rho,
144  const volVectorField& U,
145  const surfaceScalarField& alphaRhoPhi,
146  const surfaceScalarField& phi,
147  const transportModel& transport,
148  const word& propertiesName = turbulenceModel::propertiesName
149  );
150 
151 
152  //- Destructor
153  virtual ~RASModel()
154  {}
155 
156 
157  // Member Functions
158 
159  //- Read model coefficients if they have changed
160  virtual bool read();
161 
162 
163  // Access
164 
165  //- Return the lower allowable limit for k (default: small)
166  const dimensionedScalar& kMin() const
167  {
168  return kMin_;
169  }
170 
171  //- Return the lower allowable limit for epsilon (default: small)
172  const dimensionedScalar& epsilonMin() const
173  {
174  return epsilonMin_;
175  }
176 
177  //- Return the lower allowable limit for omega (default: small)
178  const dimensionedScalar& omegaMin() const
179  {
180  return omegaMin_;
181  }
182 
183  //- Allow kMin to be changed
185  {
186  return kMin_;
187  }
188 
189  //- Allow epsilonMin to be changed
191  {
192  return epsilonMin_;
193  }
194 
195  //- Allow omegaMin to be changed
197  {
198  return omegaMin_;
199  }
200 
201  //- Const access to the coefficients dictionary
202  virtual const dictionary& coeffDict() const
203  {
204  return coeffDict_;
205  }
206 
207 
208  //- Return the effective viscosity
209  virtual tmp<volScalarField> nuEff() const
210  {
211  return volScalarField::New
212  (
213  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
214  this->nut() + this->nu()
215  );
216  }
217 
218  //- Return the effective viscosity on patch
219  virtual tmp<scalarField> nuEff(const label patchi) const
220  {
221  return this->nut(patchi) + this->nu(patchi);
222  }
223 
224  //- Solve the turbulence equations and correct the turbulence viscosity
225  virtual void correct();
226 
227 
228  // Member Operators
229 
230  //- Disallow default bitwise assignment
231  void operator=(const RASModel&) = delete;
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace Foam
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #ifdef NoRepository
242  #include "RASModel.C"
243 #endif
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 #endif
248 
249 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void operator=(const RASModel &)=delete
Disallow default bitwise assignment.
surfaceScalarField & phi
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:71
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dictionary RASDict_
RAS coefficients dictionary.
Definition: RASModel.H:59
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:89
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
declareRunTimeSelectionTable(autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:169
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:90
virtual ~RASModel()
Destructor.
Definition: RASModel.H:152
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:208
dictionary coeffDict_
Model coefficients dictionary.
Definition: RASModel.H:68
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:201
static const word propertiesName
Default name of the turbulence properties dictionary.
RASModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
Definition: RASModel.C:44
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: small)
Definition: RASModel.H:165
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:31
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: RASModel.H:77
label patchi
U
Definition: pEqn.H:72
TypeName("RAS")
Runtime type information.
static autoPtr< RASModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected RAS model.
Definition: RASModel.C:116
const dimensionedScalar & epsilonMin() const
Return the lower allowable limit for epsilon (default: small)
Definition: RASModel.H:171
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: RASModel.H:65
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
volScalarField & nu
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:192
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:74
const dimensionedScalar & omegaMin() const
Return the lower allowable limit for omega (default: small)
Definition: RASModel.H:177
scalar nut
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:88
Switch turbulence_
Turbulence on/off flag.
Definition: RASModel.H:62