LESModel.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-2018 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 "LESModel.H"
27 
28 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
29 
30 template<class BasicTurbulenceModel>
32 {
33  if (printCoeffs_)
34  {
35  Info<< coeffDict_.dictName() << coeffDict_ << endl;
36  }
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
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  const word& propertiesName
53 )
54 :
55  BasicTurbulenceModel
56  (
57  type,
58  alpha,
59  rho,
60  U,
61  alphaRhoPhi,
62  phi,
63  transport,
64  propertiesName
65  ),
66 
67  LESDict_(this->subOrEmptyDict("LES")),
68  turbulence_(LESDict_.lookup("turbulence")),
69  printCoeffs_(LESDict_.lookupOrDefault<Switch>("printCoeffs", false)),
70  coeffDict_(LESDict_.optionalSubDict(type + "Coeffs")),
71 
72  kMin_
73  (
75  (
76  "kMin",
77  LESDict_,
79  small
80  )
81  ),
82 
83  epsilonMin_
84  (
86  (
87  "epsilonMin",
88  LESDict_,
89  kMin_.dimensions()/dimTime,
90  small
91  )
92  ),
93 
94  omegaMin_
95  (
97  (
98  "omegaMin",
99  LESDict_,
101  small
102  )
103  ),
104 
105  delta_
106  (
108  (
109  IOobject::groupName("delta", alphaRhoPhi.group()),
110  *this,
111  LESDict_
112  )
113  )
114 {
115  // Force the construction of the mesh deltaCoeffs which may be needed
116  // for the construction of the derived models and BCs
117  this->mesh_.deltaCoeffs();
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
122 
123 template<class BasicTurbulenceModel>
126 (
127  const alphaField& alpha,
128  const rhoField& rho,
129  const volVectorField& U,
130  const surfaceScalarField& alphaRhoPhi,
131  const surfaceScalarField& phi,
132  const transportModel& transport,
133  const word& propertiesName
134 )
135 {
136  // get model name, but do not register the dictionary
137  // otherwise it is registered in the database twice
138  const word modelType
139  (
141  (
142  IOobject
143  (
144  IOobject::groupName(propertiesName, alphaRhoPhi.group()),
145  U.time().constant(),
146  U.db(),
147  IOobject::MUST_READ_IF_MODIFIED,
148  IOobject::NO_WRITE,
149  false
150  )
151  ).subDict("LES").lookup("LESModel")
152  );
153 
154  Info<< "Selecting LES turbulence model " << modelType << endl;
155 
156  typename dictionaryConstructorTable::iterator cstrIter =
157  dictionaryConstructorTablePtr_->find(modelType);
158 
159  if (cstrIter == dictionaryConstructorTablePtr_->end())
160  {
162  << "Unknown LESModel type "
163  << modelType << nl << nl
164  << "Valid LESModel types:" << endl
165  << dictionaryConstructorTablePtr_->sortedToc()
166  << exit(FatalError);
167  }
168 
169  return autoPtr<LESModel>
170  (
171  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
172  );
173 }
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
178 template<class BasicTurbulenceModel>
180 {
182  {
183  LESDict_ <<= this->subDict("LES");
184  LESDict_.lookup("turbulence") >> turbulence_;
185 
186  coeffDict_ <<= LESDict_.optionalSubDict(type() + "Coeffs");
187 
188  delta_().read(LESDict_);
189 
190  kMin_.readIfPresent(LESDict_);
191 
192  return true;
193  }
194  else
195  {
196  return false;
197  }
198 }
199 
200 
201 template<class BasicTurbulenceModel>
203 {
204  delta_().correct();
206 }
207 
208 
209 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:111
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:176
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
type
Types of root.
Definition: Roots.H:52
surfaceScalarField & phi
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:256
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:109
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:31
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.
Definition: Switch.H:60
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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:692
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:110
static const char nl
Definition: Ostream.H:265
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:202
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/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:72
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:179
const Time & time() const
Return time.
Definition: IOobject.C:367
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
static autoPtr< LESModel > 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 LES model.
Definition: LESModel.C:126
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:361
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const dimensionSet dimVelocity