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