RASModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::RASModel
28 
29 Description
30  Templated abstract base class for RAS turbulence models
31 
32 SourceFiles
33  RASModel.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef RASModel_H
38 #define RASModel_H
39 
40 #include "TurbulenceModel.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class RASModel Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 template<class BasicTurbulenceModel>
52 class RASModel
53 :
54  public BasicTurbulenceModel
55 {
56 
57 protected:
58 
59  // Protected data
60 
61  //- RAS coefficients dictionary
63 
64  //- Turbulence on/off flag
66 
67  //- Flag to print the model coeffs at run-time
69 
70  //- Model coefficients dictionary
72 
73  //- Lower limit of k
75 
76  //- Lower limit of epsilon
78 
79  //- Lower limit for omega
81 
82 
83  // Protected Member Functions
84 
85  //- Print model coefficients
86  virtual void printCoeffs(const word& type);
87 
88  //- No copy construct
89  RASModel(const RASModel&) = delete;
90 
91  //- No copy assignment
92  void operator=(const RASModel&) = delete;
93 
94 
95 public:
96 
97  typedef typename BasicTurbulenceModel::alphaField alphaField;
98  typedef typename BasicTurbulenceModel::rhoField rhoField;
99  typedef typename BasicTurbulenceModel::transportModel transportModel;
100 
101 
102  //- Runtime type information
103  TypeName("RAS");
104 
105 
106  // Declare run-time constructor selection table
107 
109  (
110  autoPtr,
111  RASModel,
112  dictionary,
113  (
114  const alphaField& alpha,
115  const rhoField& rho,
116  const volVectorField& U,
117  const surfaceScalarField& alphaRhoPhi,
118  const surfaceScalarField& phi,
119  const transportModel& transport,
120  const word& propertiesName
121  ),
122  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
123  );
124 
125 
126  // Constructors
127 
128  //- Construct from components
129  RASModel
130  (
131  const word& type,
132  const alphaField& alpha,
133  const rhoField& rho,
134  const volVectorField& U,
135  const surfaceScalarField& alphaRhoPhi,
136  const surfaceScalarField& phi,
137  const transportModel& transport,
138  const word& propertiesName
139  );
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 transportModel& transport,
153  const word& propertiesName = turbulenceModel::propertiesName
154  );
155 
156 
157  //- Destructor
158  virtual ~RASModel() = default;
159 
160 
161  // Member Functions
162 
163  //- Read model coefficients if they have changed
164  virtual bool read();
165 
166 
167  // Access
168 
169  //- Return the lower allowable limit for k (default: SMALL)
170  const dimensionedScalar& kMin() const
171  {
172  return kMin_;
173  }
174 
175  //- Return the lower allowable limit for epsilon (default: SMALL)
176  const dimensionedScalar& epsilonMin() const
177  {
178  return epsilonMin_;
179  }
180 
181  //- Return the lower allowable limit for omega (default: SMALL)
182  const dimensionedScalar& omegaMin() const
183  {
184  return omegaMin_;
185  }
186 
187  //- Allow kMin to be changed
189  {
190  return kMin_;
191  }
192 
193  //- Allow epsilonMin to be changed
195  {
196  return epsilonMin_;
197  }
198 
199  //- Allow omegaMin to be changed
201  {
202  return omegaMin_;
203  }
204 
205  //- Const access to the coefficients dictionary
206  virtual const dictionary& coeffDict() const
207  {
208  return coeffDict_;
209  }
210 
211 
212  //- Return the effective viscosity
213  virtual tmp<volScalarField> nuEff() const
214  {
215  return tmp<volScalarField>
216  (
217  new volScalarField
218  (
219  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
220  this->nut() + this->nu()
221  )
222  );
223  }
224 
225  //- Return the effective viscosity on patch
226  virtual tmp<scalarField> nuEff(const label patchi) const
227  {
228  return this->nut(patchi) + this->nu(patchi);
229  }
230 
231  //- Solve the turbulence equations and correct the turbulence viscosity
232  virtual void correct();
233 };
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #ifdef NoRepository
243  #include "RASModel.C"
244 #endif
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 #endif
249 
250 // ************************************************************************* //
Foam::RASModel::coeffDict_
dictionary coeffDict_
Model coefficients dictionary.
Definition: RASModel.H:70
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::RASModel::omegaMin
dimensionedScalar & omegaMin()
Allow omegaMin to be changed.
Definition: RASModel.H:199
Foam::RASModel::TypeName
TypeName("RAS")
Runtime type information.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::RASModel::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:212
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModel::printCoeffs_
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: RASModel.H:67
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
RASModel.C
Foam::RASModel::RASModel
RASModel(const RASModel &)=delete
No copy construct.
Foam::RASModel::epsilonMin_
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:76
Foam::RASModel::omegaMin_
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: RASModel.H:79
rho
rho
Definition: readInitialConditions.H:88
Foam::RASModel::~RASModel
virtual ~RASModel()=default
Destructor.
Foam::ThermalDiffusivity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: ThermalDiffusivity.H:58
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::RASModel::New
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:119
Foam::RASModel::epsilonMin
const dimensionedScalar & epsilonMin() const
Return the lower allowable limit for epsilon (default: SMALL)
Definition: RASModel.H:175
Foam::RASModel::RASDict_
dictionary RASDict_
RAS coefficients dictionary.
Definition: RASModel.H:61
Foam::RASModel::operator=
void operator=(const RASModel &)=delete
No copy assignment.
Foam::RASModel::epsilonMin
dimensionedScalar & epsilonMin()
Allow epsilonMin to be changed.
Definition: RASModel.H:193
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:97
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:96
Foam::RASModel::nuEff
virtual tmp< scalarField > nuEff(const label patchi) const
Return the effective viscosity on patch.
Definition: RASModel.H:225
Foam::RASModel::declareRunTimeSelectionTable
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))
Foam::ThermalDiffusivity::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: ThermalDiffusivity.H:60
Foam::RASModel::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:196
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:98
Foam::RASModel::coeffDict
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:205
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::RASModel::omegaMin
const dimensionedScalar & omegaMin() const
Return the lower allowable limit for omega (default: SMALL)
Definition: RASModel.H:181
Foam::RASModel::turbulence_
Switch turbulence_
Turbulence on/off flag.
Definition: RASModel.H:64
Foam::RASModel::kMin_
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:73
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
TurbulenceModel.H
Foam::RASModel::kMin
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: SMALL)
Definition: RASModel.H:169
Foam::RASModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::RASModel::kMin
dimensionedScalar & kMin()
Allow kMin to be changed.
Definition: RASModel.H:187
Foam::RASModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:175