RASModel.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "RASModel.H"
28 
29 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
30 
31 template<class BasicMomentumTransportModel>
33 {
34  if (printCoeffs_)
35  {
36  Info<< coeffDict_.dictName() << coeffDict_ << endl;
37  }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 template<class BasicMomentumTransportModel>
45 (
46  const word& type,
47  const alphaField& alpha,
48  const rhoField& rho,
49  const volVectorField& U,
50  const surfaceScalarField& alphaRhoPhi,
51  const surfaceScalarField& phi,
52  const viscosity& viscosity
53 )
54 :
55  BasicMomentumTransportModel
56  (
57  type,
58  alpha,
59  rho,
60  U,
61  alphaRhoPhi,
62  phi,
63  viscosity
64  ),
65 
66  RASDict_(this->subOrEmptyDict("RAS")),
67  turbulence_(RASDict_.lookup("turbulence")),
68  printCoeffs_(RASDict_.lookupOrDefault<Switch>("printCoeffs", false)),
69  coeffDict_(RASDict_.optionalSubDict(type + "Coeffs")),
70 
71  kMin_
72  (
73  dimensioned<scalar>::lookupOrAddToDict
74  (
75  "kMin",
76  RASDict_,
78  small
79  )
80  ),
81 
82  nutMaxCoeff_
83  (
84  dimensioned<scalar>::lookupOrAddToDict
85  (
86  "nutMaxCoeff",
87  RASDict_,
88  dimless,
89  1e5
90  )
91  ),
92 
93  viscosityModel_
94  (
95  coeffDict_.found("viscosityModel")
96  ? laminarModels::generalisedNewtonianViscosityModel::New
97  (
98  coeffDict_,
99  viscosity,
100  U
101  )
102  : autoPtr<laminarModels::generalisedNewtonianViscosityModel>
103  (
104  new laminarModels::generalisedNewtonianViscosityModels::Newtonian
105  (
106  coeffDict_,
107  viscosity,
108  U
109  )
110  )
111  )
112 {
113  // Force the construction of the mesh deltaCoeffs which may be needed
114  // for the construction of the derived models and BCs
115  this->mesh_.deltaCoeffs();
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
120 
121 template<class BasicMomentumTransportModel>
124 (
125  const alphaField& alpha,
126  const rhoField& rho,
127  const volVectorField& U,
128  const surfaceScalarField& alphaRhoPhi,
129  const surfaceScalarField& phi,
130  const viscosity& viscosity
131 )
132 {
133  const IOdictionary modelDict
134  (
136  (
137  U.db(),
138  alphaRhoPhi.group()
139  )
140  );
141 
142  const word modelType =
143  modelDict.subDict("RAS").lookupBackwardsCompatible<word>
144  (
145  {"model", "RASModel"}
146  );
147 
148  Info<< "Selecting RAS turbulence model " << modelType << endl;
149 
150  typename dictionaryConstructorTable::iterator cstrIter =
151  dictionaryConstructorTablePtr_->find(modelType);
152 
153  if (cstrIter == dictionaryConstructorTablePtr_->end())
154  {
156  << "Unknown RASModel type "
157  << modelType << nl << nl
158  << "Valid RASModel types:" << endl
159  << dictionaryConstructorTablePtr_->sortedToc()
160  << exit(FatalError);
161  }
162 
163  return autoPtr<RASModel>
164  (
165  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, viscosity)
166  );
167 }
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
172 template<class BasicMomentumTransportModel>
174 {
176  {
177  RASDict_ <<= this->subDict("RAS");
178  RASDict_.lookup("turbulence") >> turbulence_;
179 
180  coeffDict_ <<= RASDict_.optionalSubDict(type() + "Coeffs");
181 
182  kMin_.readIfPresent(RASDict_);
183  nutMaxCoeff_.readIfPresent(RASDict_);
184 
185  return true;
186  }
187  else
188  {
189  return false;
190  }
191 }
192 
193 
194 template<class BasicMomentumTransportModel>
196 {
197  viscosityModel_->correct();
199 }
200 
201 
202 // ************************************************************************* //
bool found
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:93
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 void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:195
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
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
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:721
Generic dimensioned Type class.
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
messageStream Info
const dimensionSet dimVelocity
error FatalError
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488