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-2020 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"
27 
28 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
29 
30 template<class BasicMomentumTransportModel>
32 {
33  if (printCoeffs_)
34  {
35  Info<< coeffDict_.dictName() << coeffDict_ << endl;
36  }
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicMomentumTransportModel>
44 (
45  const word& type,
46  const alphaField& alpha,
47  const rhoField& rho,
48  const volVectorField& U,
49  const surfaceScalarField& alphaRhoPhi,
50  const surfaceScalarField& phi,
51  const transportModel& transport
52 )
53 :
54  BasicMomentumTransportModel
55  (
56  type,
57  alpha,
58  rho,
59  U,
60  alphaRhoPhi,
61  phi,
62  transport
63  ),
64 
65  RASDict_(this->subOrEmptyDict("RAS")),
66  turbulence_(RASDict_.lookup("turbulence")),
67  printCoeffs_(RASDict_.lookupOrDefault<Switch>("printCoeffs", false)),
68  coeffDict_(RASDict_.optionalSubDict(type + "Coeffs")),
69 
70  kMin_
71  (
73  (
74  "kMin",
75  RASDict_,
77  small
78  )
79  ),
80 
81  epsilonMin_
82  (
84  (
85  "epsilonMin",
86  RASDict_,
87  kMin_.dimensions()/dimTime,
88  small
89  )
90  ),
91 
92  omegaMin_
93  (
95  (
96  "omegaMin",
97  RASDict_,
99  small
100  )
101  )
102 {
103  // Force the construction of the mesh deltaCoeffs which may be needed
104  // for the construction of the derived models and BCs
105  this->mesh_.deltaCoeffs();
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
110 
111 template<class BasicMomentumTransportModel>
114 (
115  const alphaField& alpha,
116  const rhoField& rho,
117  const volVectorField& U,
118  const surfaceScalarField& alphaRhoPhi,
119  const surfaceScalarField& phi,
120  const transportModel& transport
121 )
122 {
123  const IOdictionary modelDict
124  (
125  momentumTransportModel::readModelDict
126  (
127  U.db(),
128  alphaRhoPhi.group()
129  )
130  );
131 
132  const word modelType
133  (
134  modelDict.subDict("RAS").found("model")
135  ? modelDict.subDict("RAS").lookup("model")
136  : modelDict.subDict("RAS").lookup("RASModel")
137  );
138 
139  Info<< "Selecting RAS turbulence model " << modelType << endl;
140 
141  typename dictionaryConstructorTable::iterator cstrIter =
142  dictionaryConstructorTablePtr_->find(modelType);
143 
144  if (cstrIter == dictionaryConstructorTablePtr_->end())
145  {
147  << "Unknown RASModel type "
148  << modelType << nl << nl
149  << "Valid RASModel types:" << endl
150  << dictionaryConstructorTablePtr_->sortedToc()
151  << exit(FatalError);
152  }
153 
154  return autoPtr<RASModel>
155  (
156  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, transport)
157  );
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
163 template<class BasicMomentumTransportModel>
165 {
167  {
168  RASDict_ <<= this->subDict("RAS");
169  RASDict_.lookup("turbulence") >> turbulence_;
170 
171  coeffDict_ <<= RASDict_.optionalSubDict(type() + "Coeffs");
172 
173  kMin_.readIfPresent(RASDict_);
174  epsilonMin_.readIfPresent(RASDict_);
175  omegaMin_.readIfPresent(RASDict_);
176 
177  return true;
178  }
179  else
180  {
181  return false;
182  }
183 }
184 
185 
186 template<class BasicMomentumTransportModel>
188 {
190 }
191 
192 
193 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:144
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
static autoPtr< RASModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Return a reference to the selected RAS model.
Definition: RASModel.C:114
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:164
phi
Definition: pEqn.H:104
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
A class for handling words, derived from string.
Definition: word.H:59
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:31
static const char nl
Definition: Ostream.H:260
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
RASModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
Definition: RASModel.C:44
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:187
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
const dimensionSet dimVelocity