EBRSM.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) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "EBRSM.H"
29#include "fvOptions.H"
30#include "bound.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace RASModels
37{
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
42tmp<volScalarField> EBRSM<BasicTurbulenceModel>::calcL() const
43{
44 // (M:Eq. C.13)
45 return
46 Cl_*max
47 (
48 pow(k_, 1.5)/epsilon_,
49 Ceta_*pow025
50 (
51 pow3
52 (
53 max
54 (
55 this->nu(),
56 dimensionedScalar(this->nu()().dimensions(), Zero)
57 )
58 )/epsilon_
59 )
60 );
61}
62
63
64template<class BasicTurbulenceModel>
65tmp<volVectorField> EBRSM<BasicTurbulenceModel>::calcN() const
66{
67 const volVectorField gradf(fvc::grad(f_));
68
69 // (M:Eq. C.9)
70 return gradf/max(mag(gradf), dimensionedScalar(dimless/dimLength, SMALL));
71}
72
73
74template<class BasicTurbulenceModel>
75tmp<volScalarField> EBRSM<BasicTurbulenceModel>::calcTau() const
76{
77 // (M:Eq. C.12)
78 return
79 max
80 (
81 k_/epsilon_,
82 Ct_*sqrt
83 (
84 max
85 (
86 this->nu(),
87 dimensionedScalar(this->nu()().dimensions(), Zero)
88 )/epsilon_
89 )
90 );
91}
92
93
94template<class BasicTurbulenceModel>
95tmp<volSymmTensorField> EBRSM<BasicTurbulenceModel>::D
96(
97 const volScalarField& tau,
98 const dimensionedScalar& sigma
99) const
100{
101 // (M:Eq. C.10, C.14)
102 return (Cmu_/sigma*tau)*this->R_ + this->nu()*I;
103}
104
105
106template<class BasicTurbulenceModel>
107tmp<volScalarField> EBRSM<BasicTurbulenceModel>::D
108(
109 const dimensionedScalar& sigma
110) const
111{
112 // (LM:p. 2)
113 return this->nut_/sigma + this->nu();
114}
115
116
117template<class BasicTurbulenceModel>
118tmp<volScalarField> EBRSM<BasicTurbulenceModel>::Ceps1Prime
119(
120 const volScalarField& G
121) const
122{
123 // (M:Eq. C.15)
124 return Ceps1_*(scalar(1) + A1_*(scalar(1) - pow3(f_))*G/this->epsilon_);
125}
126
127
128template<class BasicTurbulenceModel>
129void EBRSM<BasicTurbulenceModel>::correctNut()
130{
131 this->nut_ = Cmu_*k_*calcTau();
132 this->nut_.correctBoundaryConditions();
133 fv::options::New(this->mesh_).correct(this->nut_);
134
135 BasicTurbulenceModel::correctNut();
136}
137
138
139// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140
141template<class BasicTurbulenceModel>
143(
144 const alphaField& alpha,
145 const rhoField& rho,
146 const volVectorField& U,
147 const surfaceScalarField& alphaRhoPhi,
148 const surfaceScalarField& phi,
149 const transportModel& transport,
150 const word& propertiesName,
151 const word& type
152)
153:
154 ReynoldsStress<RASModel<BasicTurbulenceModel>>
155 (
156 type,
157 alpha,
158 rho,
159 U,
160 alphaRhoPhi,
161 phi,
162 transport,
163 propertiesName
164 ),
165
166 simpleGradientDiffusion_
167 (
168 Switch::getOrAddToDict
169 (
170 "simpleGradientDiffusion",
171 this->coeffDict_,
172 false
173 )
174 ),
175 g1_
176 (
177 dimensioned<scalar>::getOrAddToDict
178 (
179 "g1",
180 this->coeffDict_,
181 3.4
182 )
183 ),
184 g1star_
185 (
186 dimensioned<scalar>::getOrAddToDict
187 (
188 "g1star",
189 this->coeffDict_,
190 1.8
191 )
192 ),
193 g3_
194 (
195 dimensioned<scalar>::getOrAddToDict
196 (
197 "g3",
198 this->coeffDict_,
199 0.8
200 )
201 ),
202 g3star_
203 (
204 dimensioned<scalar>::getOrAddToDict
205 (
206 "g3star",
207 this->coeffDict_,
208 1.3
209 )
210 ),
211 g4_
212 (
213 dimensioned<scalar>::getOrAddToDict
214 (
215 "g4",
216 this->coeffDict_,
217 1.25
218 )
219 ),
220 g5_
221 (
222 dimensioned<scalar>::getOrAddToDict
223 (
224 "g5",
225 this->coeffDict_,
226 0.2
227 )
228 ),
229 Cmu_
230 (
231 dimensioned<scalar>::getOrAddToDict
232 (
233 "Cmu",
234 this->coeffDict_,
235 0.21
236 )
237 ),
238 Ceps1_
239 (
240 dimensioned<scalar>::getOrAddToDict
241 (
242 "Ceps1",
243 this->coeffDict_,
244 1.44
245 )
246 ),
247 Ceps2_
248 (
249 dimensioned<scalar>::getOrAddToDict
250 (
251 "Ceps2",
252 this->coeffDict_,
253 1.83
254 )
255 ),
256 sigmaK_
257 (
258 dimensioned<scalar>::getOrAddToDict
259 (
260 "sigmaK",
261 this->coeffDict_,
262 1.0
263 )
264 ),
265 sigmaEps_
266 (
267 dimensioned<scalar>::getOrAddToDict
268 (
269 "sigmaEps",
270 this->coeffDict_,
271 1.15
272 )
273 ),
274 A1_
275 (
276 dimensioned<scalar>::getOrAddToDict
277 (
278 "A1",
279 this->coeffDict_,
280 0.065
281 )
282 ),
283 Ct_
284 (
285 dimensioned<scalar>::getOrAddToDict
286 (
287 "Ct",
288 this->coeffDict_,
289 6.0
290 )
291 ),
292 Cl_
293 (
294 dimensioned<scalar>::getOrAddToDict
295 (
296 "Cl",
297 this->coeffDict_,
298 0.133
299 )
300 ),
301 Ceta_
302 (
303 dimensioned<scalar>::getOrAddToDict
304 (
305 "Ceta",
306 this->coeffDict_,
307 80.0
308 )
309 ),
310 Cstability_
311 (
312 dimensioned<scalar>::getOrAddToDict
313 (
314 "Cstability",
315 this->coeffDict_,
317 10.0
318 )
319 ),
320
321 f_
322 (
324 (
325 IOobject::groupName("f", alphaRhoPhi.group()),
326 this->runTime_.timeName(),
327 this->mesh_,
328 IOobject::MUST_READ,
329 IOobject::AUTO_WRITE
330 ),
331 this->mesh_
332 ),
333 k_
334 (
336 (
337 IOobject::groupName("k", alphaRhoPhi.group()),
338 this->runTime_.timeName(),
339 this->mesh_,
340 IOobject::NO_READ,
341 IOobject::AUTO_WRITE
342 ),
343 0.5*tr(this->R_)
344 ),
345 epsilon_
346 (
348 (
349 IOobject::groupName("epsilon", alphaRhoPhi.group()),
350 this->runTime_.timeName(),
351 this->mesh_,
352 IOobject::MUST_READ,
353 IOobject::AUTO_WRITE
354 ),
355 this->mesh_
356 )
357{
358 this->boundNormalStress(this->R_);
359 bound(epsilon_, this->epsilonMin_);
360 bound(k_, this->kMin_);
361
362 if (type == typeName)
363 {
364 this->printCoeffs(type);
365 }
366}
367
368
369// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
370
371template<class BasicTurbulenceModel>
373{
375 {
376 simpleGradientDiffusion_.readIfPresent
377 (
378 "simpleGradientDiffusion",
379 this->coeffDict()
380 );
381 g1_.readIfPresent(this->coeffDict());
382 g1star_.readIfPresent(this->coeffDict());
383 g3_.readIfPresent(this->coeffDict());
384 g3star_.readIfPresent(this->coeffDict());
385 g4_.readIfPresent(this->coeffDict());
386 g5_.readIfPresent(this->coeffDict());
387 Cmu_.readIfPresent(this->coeffDict());
388 Ceps1_.readIfPresent(this->coeffDict());
389 Ceps2_.readIfPresent(this->coeffDict());
390 sigmaK_.readIfPresent(this->coeffDict());
391 sigmaEps_.readIfPresent(this->coeffDict());
392 A1_.readIfPresent(this->coeffDict());
393 Ct_.readIfPresent(this->coeffDict());
394 Cl_.readIfPresent(this->coeffDict());
395 Ceta_.readIfPresent(this->coeffDict());
396 Cstability_.readIfPresent(this->coeffDict());
397
398 return true;
399 }
400
401 return false;
402}
403
404
405template<class BasicTurbulenceModel>
407{
408 if (!this->turbulence_)
409 {
410 return;
411 }
412
413 // Construct local convenience references
414 const alphaField& alpha = this->alpha_;
415 const rhoField& rho = this->rho_;
416 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
417 const volVectorField& U = this->U_;
418 volSymmTensorField& R = this->R_;
420
422
423
424 // Calculate the velocity gradient tensor in Hessian form (delta_i u_j)
425 // Tranpose of the classical Jacobian form (delta_j u_i)
427 const volTensorField& gradU = tgradU.cref();
428
429 // Calculate the production tensor
430 tmp<volSymmTensorField> tP = -twoSymm(R & gradU);
431 const volSymmTensorField& P = tP.cref();
432
433 // Calculate turbulent kinetic energy production rate
434 const volScalarField G(this->GName(), 0.5*mag(tr(P)));
435
436
437 // Calculate elliptic blending function
438 // between near-wall and weakly-inhomogeneous regions
439 {
440 // (M:Eq. C.13)
441 tmp<volScalarField> tinvLsqr = scalar(1)/sqr(calcL());
442 const volScalarField& invLsqr = tinvLsqr.cref();
443
444 // (M:Eq. C.8)
446 (
447 fvm::Sp(invLsqr, f_)
448 - fvm::laplacian(f_)
449 ==
450 invLsqr
451 );
452
453 tinvLsqr.clear();
454
455 fEqn.ref().relax();
456 solve(fEqn);
457 }
458
459
460 // Calculate approximate wall-normal vector field (M:Eq. C.9)
461 tmp<volVectorField> tn = calcN();
462 const volVectorField& n = tn.cref();
463
464 // Calculate turbulent time scale (M:Eq. C.12)
465 tmp<volScalarField> ttau = calcTau();
466 const volScalarField& tau = ttau.cref();
467
468
469 // Calculate turbulent dissipation rate field
470 {
471 // Dissipation-production stimulator in the buffer layer (M:Eq. C.15)
472 tmp<volScalarField> tCeps1Prime = Ceps1Prime(G);
473 const volScalarField& Ceps1Prime = tCeps1Prime.cref();
474
475 // Update epsilon and G at the wall
476 epsilon_.boundaryFieldRef().updateCoeffs();
477
478 // (M:Eq. C.14)
480 (
481 fvm::ddt(alpha, rho, epsilon_)
482 + fvm::div(alphaRhoPhi, epsilon_)
483 - (
484 simpleGradientDiffusion_
485 ? fvm::laplacian(alpha*rho*D(sigmaEps_), epsilon_)
486 : fvm::laplacian(alpha*rho*D(tau, sigmaEps_), epsilon_)
487 )
488 + fvm::Sp(Cstability_/k_, epsilon_)
489 ==
490 alpha()*rho()*Ceps1Prime()*G()/tau()
491 - fvm::Sp(alpha*rho*Ceps2_/tau, epsilon_)
492 + Cstability_*epsilon_()/k_()
493 + fvOptions(alpha, rho, epsilon_)
494 );
495
496 tCeps1Prime.clear();
497
498 epsEqn.ref().relax();
499 fvOptions.constrain(epsEqn.ref());
500 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
501 solve(epsEqn);
502 fvOptions.correct(epsilon_);
503
504 bound(epsilon_, this->epsilonMin_);
505 }
506
507
508 // Calculate Reynolds-stress field
509 {
510 // Homogeneous form of the redistribution term (M:Eq. C.3)
512 {
513 // Reynolds stress anisotropy tensor (M:Eq. C.4)
514 const volSymmTensorField B(R/(scalar(2)*k_) - oneThirdI);
515
516 // Rate-of-strain tensor (M:Eq. C.5)
517 const volSymmTensorField S(symm(gradU));
518
519 // Rate-of-rotation tensor (M:Eq. C.6)
520 // Note the Hessian form of gradient
521 const volTensorField W(gradU.T() - gradU);
522
523 tPhiH =
524 k_
525 *(
526 (g3_ - g3star_*mag(B))*S
527 + g4_*dev(twoSymm(B & S))
528 + g5_*twoSymm(B & W.T())
529 );
530 }
531
532
533 // Near-wall form of the redistribution model (M:Eq. C.7)
535 {
537 const volSymmTensorField& nn = tnn.cref();
538
539 tn.clear();
540
541 tPhiW =
542 - scalar(5)*epsilon_/k_*
543 (
544 twoSymm(R & nn)
545 - 0.5*(R && nn)*(nn + I)
546 );
547 }
548
549
551 {
552 const volScalarField fCube(pow3(f_));
553
554 // Velocity-pressure gradient correlation (M:Eq. C.2)
555 const volSymmTensorField Phi
556 (
557 (scalar(1) - fCube)*tPhiW + fCube*tPhiH
558 );
559
560 // Near-wall part of the dissipation tensor (M:Eq. C.11)
561 const volScalarField epsilonW
562 (
563 (scalar(1) - fCube)*epsilon_/k_
564 );
565
566 // Homogeneous part of the dissipation tensor (M:Eq. C.11)
567 const volSphericalTensorField epsilonH
568 (
569 fCube*epsilon_*twoThirdsI
570 );
571
572 // Implicit part of the g1-term (M:Eq. C.3)
573 const volScalarField Phi1Implicit
574 (
575 fCube*(g1_*epsilon_ + g1star_*G)/(scalar(2)*k_)
576 );
577
578 // Explicit part of the g1-term (M:Eq. C.3)
579 const volSphericalTensorField Phi1Explicit
580 (
581 fCube*(g1_*epsilon_ + g1star_*G)*oneThirdI
582 );
583
584
585 // (M:Eq. C.1)
586 REqn =
587 (
589 + fvm::div(alphaRhoPhi, R)
590 - (
591 simpleGradientDiffusion_
592 ? fvm::laplacian(alpha*rho*D(sigmaK_), R)
593 : fvm::laplacian(alpha*rho*D(tau, sigmaK_), R)
594 )
595 + fvm::Sp(alpha*rho*Phi1Implicit, R)
596 + fvm::Sp(alpha*rho*epsilonW, R)
597 ==
598 alpha()*rho()*
599 (
600 P()
601 + Phi()
602 + Phi1Explicit()
603 - epsilonH()
604 )
605 + fvOptions(alpha, rho, R)
606 );
607
608 ttau.clear();
609 tP.clear();
610 }
611
612
613 REqn.ref().relax();
614 fvOptions.constrain(REqn.ref());
615 solve(REqn);
616 fvOptions.correct(R);
617
618 this->boundNormalStress(R);
619
620 #ifdef FULLDEBUG
621 this->checkRealizabilityConditions(R);
622 #endif
623
624 k_ == 0.5*tr(R);
625 bound(k_, this->kMin_);
626 }
627
628 correctNut();
629}
630
631
632// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
633
634} // End namespace RASModels
635} // End namespace Foam
636
637// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
#define R(A, B, C, D, E, F, K, M)
label n
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
surfaceScalarField & phi
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
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
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:77
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:74
Manceau (2015)'s elliptic-blending Reynolds-stress turbulence model for incompressible and compressib...
Definition: EBRSM.H:108
BasicTurbulenceModel::alphaField alphaField
Definition: EBRSM.H:187
BasicTurbulenceModel::rhoField rhoField
Definition: EBRSM.H:188
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
Definition: EBRSM.C:406
BasicTurbulenceModel::transportModel transportModel
Definition: EBRSM.H:189
virtual bool read()
Re-read model coefficients if they have changed.
Definition: EBRSM.C:372
Reynolds-stress turbulence model base class.
void boundNormalStress(volSymmTensorField &R) const
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
Generic dimensioned Type class.
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
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()
word timeName
Definition: getTimeIndex.H:3
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
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 Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:94
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
static const sphericalTensor twoThirdsI(2.0/3.0)
static const sphericalTensor oneThirdI(1.0/3.0)
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & nu
volScalarField & alpha
const dimensionedScalar & D
CEqn solve()