RASModel.C
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  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "RASModel.H"
30 
31 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32 
33 template<class BasicTurbulenceModel>
35 {
36  if (printCoeffs_)
37  {
38  Info<< coeffDict_.dictName() << coeffDict_ << endl;
39  }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class BasicTurbulenceModel>
47 (
48  const word& type,
49  const alphaField& alpha,
50  const rhoField& rho,
51  const volVectorField& U,
52  const surfaceScalarField& alphaRhoPhi,
53  const surfaceScalarField& phi,
54  const transportModel& transport,
55  const word& propertiesName
56 )
57 :
58  BasicTurbulenceModel
59  (
60  type,
61  alpha,
62  rho,
63  U,
64  alphaRhoPhi,
65  phi,
66  transport,
67  propertiesName
68  ),
69  RASDict_(this->subOrEmptyDict("RAS")),
70  turbulence_(RASDict_.getOrDefault<Switch>("turbulence", true)),
71  printCoeffs_(RASDict_.getOrDefault<Switch>("printCoeffs", false)),
72  coeffDict_(RASDict_.optionalSubDict(type + "Coeffs")),
73  kMin_
74  (
76  (
77  "kMin",
78  RASDict_,
80  SMALL
81  )
82  ),
83  epsilonMin_
84  (
86  (
87  "epsilonMin",
88  RASDict_,
89  kMin_.dimensions()/dimTime,
90  SMALL
91  )
92  ),
93  omegaMin_
94  (
96  (
97  "omegaMin",
98  RASDict_,
100  SMALL
101  )
102  )
103 {
104  // Force the construction of the mesh deltaCoeffs which may be needed
105  // for the construction of the derived models and BCs
106  this->mesh_.deltaCoeffs();
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
111 
112 template<class BasicTurbulenceModel>
115 (
116  const alphaField& alpha,
117  const rhoField& rho,
118  const volVectorField& U,
119  const surfaceScalarField& alphaRhoPhi,
120  const surfaceScalarField& phi,
121  const transportModel& transport,
122  const word& propertiesName
123 )
124 {
125  const IOdictionary modelDict
126  (
127  IOobject
128  (
129  IOobject::groupName(propertiesName, alphaRhoPhi.group()),
130  U.time().constant(),
131  U.db(),
132  IOobject::MUST_READ_IF_MODIFIED,
133  IOobject::NO_WRITE,
134  false // Do not register
135  )
136  );
137 
138  const dictionary& dict = modelDict.subDict("RAS");
139 
140  const word modelType
141  (
142  // "RASModel" for v2006 and earlier
143  dict.getCompat<word>("model", {{"RASModel", -2006}})
144  );
145 
146  Info<< "Selecting RAS turbulence model " << modelType << endl;
147 
148  auto* ctorPtr = dictionaryConstructorTable(modelType);
149 
150  if (!ctorPtr)
151  {
153  (
154  dict,
155  "RAS model",
156  modelType,
157  *dictionaryConstructorTablePtr_
158  ) << exit(FatalIOError);
159  }
160 
161  return autoPtr<RASModel>
162  (
163  ctorPtr(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
164  );
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
169 
170 template<class BasicTurbulenceModel>
172 {
174  {
175  RASDict_ <<= this->subDict("RAS");
176  RASDict_.readEntry("turbulence", turbulence_);
177 
178  coeffDict_ <<= RASDict_.optionalSubDict(type() + "Coeffs");
179 
180  kMin_.readIfPresent(RASDict_);
181  epsilonMin_.readIfPresent(RASDict_);
182  omegaMin_.readIfPresent(RASDict_);
183 
184  return true;
185  }
186 
187  return false;
188 }
189 
190 template<class BasicTurbulenceModel>
193 {
194  const scalar Cmu = 0.09;
195 
197  (
198  IOobject
199  (
200  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
201  this->mesh_.time().timeName(),
202  this->mesh_
203  ),
204  Cmu*this->k()*this->omega()
205  );
206 }
207 
208 
209 template<class BasicTurbulenceModel>
212 {
213  const scalar betaStar = 0.09;
214  const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
215 
217  (
218  IOobject
219  (
220  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
221  this->mesh_.time().timeName(),
222  this->mesh_
223  ),
224  this->epsilon()/(betaStar*(this->k() + k0))
225  );
226 }
227 
228 
229 template<class BasicTurbulenceModel>
231 {
233 }
234 
235 
236 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
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::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::RASModel::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: RASModel.C:211
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModel::RASModel
RASModel(const RASModel &)=delete
No copy construct.
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rho
rho
Definition: readInitialConditions.H:88
Foam::ThermalDiffusivity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: ThermalDiffusivity.H:58
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
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:115
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
correct
fvOptions correct(rho)
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
RASModel.H
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
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:230
Foam::RASModel::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: RASModel.C:192
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::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::RASModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::RASModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:171
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189