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-2023 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  with support for generalised Newtonian viscosity models including
31  strain-rate dependency.
32 
33 SourceFiles
34  RASModel.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef RASModel_H
39 #define RASModel_H
40 
41 #include "momentumTransportModel.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class RASModel Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class BasicMomentumTransportModel>
54 class RASModel
55 :
56  public BasicMomentumTransportModel
57 {
58 
59 protected:
60 
61  // Protected data
62 
63  //- RAS coefficients dictionary
65 
66  //- Turbulence on/off flag
68 
69  //- Flag to print the model coeffs at run-time
71 
72  //- Model coefficients dictionary
74 
75  //- Lower limit of k
77 
78  //- Lower limit of epsilon
80 
81  //- Lower limit for omega
83 
84  //- Run-time selectable generalised Newtonian viscosity model
87 
88 
89  // Protected Member Functions
90 
91  //- Print model coefficients
92  virtual void printCoeffs(const word& type);
93 
94 
95 public:
96 
97  typedef typename BasicMomentumTransportModel::alphaField alphaField;
98  typedef typename BasicMomentumTransportModel::rhoField rhoField;
99 
100 
101  //- Runtime type information
102  TypeName("RAS");
103 
104 
105  // Declare run-time constructor selection table
106 
108  (
109  autoPtr,
110  RASModel,
111  dictionary,
112  (
113  const alphaField& alpha,
114  const rhoField& rho,
115  const volVectorField& U,
116  const surfaceScalarField& alphaRhoPhi,
117  const surfaceScalarField& phi,
118  const viscosity& viscosity
119  ),
120  (alpha, rho, U, alphaRhoPhi, phi, viscosity)
121  );
122 
123 
124  // Constructors
125 
126  //- Construct from components
127  RASModel
128  (
129  const word& type,
130  const alphaField& alpha,
131  const rhoField& rho,
132  const volVectorField& U,
133  const surfaceScalarField& alphaRhoPhi,
134  const surfaceScalarField& phi,
135  const viscosity& viscosity
136  );
137 
138  //- Disallow default bitwise copy construction
139  RASModel(const RASModel&) = delete;
140 
141 
142  // Selectors
143 
144  //- Return a reference to the selected RAS model
145  static autoPtr<RASModel> New
146  (
147  const alphaField& alpha,
148  const rhoField& rho,
149  const volVectorField& U,
150  const surfaceScalarField& alphaRhoPhi,
151  const surfaceScalarField& phi,
152  const viscosity& viscosity
153  );
154 
155 
156  //- Destructor
157  virtual ~RASModel()
158  {}
159 
160 
161  // Member Functions
162 
163  //- Read model coefficients if they have changed
164  virtual bool read();
165 
166  //- Return the lower allowable limit for k (default: small)
167  const dimensionedScalar& kMin() const
168  {
169  return kMin_;
170  }
171 
172  //- Return the lower allowable limit for epsilon (default: small)
173  const dimensionedScalar& epsilonMin() const
174  {
175  return epsilonMin_;
176  }
177 
178  //- Return the lower allowable limit for omega (default: small)
179  const dimensionedScalar& omegaMin() const
180  {
181  return omegaMin_;
182  }
183 
184  //- Allow kMin to be changed
186  {
187  return kMin_;
188  }
189 
190  //- Allow epsilonMin to be changed
192  {
193  return epsilonMin_;
194  }
195 
196  //- Allow omegaMin to be changed
198  {
199  return omegaMin_;
200  }
201 
202  //- Const access to the coefficients dictionary
203  virtual const dictionary& coeffDict() const
204  {
205  return coeffDict_;
206  }
207 
208  //- Return the laminar viscosity
209  virtual tmp<volScalarField> nu() const
210  {
211  return viscosityModel_->nu();
212  }
213 
214  //- Return the laminar viscosity on patchi
215  virtual tmp<scalarField> nu(const label patchi) const
216  {
217  return viscosityModel_->nu(patchi);
218  }
219 
220  //- Return the effective viscosity
221  virtual tmp<volScalarField> nuEff() const
222  {
223  return volScalarField::New
224  (
225  this->groupName("nuEff"),
226  this->nut() + this->nu()
227  );
228  }
229 
230  //- Return the effective viscosity on patch
231  virtual tmp<scalarField> nuEff(const label patchi) const
232  {
233  return this->nut(patchi) + this->nu(patchi);
234  }
235 
236  //- Predict the turbulence transport coefficients if possible
237  // without solving turbulence transport model equations
238  virtual void predict()
239  {}
240 
241  //- Solve the turbulence equations and correct the turbulence viscosity
242  virtual void correct();
243 
244 
245  // Member Operators
246 
247  //- Disallow default bitwise assignment
248  void operator=(const RASModel&) = delete;
249 };
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 } // End namespace Foam
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #ifdef NoRepository
259  #include "RASModel.C"
260 #endif
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 #endif
265 
266 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:96
dictionary coeffDict_
Model coefficients dictionary.
Definition: RASModel.H:72
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:32
static autoPtr< RASModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected RAS model.
Definition: RASModel.C:135
virtual ~RASModel()
Destructor.
Definition: RASModel.H:156
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:202
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:220
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:207
dictionary RASDict_
RAS coefficients dictionary.
Definition: RASModel.H:63
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: RASModel.H:81
const dimensionedScalar & omegaMin() const
Return the lower allowable limit for omega (default: small)
Definition: RASModel.H:178
autoPtr< laminarModels::generalisedNewtonianViscosityModel > viscosityModel_
Run-time selectable generalised Newtonian viscosity model.
Definition: RASModel.H:85
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: RASModel.H:208
virtual void predict()
Predict the turbulence transport coefficients if possible.
Definition: RASModel.H:237
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:78
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: small)
Definition: RASModel.H:166
const dimensionedScalar & epsilonMin() const
Return the lower allowable limit for epsilon (default: small)
Definition: RASModel.H:172
RASModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
Definition: RASModel.C:45
void operator=(const RASModel &)=delete
Disallow default bitwise assignment.
Switch turbulence_
Turbulence on/off flag.
Definition: RASModel.H:66
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:75
declareRunTimeSelectionTable(autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity))
TypeName("RAS")
Runtime type information.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: RASModel.H:69
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:184
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:97
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
const scalar nut
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488