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-2022 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  (
74  (
75  "kMin",
76  RASDict_,
78  small
79  )
80  ),
81 
82  epsilonMin_
83  (
85  (
86  "epsilonMin",
87  RASDict_,
88  kMin_.dimensions()/dimTime,
89  small
90  )
91  ),
92 
93  omegaMin_
94  (
96  (
97  "omegaMin",
98  RASDict_,
100  small
101  )
102  ),
103 
104  viscosityModel_
105  (
106  coeffDict_.found("viscosityModel")
108  (
109  coeffDict_,
110  viscosity,
111  U
112  )
113  : autoPtr<laminarModels::generalisedNewtonianViscosityModel>
114  (
115  new laminarModels::generalisedNewtonianViscosityModels::Newtonian
116  (
117  coeffDict_,
118  viscosity,
119  U
120  )
121  )
122  )
123 {
124  // Force the construction of the mesh deltaCoeffs which may be needed
125  // for the construction of the derived models and BCs
126  this->mesh_.deltaCoeffs();
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
131 
132 template<class BasicMomentumTransportModel>
135 (
136  const alphaField& alpha,
137  const rhoField& rho,
138  const volVectorField& U,
139  const surfaceScalarField& alphaRhoPhi,
140  const surfaceScalarField& phi,
141  const viscosity& viscosity
142 )
143 {
144  const IOdictionary modelDict
145  (
146  momentumTransportModel::readModelDict
147  (
148  U.db(),
149  alphaRhoPhi.group()
150  )
151  );
152 
153  const word modelType =
154  modelDict.subDict("RAS").lookupBackwardsCompatible<word>
155  (
156  {"model", "RASModel"}
157  );
158 
159  Info<< "Selecting RAS turbulence model " << modelType << endl;
160 
161  typename dictionaryConstructorTable::iterator cstrIter =
162  dictionaryConstructorTablePtr_->find(modelType);
163 
164  if (cstrIter == dictionaryConstructorTablePtr_->end())
165  {
167  << "Unknown RASModel type "
168  << modelType << nl << nl
169  << "Valid RASModel types:" << endl
170  << dictionaryConstructorTablePtr_->sortedToc()
171  << exit(FatalError);
172  }
173 
174  return autoPtr<RASModel>
175  (
176  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, viscosity)
177  );
178 }
179 
180 
181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 
183 template<class BasicMomentumTransportModel>
185 {
187  {
188  RASDict_ <<= this->subDict("RAS");
189  RASDict_.lookup("turbulence") >> turbulence_;
190 
191  coeffDict_ <<= RASDict_.optionalSubDict(type() + "Coeffs");
192 
193  kMin_.readIfPresent(RASDict_);
194  epsilonMin_.readIfPresent(RASDict_);
195  omegaMin_.readIfPresent(RASDict_);
196 
197  return true;
198  }
199  else
200  {
201  return false;
202  }
203 }
204 
205 
206 template<class BasicMomentumTransportModel>
208 {
209  viscosityModel_->correct();
211 }
212 
213 
214 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
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
U
Definition: pEqn.H:72
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
const dimensionSet dimless
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:184
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
const dimensionSet dimTime
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
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:135
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
phi
Definition: correctPhi.H:3
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:32
static const char nl
Definition: Ostream.H:260
const dimensionSet dimVelocity
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
messageStream Info
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:207
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:875