kOmegaSSTLM.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) 2016-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 "kOmegaSSTLM.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 (
45  const volScalarField& CDkOmega
46 ) const
47 {
48  const volScalarField Ry(this->y_*sqrt(this->k_)/this->nu());
49  const volScalarField F3(exp(-pow(Ry/120.0, 8)));
50 
51  return max(kOmegaSST<BasicTurbulenceModel>::F1(CDkOmega), F3);
52 }
53 
54 
55 template<class BasicTurbulenceModel>
57 (
59 ) const
60 {
61  return gammaIntEff_*kOmegaSST<BasicTurbulenceModel>::Pk(G);
62 }
63 
64 
65 template<class BasicTurbulenceModel>
67 (
68  const volScalarField& F1,
69  const volTensorField& gradU
70 ) const
71 {
72  return
73  clamp(gammaIntEff_, scalarMinMax(0.1, 1))
75 }
76 
77 
78 template<class BasicTurbulenceModel>
80 (
82  const volScalarField::Internal& Omega,
84 ) const
85 {
86  const volScalarField::Internal& omega = this->omega_();
87  const volScalarField::Internal& y = this->y_();
88 
89  dimensionedScalar deltaMin("deltaMin", dimLength, SMALL);
91  (
92  max(375*Omega*nu*ReThetat_()*y/sqr(Us), deltaMin)
93  );
94 
95  const volScalarField::Internal ReOmega(sqr(y)*omega/nu);
96  const volScalarField::Internal Fwake(exp(-sqr(ReOmega/1e5)));
97 
99  (
101  (
102  IOobject::groupName("Fthetat", this->alphaRhoPhi_.group()),
103  min
104  (
105  max
106  (
107  Fwake*exp(-pow4((y/delta))),
108  (1 - sqr((gammaInt_() - 1.0/ce2_)/(1 - 1.0/ce2_)))
109  ),
110  scalar(1)
111  )
112  )
113  );
114 }
115 
116 
117 template<class BasicTurbulenceModel>
120 {
122  (
124  (
125  IOobject
126  (
127  IOobject::groupName("ReThetac", this->alphaRhoPhi_.group()),
128  this->runTime_.timeName(),
129  this->mesh_
130  ),
131  this->mesh_,
132  dimless
133  )
134  );
135  volScalarField::Internal& ReThetac = tReThetac.ref();
136 
137  forAll(ReThetac, celli)
138  {
139  const scalar ReThetat = ReThetat_[celli];
140 
141  ReThetac[celli] =
142  ReThetat <= 1870
143  ?
144  ReThetat
145  - 396.035e-2
146  + 120.656e-4*ReThetat
147  - 868.230e-6*sqr(ReThetat)
148  + 696.506e-9*pow3(ReThetat)
149  - 174.105e-12*pow4(ReThetat)
150  :
151  ReThetat - 593.11 - 0.482*(ReThetat - 1870);
152  }
154  return tReThetac;
155 }
156 
157 
158 template<class BasicTurbulenceModel>
160 (
162 ) const
163 {
165  (
167  (
168  IOobject
169  (
170  IOobject::groupName("Flength", this->alphaRhoPhi_.group()),
171  this->runTime_.timeName(),
172  this->mesh_
173  ),
174  this->mesh_,
175  dimless
176  )
177  );
178  volScalarField::Internal& Flength = tFlength.ref();
179 
180  const volScalarField::Internal& omega = this->omega_();
181  const volScalarField::Internal& y = this->y_();
182 
183  forAll(ReThetat_, celli)
184  {
185  const scalar ReThetat = ReThetat_[celli];
186 
187  if (ReThetat < 400)
188  {
189  Flength[celli] =
190  398.189e-1
191  - 119.270e-4*ReThetat
192  - 132.567e-6*sqr(ReThetat);
193  }
194  else if (ReThetat < 596)
195  {
196  Flength[celli] =
197  263.404
198  - 123.939e-2*ReThetat
199  + 194.548e-5*sqr(ReThetat)
200  - 101.695e-8*pow3(ReThetat);
201  }
202  else if (ReThetat < 1200)
203  {
204  Flength[celli] = 0.5 - 3e-4*(ReThetat - 596);
205  }
206  else
207  {
208  Flength[celli] = 0.3188;
209  }
210 
211  const scalar Fsublayer =
212  exp(-sqr(sqr(y[celli])*omega[celli]/(200*nu[celli])));
213 
214  Flength[celli] = Flength[celli]*(1 - Fsublayer) + 40*Fsublayer;
215  }
217  return tFlength;
218 }
219 
220 
221 template<class BasicTurbulenceModel>
223 (
225  const volScalarField::Internal& dUsds,
227 ) const
228 {
230  (
232  (
233  IOobject
234  (
235  IOobject::groupName("ReThetat0", this->alphaRhoPhi_.group()),
236  this->runTime_.timeName(),
237  this->mesh_
238  ),
239  this->mesh_,
240  dimless
241  )
242  );
243  volScalarField::Internal& ReThetat0 = tReThetat0.ref();
244 
245  const volScalarField& k = this->k_;
246 
247  label maxIter = 0;
248 
249  forAll(ReThetat0, celli)
250  {
251  const scalar Tu
252  (
253  max(100*sqrt((2.0/3.0)*k[celli])/Us[celli], scalar(0.027))
254  );
255 
256  // Initialize lambda to zero.
257  // If lambda were cached between time-steps convergence would be faster
258  // starting from the previous time-step value.
259  scalar lambda = 0;
260 
261  scalar lambdaErr;
262  scalar thetat;
263  label iter = 0;
264 
265  do
266  {
267  // Previous iteration lambda for convergence test
268  const scalar lambda0 = lambda;
269 
270  if (Tu <= 1.3)
271  {
272  const scalar Flambda =
273  dUsds[celli] <= 0
274  ?
275  1
276  - (
277  - 12.986*lambda
278  - 123.66*sqr(lambda)
279  - 405.689*pow3(lambda)
280  )*exp(-pow(Tu/1.5, 1.5))
281  :
282  1
283  + 0.275*(1 - exp(-35*lambda))
284  *exp(-Tu/0.5);
285 
286  thetat =
287  (1173.51 - 589.428*Tu + 0.2196/sqr(Tu))
288  *Flambda*nu[celli]
289  /Us[celli];
290  }
291  else
292  {
293  const scalar Flambda =
294  dUsds[celli] <= 0
295  ?
296  1
297  - (
298  -12.986*lambda
299  -123.66*sqr(lambda)
300  -405.689*pow3(lambda)
301  )*exp(-pow(Tu/1.5, 1.5))
302  :
303  1
304  + 0.275*(1 - exp(-35*lambda))
305  *exp(-2*Tu);
306 
307  thetat =
308  331.50*pow((Tu - 0.5658), -0.671)
309  *Flambda*nu[celli]/Us[celli];
310  }
311 
312  lambda = sqr(thetat)/nu[celli]*dUsds[celli];
313  lambda = clamp(lambda, scalar(-0.1), scalar(0.1));
314 
315  lambdaErr = mag(lambda - lambda0);
316 
317  maxIter = max(maxIter, ++iter);
318 
319  } while (lambdaErr > lambdaErr_);
320 
321  ReThetat0[celli] = max(thetat*Us[celli]/nu[celli], scalar(20));
322  }
323 
324  if (maxIter > maxLambdaIter_)
325  {
327  << "Number of lambda iterations exceeds maxLambdaIter("
328  << maxLambdaIter_ << ')'<< endl;
329  }
331  return tReThetat0;
332 }
333 
334 
335 template<class BasicTurbulenceModel>
337 (
338  const volScalarField::Internal& Rev,
339  const volScalarField::Internal& ReThetac,
340  const volScalarField::Internal& RT
341 ) const
342 {
343  const volScalarField::Internal Fonset1(Rev/(2.193*ReThetac));
344 
345  const volScalarField::Internal Fonset2
346  (
347  min(max(Fonset1, pow4(Fonset1)), scalar(2))
348  );
349 
350  const volScalarField::Internal Fonset3(max(1 - pow3(RT/2.5), scalar(0)));
351 
353  (
355  (
356  IOobject::groupName("Fonset", this->alphaRhoPhi_.group()),
357  max(Fonset2 - Fonset3, scalar(0))
358  )
359  );
360 }
361 
362 
363 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
364 
365 template<class BasicTurbulenceModel>
367 (
368  const alphaField& alpha,
369  const rhoField& rho,
370  const volVectorField& U,
371  const surfaceScalarField& alphaRhoPhi,
372  const surfaceScalarField& phi,
373  const transportModel& transport,
374  const word& propertiesName,
375  const word& type
376 )
377 :
378  kOmegaSST<BasicTurbulenceModel>
379  (
380  alpha,
381  rho,
382  U,
383  alphaRhoPhi,
384  phi,
385  transport,
386  propertiesName,
387  typeName
388  ),
389 
390  ca1_
391  (
392  dimensionedScalar::getOrAddToDict
393  (
394  "ca1",
395  this->coeffDict_,
396  2
397  )
398  ),
399  ca2_
400  (
401  dimensionedScalar::getOrAddToDict
402  (
403  "ca2",
404  this->coeffDict_,
405  0.06
406  )
407  ),
408  ce1_
409  (
410  dimensionedScalar::getOrAddToDict
411  (
412  "ce1",
413  this->coeffDict_,
414  1
415  )
416  ),
417  ce2_
418  (
419  dimensionedScalar::getOrAddToDict
420  (
421  "ce2",
422  this->coeffDict_,
423  50
424  )
425  ),
426  cThetat_
427  (
428  dimensionedScalar::getOrAddToDict
429  (
430  "cThetat",
431  this->coeffDict_,
432  0.03
433  )
434  ),
435  sigmaThetat_
436  (
437  dimensionedScalar::getOrAddToDict
438  (
439  "sigmaThetat",
440  this->coeffDict_,
441  2
442  )
443  ),
444  lambdaErr_
445  (
446  this->coeffDict_.template getOrDefault<scalar>("lambdaErr", 1e-6)
447  ),
448  maxLambdaIter_
449  (
450  this->coeffDict_.template getOrDefault<label>("maxLambdaIter", 10)
451  ),
452  deltaU_("deltaU", dimVelocity, SMALL),
453 
454  ReThetat_
455  (
456  IOobject
457  (
458  IOobject::groupName("ReThetat", alphaRhoPhi.group()),
459  this->runTime_.timeName(),
460  this->mesh_,
461  IOobject::MUST_READ,
462  IOobject::AUTO_WRITE
463  ),
464  this->mesh_
465  ),
466 
467  gammaInt_
468  (
469  IOobject
470  (
471  IOobject::groupName("gammaInt", alphaRhoPhi.group()),
472  this->runTime_.timeName(),
473  this->mesh_,
474  IOobject::MUST_READ,
475  IOobject::AUTO_WRITE
476  ),
477  this->mesh_
478  ),
479 
480  gammaIntEff_
481  (
482  IOobject
483  (
484  IOobject::groupName("gammaIntEff", alphaRhoPhi.group()),
485  this->runTime_.timeName(),
486  this->mesh_
487  ),
488  this->mesh_,
490  )
491 {
492  if (type == typeName)
493  {
494  this->printCoeffs(type);
495  }
496 }
497 
498 
499 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
500 
501 template<class BasicTurbulenceModel>
503 {
505  {
506  ca1_.readIfPresent(this->coeffDict());
507  ca2_.readIfPresent(this->coeffDict());
508  ce1_.readIfPresent(this->coeffDict());
509  ce2_.readIfPresent(this->coeffDict());
510  sigmaThetat_.readIfPresent(this->coeffDict());
511  cThetat_.readIfPresent(this->coeffDict());
512  this->coeffDict().readIfPresent("lambdaErr", lambdaErr_);
513  this->coeffDict().readIfPresent("maxLambdaIter", maxLambdaIter_);
514 
515  return true;
516  }
517 
518  return false;
519 }
520 
521 
522 template<class BasicTurbulenceModel>
524 {
525  // Local references
526  const alphaField& alpha = this->alpha_;
527  const rhoField& rho = this->rho_;
528  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
529  const volVectorField& U = this->U_;
530  const volScalarField& k = this->k_;
531  const volScalarField& omega = this->omega_;
532  const tmp<volScalarField> tnu = this->nu();
533  const volScalarField::Internal& nu = tnu()();
534  const volScalarField::Internal& y = this->y_();
535  fv::options& fvOptions(fv::options::New(this->mesh_));
536 
537  // Fields derived from the velocity gradient
538  tmp<volTensorField> tgradU = fvc::grad(U);
539  const volScalarField::Internal Omega(sqrt(2*magSqr(skew(tgradU()()))));
540  const volScalarField::Internal S(sqrt(2*magSqr(symm(tgradU()()))));
541  const volScalarField::Internal Us(max(mag(U()), deltaU_));
542  const volScalarField::Internal dUsds((U() & (U() & tgradU()()))/sqr(Us));
543  tgradU.clear();
544 
545  const volScalarField::Internal Fthetat(this->Fthetat(Us, Omega, nu));
546 
547  {
548  const volScalarField::Internal t(500*nu/sqr(Us));
549  const volScalarField::Internal Pthetat
550  (
551  alpha()*rho()*(cThetat_/t)*(1 - Fthetat)
552  );
553 
554  // Transition onset momentum-thickness Reynolds number equation
555  tmp<fvScalarMatrix> ReThetatEqn
556  (
557  fvm::ddt(alpha, rho, ReThetat_)
558  + fvm::div(alphaRhoPhi, ReThetat_)
559  - fvm::laplacian(alpha*rho*DReThetatEff(), ReThetat_)
560  ==
561  Pthetat*ReThetat0(Us, dUsds, nu) - fvm::Sp(Pthetat, ReThetat_)
562  + fvOptions(alpha, rho, ReThetat_)
563  );
564 
565  ReThetatEqn.ref().relax();
566  fvOptions.constrain(ReThetatEqn.ref());
567  solve(ReThetatEqn);
568  fvOptions.correct(ReThetat_);
569  bound(ReThetat_, 0);
570  }
571 
572  const volScalarField::Internal ReThetac(this->ReThetac());
573  const volScalarField::Internal Rev(sqr(y)*S/nu);
574  const volScalarField::Internal RT(k()/(nu*omega()));
575 
576  {
577  const volScalarField::Internal Pgamma
578  (
579  alpha()*rho()
580  *ca1_*Flength(nu)*S*sqrt(gammaInt_()*Fonset(Rev, ReThetac, RT))
581  );
582 
583  const volScalarField::Internal Fturb(exp(-pow4(0.25*RT)));
584 
585  const volScalarField::Internal Egamma
586  (
587  alpha()*rho()*ca2_*Omega*Fturb*gammaInt_()
588  );
589 
590  // Intermittency equation
591  tmp<fvScalarMatrix> gammaIntEqn
592  (
593  fvm::ddt(alpha, rho, gammaInt_)
594  + fvm::div(alphaRhoPhi, gammaInt_)
595  - fvm::laplacian(alpha*rho*DgammaIntEff(), gammaInt_)
596  ==
597  Pgamma - fvm::Sp(ce1_*Pgamma, gammaInt_)
598  + Egamma - fvm::Sp(ce2_*Egamma, gammaInt_)
599  + fvOptions(alpha, rho, gammaInt_)
600  );
601 
602  gammaIntEqn.ref().relax();
603  fvOptions.constrain(gammaIntEqn.ref());
604  solve(gammaIntEqn);
605  fvOptions.correct(gammaInt_);
606  bound(gammaInt_, 0);
607  }
608 
609  const volScalarField::Internal Freattach(exp(-pow4(RT/20.0)));
610  const volScalarField::Internal gammaSep
611  (
612  min(2*max(Rev/(3.235*ReThetac) - 1, scalar(0))*Freattach, scalar(2))
613  *Fthetat
614  );
615 
616  gammaIntEff_ = max(gammaInt_(), gammaSep);
617 }
618 
619 
620 template<class BasicTurbulenceModel>
622 {
623  if (!this->turbulence_)
624  {
625  return;
626  }
627 
628  // Correct k and omega
630 
631  // Correct ReThetat and gammaInt
632  correctReThetatGammaInt();
633 }
634 
635 
636 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
637 
638 } // End namespace RASModels
639 } // End namespace Foam
640 
641 // ************************************************************************* //
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k which for standard RAS is betaStar*omega.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< volScalarField::Internal > ReThetat0(const volScalarField::Internal &Us, const volScalarField::Internal &dUsds, const volScalarField::Internal &nu) const
Return the transition onset momentum-thickness Reynolds number.
Definition: kOmegaSSTLM.C:216
#define F1(B, C, D)
Definition: SHA1.C:149
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
label k
Boltzmann constant.
const dimensionSet dimless
Dimensionless.
Finite-volume options.
Definition: fvOptions.H:51
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
fv::options & fvOptions
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
tmp< volScalarField::Internal > ReThetac() const
Empirical correlation for critical Reynolds number where the.
Definition: kOmegaSSTLM.C:112
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dimensionedScalar exp(const dimensionedScalar &ds)
#define F3(B, C, D)
Definition: SHA1.C:151
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Us
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:614
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
Definition: kOmegaSST.H:125
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTLM.C:495
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Bound the given scalar field if it has gone unbounded.
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Modified form of the k-omega SST epsilon/k.
Definition: kOmegaSSTLM.C:60
tmp< volScalarField::Internal > Fonset(const volScalarField::Internal &Rev, const volScalarField::Internal &ReThetac, const volScalarField::Internal &RT) const
Transition onset location control function.
Definition: kOmegaSSTLM.C:330
tmp< volScalarField::Internal > Flength(const volScalarField::Internal &nu) const
Empirical correlation that controls the length of the.
Definition: kOmegaSSTLM.C:153
tensor Ry(const scalar omega)
Rotational transformation tensor about the y-axis by omega radians.
Definition: transform.H:99
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Langtry-Menter 4-equation transitional SST model based on the k-omega-SST RAS model.
Definition: kOmegaSSTLM.H:103
Base-class for all transport models used by the incompressible turbulence models. ...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
virtual bool read()
Re-read model coefficients if they have changed.
void correctReThetatGammaInt()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTLM.C:516
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:289
tmp< volScalarField::Internal > Fthetat(const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
Freestream blending-function.
Definition: kOmegaSSTLM.C:73
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Modified form of the k-omega SST F1 function.
Definition: kOmegaSSTLM.C:37
volScalarField & nu
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Modified form of the k-omega SST k production rate.
Definition: kOmegaSSTLM.C:50