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  //- Upper limit coefficient for nut
80 
81  //- Run-time selectable generalised Newtonian viscosity model
84 
85 
86  // Protected Member Functions
87 
88  //- Print model coefficients
89  virtual void printCoeffs(const word& type);
90 
91 
92 public:
93 
94  typedef typename BasicMomentumTransportModel::alphaField alphaField;
95  typedef typename BasicMomentumTransportModel::rhoField rhoField;
96 
97 
98  //- Runtime type information
99  TypeName("RAS");
100 
101 
102  // Declare run-time constructor selection table
103 
105  (
106  autoPtr,
107  RASModel,
108  dictionary,
109  (
110  const alphaField& alpha,
111  const rhoField& rho,
112  const volVectorField& U,
113  const surfaceScalarField& alphaRhoPhi,
114  const surfaceScalarField& phi,
115  const viscosity& viscosity
116  ),
117  (alpha, rho, U, alphaRhoPhi, phi, viscosity)
118  );
119 
120 
121  // Constructors
122 
123  //- Construct from components
124  RASModel
125  (
126  const word& type,
127  const alphaField& alpha,
128  const rhoField& rho,
129  const volVectorField& U,
130  const surfaceScalarField& alphaRhoPhi,
131  const surfaceScalarField& phi,
132  const viscosity& viscosity
133  );
134 
135  //- Disallow default bitwise copy construction
136  RASModel(const RASModel&) = delete;
137 
138 
139  // Selectors
140 
141  //- Return a reference to the selected RAS model
142  static autoPtr<RASModel> New
143  (
144  const alphaField& alpha,
145  const rhoField& rho,
146  const volVectorField& U,
147  const surfaceScalarField& alphaRhoPhi,
148  const surfaceScalarField& phi,
149  const viscosity& viscosity
150  );
151 
152 
153  //- Destructor
154  virtual ~RASModel()
155  {}
156 
157 
158  // Member Functions
159 
160  //- Read model coefficients if they have changed
161  virtual bool read();
162 
163  //- Const access to the coefficients dictionary
164  virtual const dictionary& coeffDict() const
165  {
166  return coeffDict_;
167  }
168 
169  //- Return the laminar viscosity
170  virtual tmp<volScalarField> nu() const
171  {
172  return viscosityModel_->nu();
173  }
174 
175  //- Return the laminar viscosity on patchi
176  virtual tmp<scalarField> nu(const label patchi) const
177  {
178  return viscosityModel_->nu(patchi);
179  }
180 
181  //- Return the effective viscosity
182  virtual tmp<volScalarField> nuEff() const
183  {
184  return volScalarField::New
185  (
186  this->groupName("nuEff"),
187  this->nut() + this->nu()
188  );
189  }
190 
191  //- Return the effective viscosity on patch
192  virtual tmp<scalarField> nuEff(const label patchi) const
193  {
194  return this->nut(patchi) + this->nu(patchi);
195  }
196 
197  //- Predict the turbulence transport coefficients if possible
198  // without solving turbulence transport model equations
199  virtual void predict()
200  {}
201 
202  //- Solve the turbulence equations and correct the turbulence viscosity
203  virtual void correct();
204 
205 
206  // Member Operators
207 
208  //- Disallow default bitwise assignment
209  void operator=(const RASModel&) = delete;
210 };
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #ifdef NoRepository
220  #include "RASModel.C"
221 #endif
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #endif
226 
227 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
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:93
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:124
virtual ~RASModel()
Destructor.
Definition: RASModel.H:153
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:163
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:181
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:195
dictionary RASDict_
RAS coefficients dictionary.
Definition: RASModel.H:63
dimensionedScalar nutMaxCoeff_
Upper limit coefficient for nut.
Definition: RASModel.H:78
autoPtr< laminarModels::generalisedNewtonianViscosityModel > viscosityModel_
Run-time selectable generalised Newtonian viscosity model.
Definition: RASModel.H:82
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: RASModel.H:169
virtual void predict()
Predict the turbulence transport coefficients if possible.
Definition: RASModel.H:198
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:173
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:94
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:162
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