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-2020 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 "MomentumTransportModel.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 /*---------------------------------------------------------------------------*\
46  Class RASModel Declaration
47 \*---------------------------------------------------------------------------*/
48 
49 template<class BasicMomentumTransportModel>
50 class RASModel
51 :
52  public BasicMomentumTransportModel
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 BasicMomentumTransportModel::alphaField alphaField;
90  typedef typename BasicMomentumTransportModel::rhoField rhoField;
91  typedef typename BasicMomentumTransportModel::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  ),
113  (alpha, rho, U, alphaRhoPhi, phi, transport)
114  );
115 
116 
117  // Constructors
118 
119  //- Construct from components
120  RASModel
121  (
122  const word& type,
123  const alphaField& alpha,
124  const rhoField& rho,
125  const volVectorField& U,
126  const surfaceScalarField& alphaRhoPhi,
127  const surfaceScalarField& phi,
128  const transportModel& transport
129  );
130 
131  //- Disallow default bitwise copy construction
132  RASModel(const RASModel&) = delete;
133 
134 
135  // Selectors
136 
137  //- Return a reference to the selected RAS model
138  static autoPtr<RASModel> New
139  (
140  const alphaField& alpha,
141  const rhoField& rho,
142  const volVectorField& U,
143  const surfaceScalarField& alphaRhoPhi,
144  const surfaceScalarField& phi,
145  const transportModel& transport
146  );
147 
148 
149  //- Destructor
150  virtual ~RASModel()
151  {}
152 
153 
154  // Member Functions
155 
156  //- Read model coefficients if they have changed
157  virtual bool read();
158 
159 
160  // Access
161 
162  //- Return the lower allowable limit for k (default: small)
163  const dimensionedScalar& kMin() const
164  {
165  return kMin_;
166  }
167 
168  //- Return the lower allowable limit for epsilon (default: small)
169  const dimensionedScalar& epsilonMin() const
170  {
171  return epsilonMin_;
172  }
173 
174  //- Return the lower allowable limit for omega (default: small)
175  const dimensionedScalar& omegaMin() const
176  {
177  return omegaMin_;
178  }
179 
180  //- Allow kMin to be changed
182  {
183  return kMin_;
184  }
185 
186  //- Allow epsilonMin to be changed
188  {
189  return epsilonMin_;
190  }
191 
192  //- Allow omegaMin to be changed
194  {
195  return omegaMin_;
196  }
197 
198  //- Const access to the coefficients dictionary
199  virtual const dictionary& coeffDict() const
200  {
201  return coeffDict_;
202  }
203 
204 
205  //- Return the effective viscosity
206  virtual tmp<volScalarField> nuEff() const
207  {
208  return volScalarField::New
209  (
210  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
211  this->nut() + this->nu()
212  );
213  }
214 
215  //- Return the effective viscosity on patch
216  virtual tmp<scalarField> nuEff(const label patchi) const
217  {
218  return this->nut(patchi) + this->nu(patchi);
219  }
220 
221  //- Solve the turbulence equations and correct the turbulence viscosity
222  virtual void correct();
223 
224 
225  // Member Operators
226 
227  //- Disallow default bitwise assignment
228  void operator=(const RASModel&) = delete;
229 };
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #ifdef NoRepository
239  #include "RASModel.C"
240 #endif
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #endif
245 
246 // ************************************************************************* //
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.
static autoPtr< RASModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Return a reference to the selected RAS model.
Definition: RASModel.C:114
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
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
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
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:164
virtual ~RASModel()
Destructor.
Definition: RASModel.H:149
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
phi
Definition: pEqn.H:104
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:205
dictionary coeffDict_
Model coefficients dictionary.
Definition: RASModel.H:68
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:198
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
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:162
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:31
declareRunTimeSelectionTable(autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport),(alpha, rho, U, alphaRhoPhi, phi, transport))
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: RASModel.H:77
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
label patchi
U
Definition: pEqn.H:72
TypeName("RAS")
Runtime type information.
const dimensionedScalar & epsilonMin() const
Return the lower allowable limit for epsilon (default: small)
Definition: RASModel.H:168
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
RASModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
Definition: RASModel.C:44
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
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:187
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:174
scalar nut
Namespace for OpenFOAM.
Switch turbulence_
Turbulence on/off flag.
Definition: RASModel.H:62