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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33template<class BasicTurbulenceModel>
35{
36 if (printCoeffs_)
37 {
38 Info<< coeffDict_.dictName() << coeffDict_ << endl;
39 }
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45template<class BasicTurbulenceModel>
47(
48 const word& type,
49 const alphaField& alpha,
50 const rhoField& rho,
51 const volVectorField& U,
52 const surfaceScalarField& alphaRhoPhi,
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 LESDict_(this->subOrEmptyDict("LES")),
70 turbulence_(LESDict_.getOrDefault<Switch>("turbulence", true)),
71 printCoeffs_(LESDict_.getOrDefault<Switch>("printCoeffs", false)),
72 coeffDict_(LESDict_.optionalSubDict(type + "Coeffs")),
73 Ce_
74 (
75 dimensioned<scalar>::getOrAddToDict
76 (
77 "Ce",
78 LESDict_,
79 1.048
80 )
81 ),
82 kMin_
83 (
84 dimensioned<scalar>::getOrAddToDict
85 (
86 "kMin",
87 LESDict_,
89 SMALL
90 )
91 ),
92 epsilonMin_
93 (
94 dimensioned<scalar>::getOrAddToDict
95 (
96 "epsilonMin",
97 LESDict_,
98 kMin_.dimensions()/dimTime,
99 SMALL
100 )
101 ),
102 omegaMin_
103 (
104 dimensioned<scalar>::getOrAddToDict
105 (
106 "omegaMin",
107 LESDict_,
109 SMALL
110 )
111 ),
112 delta_
113 (
115 (
116 IOobject::groupName("delta", alphaRhoPhi.group()),
117 *this,
118 LESDict_
119 )
120 )
121{
122 // Force the construction of the mesh deltaCoeffs which may be needed
123 // for the construction of the derived models and BCs
124 this->mesh_.deltaCoeffs();
125}
126
127
128// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
129
130template<class BasicTurbulenceModel>
133(
134 const alphaField& alpha,
135 const rhoField& rho,
136 const volVectorField& U,
137 const surfaceScalarField& alphaRhoPhi,
138 const surfaceScalarField& phi,
139 const transportModel& transport,
140 const word& propertiesName
141)
142{
143 const IOdictionary modelDict
144 (
146 (
147 IOobject::groupName(propertiesName, alphaRhoPhi.group()),
148 U.time().constant(),
149 U.db(),
152 false // Do not register
153 )
154 );
155
156 const dictionary& dict = modelDict.subDict("LES");
157
158 const word modelType
159 (
160 // "LESModel" for v2006 and earlier
161 dict.getCompat<word>("model", {{"LESModel", -2006}})
162 );
163
164 Info<< "Selecting LES turbulence model " << modelType << endl;
165
166 auto* ctorPtr = dictionaryConstructorTable(modelType);
167
168 if (!ctorPtr)
169 {
171 (
172 dict,
173 "LES model",
174 modelType,
175 *dictionaryConstructorTablePtr_
176 ) << exit(FatalIOError);
177 }
178
179 return autoPtr<LESModel>
180 (
181 ctorPtr(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
182 );
183}
184
185
186// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187
188template<class BasicTurbulenceModel>
190{
191 if (BasicTurbulenceModel::read())
192 {
193 LESDict_ <<= this->subDict("LES");
194 LESDict_.readEntry("turbulence", turbulence_);
195
196 coeffDict_ <<= LESDict_.optionalSubDict(type() + "Coeffs");
197
198 delta_().read(LESDict_);
199
200 Ce_.readIfPresent(LESDict_);
201
202 kMin_.readIfPresent(LESDict_);
203
204 return true;
205 }
206
207 return false;
208}
209
210
211template<class BasicTurbulenceModel>
214{
216 (
218 (
219 IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
220 this->mesh_.time().timeName(),
221 this->mesh_
222 ),
223 this->Ce()*pow(this->k(), 1.5)/this->delta()
224 );
225}
226
227
228template<class BasicTurbulenceModel>
231{
232 const scalar betaStar = 0.09;
233 const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
234
236 (
238 (
239 IOobject::groupName("omega", this->alphaRhoPhi_.group()),
240 this->mesh_.time().timeName(),
241 this->mesh_
242 ),
243 this->epsilon()/(betaStar*(this->k() + k0))
244 );
245}
246
247
248template<class BasicTurbulenceModel>
250{
251 delta_().correct();
252 BasicTurbulenceModel::correct();
253}
254
255
256// ************************************************************************* //
label k
scalar delta
surfaceScalarField & phi
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:282
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:180
Templated abstract base class for LES SGS models.
Definition: LESModel.H:62
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: LESModel.C:213
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:109
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:110
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:249
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:111
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: LESModel.C:230
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:189
Abstract base class for LES deltas.
Definition: LESdelta.H:54
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
T getCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, enum keyType::option matchOpt=keyType::REGEX) const
Generic dimensioned Type class.
virtual void printCoeffs()
Print model coefficients.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
volScalarField & alpha
dictionary dict