LaunderSharmaKE.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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 "LaunderSharmaKE.H"
30#include "fvOptions.H"
31#include "bound.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace RASModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44{
45 return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
46}
47
48
49template<class BasicTurbulenceModel>
51{
52 return
53 scalar(1)
54 - 0.3*exp(-min(sqr(sqr(k_)/(this->nu()*epsilon_)), scalar(50)));
55}
56
57
58template<class BasicTurbulenceModel>
60{
61 this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
62 this->nut_.correctBoundaryConditions();
63 fv::options::New(this->mesh_).correct(this->nut_);
64
65 BasicTurbulenceModel::correctNut();
66}
67
68
69template<class BasicTurbulenceModel>
71{
73 (
75 (
76 k_,
77 dimVolume*this->rho_.dimensions()*k_.dimensions()
78 /dimTime
79 )
80 );
81}
82
83
84template<class BasicTurbulenceModel>
86{
88 (
90 (
91 epsilon_,
92 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
93 /dimTime
94 )
95 );
96}
97
98
99// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100
101template<class BasicTurbulenceModel>
103(
104 const alphaField& alpha,
105 const rhoField& rho,
106 const volVectorField& U,
107 const surfaceScalarField& alphaRhoPhi,
108 const surfaceScalarField& phi,
109 const transportModel& transport,
110 const word& propertiesName,
111 const word& type
112)
113:
114 eddyViscosity<RASModel<BasicTurbulenceModel>>
115 (
116 type,
117 alpha,
118 rho,
119 U,
120 alphaRhoPhi,
121 phi,
122 transport,
123 propertiesName
124 ),
125
126 Cmu_
127 (
128 dimensioned<scalar>::getOrAddToDict
129 (
130 "Cmu",
131 this->coeffDict_,
132 0.09
133 )
134 ),
135 C1_
136 (
137 dimensioned<scalar>::getOrAddToDict
138 (
139 "C1",
140 this->coeffDict_,
141 1.44
142 )
143 ),
144 C2_
145 (
146 dimensioned<scalar>::getOrAddToDict
147 (
148 "C2",
149 this->coeffDict_,
150 1.92
151 )
152 ),
153 C3_
154 (
155 dimensioned<scalar>::getOrAddToDict
156 (
157 "C3",
158 this->coeffDict_,
159 0
160 )
161 ),
162 sigmak_
163 (
164 dimensioned<scalar>::getOrAddToDict
165 (
166 "sigmak",
167 this->coeffDict_,
168 1.0
169 )
170 ),
171 sigmaEps_
172 (
173 dimensioned<scalar>::getOrAddToDict
174 (
175 "sigmaEps",
176 this->coeffDict_,
177 1.3
178 )
179 ),
180
181 k_
182 (
184 (
185 "k",
186 this->runTime_.timeName(),
187 this->mesh_,
188 IOobject::MUST_READ,
189 IOobject::AUTO_WRITE
190 ),
191 this->mesh_
192 ),
193
194 epsilon_
195 (
197 (
198 "epsilon",
199 this->runTime_.timeName(),
200 this->mesh_,
201 IOobject::MUST_READ,
202 IOobject::AUTO_WRITE
203 ),
204 this->mesh_
205 )
206{
207 bound(k_, this->kMin_);
208 bound(epsilon_, this->epsilonMin_);
209
210 if (type == typeName)
211 {
212 this->printCoeffs(type);
213 }
214}
215
216
217// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218
219template<class BasicTurbulenceModel>
221{
223 {
224 Cmu_.readIfPresent(this->coeffDict());
225 C1_.readIfPresent(this->coeffDict());
226 C2_.readIfPresent(this->coeffDict());
227 C3_.readIfPresent(this->coeffDict());
228 sigmak_.readIfPresent(this->coeffDict());
229 sigmaEps_.readIfPresent(this->coeffDict());
230
231 return true;
232 }
233
234 return false;
235}
236
237
238template<class BasicTurbulenceModel>
240{
241 if (!this->turbulence_)
242 {
243 return;
244 }
245
246 // Local references
247 const alphaField& alpha = this->alpha_;
248 const rhoField& rho = this->rho_;
249 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
250 const volVectorField& U = this->U_;
251 volScalarField& nut = this->nut_;
253
255
257
258 // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
259 // number model
260
261 volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
262 volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
263
265 volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
266 tgradU.clear();
267
268
269 // Dissipation equation
271 (
272 fvm::ddt(alpha, rho, epsilon_)
273 + fvm::div(alphaRhoPhi, epsilon_)
274 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
275 ==
276 C1_*alpha*rho*G*epsilon_/k_
277 - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
278 - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
279 + alpha*rho*E
280 + epsilonSource()
281 + fvOptions(alpha, rho, epsilon_)
282 );
283
284 epsEqn.ref().relax();
285 fvOptions.constrain(epsEqn.ref());
286 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
287 solve(epsEqn);
288 fvOptions.correct(epsilon_);
289 bound(epsilon_, this->epsilonMin_);
290
291
292 // Turbulent kinetic energy equation
294 (
295 fvm::ddt(alpha, rho, k_)
296 + fvm::div(alphaRhoPhi, k_)
297 - fvm::laplacian(alpha*rho*DkEff(), k_)
298 ==
299 alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
300 - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
301 + kSource()
302 + fvOptions(alpha, rho, k_)
303 );
304
305 kEqn.ref().relax();
306 fvOptions.constrain(kEqn.ref());
307 solve(kEqn);
308 fvOptions.correct(k_);
309 bound(k_, this->kMin_);
310
311 correctNut();
312}
313
314
315// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316
317} // End namespace RASModels
318} // End namespace Foam
319
320// ************************************************************************* //
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:55
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > f2() const
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
tmp< volScalarField > fMu() const
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
thermo correct()
scalar nut
zeroField divU
Definition: alphaSuSp.H:3
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
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 sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
volScalarField & nu
volScalarField & alpha
const dimensionedScalar & D
CEqn solve()