LESModel.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-2020 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 "LESModel.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 
70  LESDict_(this->subOrEmptyDict("LES")),
71  turbulence_(LESDict_.get<Switch>("turbulence")),
72  printCoeffs_(LESDict_.getOrDefault<Switch>("printCoeffs", false)),
73  coeffDict_(LESDict_.optionalSubDict(type + "Coeffs")),
74 
75  kMin_
76  (
78  (
79  "kMin",
80  LESDict_,
82  SMALL
83  )
84  ),
85 
86  epsilonMin_
87  (
89  (
90  "epsilonMin",
91  LESDict_,
92  kMin_.dimensions()/dimTime,
93  SMALL
94  )
95  ),
96 
97  omegaMin_
98  (
100  (
101  "omegaMin",
102  LESDict_,
104  SMALL
105  )
106  ),
107 
108  delta_
109  (
111  (
112  IOobject::groupName("delta", alphaRhoPhi.group()),
113  *this,
114  LESDict_
115  )
116  )
117 {
118  // Force the construction of the mesh deltaCoeffs which may be needed
119  // for the construction of the derived models and BCs
120  this->mesh_.deltaCoeffs();
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
125 
126 template<class BasicTurbulenceModel>
129 (
130  const alphaField& alpha,
131  const rhoField& rho,
132  const volVectorField& U,
133  const surfaceScalarField& alphaRhoPhi,
134  const surfaceScalarField& phi,
135  const transportModel& transport,
136  const word& propertiesName
137 )
138 {
139  const IOdictionary modelDict
140  (
141  IOobject
142  (
143  IOobject::groupName(propertiesName, alphaRhoPhi.group()),
144  U.time().constant(),
145  U.db(),
146  IOobject::MUST_READ_IF_MODIFIED,
147  IOobject::NO_WRITE,
148  false // Do not register
149  )
150  );
151 
152  const dictionary& dict = modelDict.subDict("LES");
153 
154  const word modelType
155  (
156  // "LESModel" for v2006 and earlier
157  dict.getCompat<word>("model", {{"LESModel", -2006}})
158  );
159 
160  Info<< "Selecting LES turbulence model " << modelType << endl;
161 
162  auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
163 
164  if (!cstrIter.found())
165  {
167  (
168  dict,
169  "LES model",
170  modelType,
171  *dictionaryConstructorTablePtr_
172  ) << exit(FatalIOError);
173  }
174 
175  return autoPtr<LESModel>
176  (
177  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
178  );
179 }
180 
181 
182 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183 
184 template<class BasicTurbulenceModel>
186 {
188  {
189  LESDict_ <<= this->subDict("LES");
190  LESDict_.readEntry("turbulence", turbulence_);
191 
192  coeffDict_ <<= LESDict_.optionalSubDict(type() + "Coeffs");
193 
194  delta_().read(LESDict_);
195 
196  kMin_.readIfPresent(LESDict_);
197 
198  return true;
199  }
200 
201  return false;
202 }
203 
204 
205 template<class BasicTurbulenceModel>
207 {
208  delta_().correct();
210 }
211 
212 
213 // ************************************************************************* //
Foam::LESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:107
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:104
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::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::LESModel::New
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:129
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
rho
rho
Definition: readInitialConditions.H:88
Foam::LESModel::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:206
Foam::LESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:108
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:406
Foam::LESModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:185
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::LESModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:34
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:528
LESModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::LESModel::LESModel
LESModel(const LESModel &)=delete
No copy construct.
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
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:121
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
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::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
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::LESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:106
Foam::GeometricField< vector, fvPatchField, volMesh >