kOmegaSSTBase.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-2015 OpenFOAM Foundation
9  Copyright (C) 2016-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 "kOmegaSSTBase.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 #include "wallDist.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
40 
41 template<class BasicEddyViscosityModel>
43 (
44  const volScalarField& CDkOmega
45 ) const
46 {
47  tmp<volScalarField> CDkOmegaPlus = max
48  (
49  CDkOmega,
50  dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
51  );
52 
54  (
55  min
56  (
57  max
58  (
59  (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
60  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
61  ),
62  (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
63  ),
64  scalar(10)
65  );
66 
67  return tanh(pow4(arg1));
68 }
69 
70 
71 template<class BasicEddyViscosityModel>
73 {
75  (
76  max
77  (
78  (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
79  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
80  ),
81  scalar(100)
82  );
83 
84  return tanh(sqr(arg2));
85 }
86 
87 
88 template<class BasicEddyViscosityModel>
90 {
92  (
93  150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
94  scalar(10)
95  );
96 
97  return 1 - tanh(pow4(arg3));
98 }
99 
100 
101 template<class BasicEddyViscosityModel>
103 {
104  tmp<volScalarField> f23(F2());
105 
106  if (F3_)
107  {
108  f23.ref() *= F3();
109  }
110 
111  return f23;
112 }
113 
114 
115 template<class BasicEddyViscosityModel>
117 (
118  const volScalarField& S2
119 )
120 {
121  // Correct the turbulence viscosity
122  this->nut_ = a1_*k_/max(a1_*omega_, b1_*F23()*sqrt(S2));
123  this->nut_.correctBoundaryConditions();
124  fv::options::New(this->mesh_).correct(this->nut_);
125 }
126 
127 
128 template<class BasicEddyViscosityModel>
130 {
131  correctNut(2*magSqr(symm(fvc::grad(this->U_))));
132 }
133 
134 
135 template<class BasicEddyViscosityModel>
137 (
139 ) const
140 {
141  return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
142 }
143 
144 
145 template<class BasicEddyViscosityModel>
148 (
149  const volScalarField& F1,
150  const volTensorField& gradU
151 ) const
152 {
153  return betaStar_*omega_();
154 }
155 
156 
157 template<class BasicEddyViscosityModel>
159 (
160  const volScalarField::Internal& GbyNu0,
162  const volScalarField::Internal& S2
163 ) const
164 {
165  return min
166  (
167  GbyNu0,
168  (c1_/a1_)*betaStar_*omega_()
169  *max(a1_*omega_(), b1_*F2*sqrt(S2))
170  );
171 }
172 
173 
174 template<class BasicEddyViscosityModel>
176 {
177  return tmp<fvScalarMatrix>
178  (
179  new fvScalarMatrix
180  (
181  k_,
182  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
183  )
184  );
185 }
186 
187 
188 template<class BasicEddyViscosityModel>
190 {
191  return tmp<fvScalarMatrix>
192  (
193  new fvScalarMatrix
194  (
195  omega_,
196  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
197  )
198  );
199 }
200 
201 
202 template<class BasicEddyViscosityModel>
204 (
205  const volScalarField::Internal& S2,
208 ) const
209 {
210  return tmp<fvScalarMatrix>
211  (
212  new fvScalarMatrix
213  (
214  omega_,
215  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
216  )
217  );
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
222 
223 template<class BasicEddyViscosityModel>
225 (
226  const word& type,
227  const alphaField& alpha,
228  const rhoField& rho,
229  const volVectorField& U,
230  const surfaceScalarField& alphaRhoPhi,
231  const surfaceScalarField& phi,
232  const transportModel& transport,
233  const word& propertiesName
234 )
235 :
236  BasicEddyViscosityModel
237  (
238  type,
239  alpha,
240  rho,
241  U,
242  alphaRhoPhi,
243  phi,
244  transport,
245  propertiesName
246  ),
247 
248  alphaK1_
249  (
251  (
252  "alphaK1",
253  this->coeffDict_,
254  0.85
255  )
256  ),
257  alphaK2_
258  (
260  (
261  "alphaK2",
262  this->coeffDict_,
263  1.0
264  )
265  ),
266  alphaOmega1_
267  (
269  (
270  "alphaOmega1",
271  this->coeffDict_,
272  0.5
273  )
274  ),
275  alphaOmega2_
276  (
278  (
279  "alphaOmega2",
280  this->coeffDict_,
281  0.856
282  )
283  ),
284  gamma1_
285  (
287  (
288  "gamma1",
289  this->coeffDict_,
290  5.0/9.0
291  )
292  ),
293  gamma2_
294  (
296  (
297  "gamma2",
298  this->coeffDict_,
299  0.44
300  )
301  ),
302  beta1_
303  (
305  (
306  "beta1",
307  this->coeffDict_,
308  0.075
309  )
310  ),
311  beta2_
312  (
314  (
315  "beta2",
316  this->coeffDict_,
317  0.0828
318  )
319  ),
320  betaStar_
321  (
323  (
324  "betaStar",
325  this->coeffDict_,
326  0.09
327  )
328  ),
329  a1_
330  (
332  (
333  "a1",
334  this->coeffDict_,
335  0.31
336  )
337  ),
338  b1_
339  (
341  (
342  "b1",
343  this->coeffDict_,
344  1.0
345  )
346  ),
347  c1_
348  (
350  (
351  "c1",
352  this->coeffDict_,
353  10.0
354  )
355  ),
356  F3_
357  (
359  (
360  "F3",
361  this->coeffDict_,
362  false
363  )
364  ),
365 
366  y_(wallDist::New(this->mesh_).y()),
367 
368  k_
369  (
370  IOobject
371  (
372  IOobject::groupName("k", alphaRhoPhi.group()),
373  this->runTime_.timeName(),
374  this->mesh_,
377  ),
378  this->mesh_
379  ),
380  omega_
381  (
382  IOobject
383  (
384  IOobject::groupName("omega", alphaRhoPhi.group()),
385  this->runTime_.timeName(),
386  this->mesh_,
389  ),
390  this->mesh_
391  ),
392  decayControl_
393  (
395  (
396  "decayControl",
397  this->coeffDict_,
398  false
399  )
400  ),
401  kInf_
402  (
404  (
405  "kInf",
406  this->coeffDict_,
407  k_.dimensions(),
408  0
409  )
410  ),
411  omegaInf_
412  (
414  (
415  "omegaInf",
416  this->coeffDict_,
417  omega_.dimensions(),
418  0
419  )
420  )
421 {
422  bound(k_, this->kMin_);
423  bound(omega_, this->omegaMin_);
424 
425  setDecayControl(this->coeffDict_);
426 }
427 
428 
429 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
430 
431 template<class BasicEddyViscosityModel>
433 (
434  const dictionary& dict
435 )
436 {
437  decayControl_.readIfPresent("decayControl", dict);
438 
439  if (decayControl_)
440  {
441  kInf_.read(dict);
442  omegaInf_.read(dict);
443 
444  Info<< " Employing decay control with kInf:" << kInf_
445  << " and omegaInf:" << omegaInf_ << endl;
446  }
447  else
448  {
449  kInf_.value() = 0;
450  omegaInf_.value() = 0;
451  }
452 }
453 
454 
455 template<class BasicEddyViscosityModel>
457 {
459  {
460  alphaK1_.readIfPresent(this->coeffDict());
461  alphaK2_.readIfPresent(this->coeffDict());
462  alphaOmega1_.readIfPresent(this->coeffDict());
463  alphaOmega2_.readIfPresent(this->coeffDict());
464  gamma1_.readIfPresent(this->coeffDict());
465  gamma2_.readIfPresent(this->coeffDict());
466  beta1_.readIfPresent(this->coeffDict());
467  beta2_.readIfPresent(this->coeffDict());
468  betaStar_.readIfPresent(this->coeffDict());
469  a1_.readIfPresent(this->coeffDict());
470  b1_.readIfPresent(this->coeffDict());
471  c1_.readIfPresent(this->coeffDict());
472  F3_.readIfPresent("F3", this->coeffDict());
473 
474  setDecayControl(this->coeffDict());
475 
476  return true;
477  }
478 
479  return false;
480 }
481 
482 
483 template<class BasicEddyViscosityModel>
485 {
486  if (!this->turbulence_)
487  {
488  return;
489  }
490 
491  // Local references
492  const alphaField& alpha = this->alpha_;
493  const rhoField& rho = this->rho_;
494  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
495  const volVectorField& U = this->U_;
496  volScalarField& nut = this->nut_;
497  fv::options& fvOptions(fv::options::New(this->mesh_));
498 
500 
502 
503  tmp<volTensorField> tgradU = fvc::grad(U);
504  volScalarField S2(2*magSqr(symm(tgradU())));
506  (
507  this->type() + ":GbyNu",
508  (tgradU() && dev(twoSymm(tgradU())))
509  );
510  volScalarField::Internal G(this->GName(), nut*GbyNu0);
511 
512  // Update omega and G at the wall
513  omega_.boundaryFieldRef().updateCoeffs();
514 
515  volScalarField CDkOmega
516  (
517  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
518  );
519 
520  volScalarField F1(this->F1(CDkOmega));
521  volScalarField F23(this->F23());
522 
523  {
526 
527  GbyNu0 = GbyNu(GbyNu0, F23(), S2());
528 
529  // Turbulent frequency equation
530  tmp<fvScalarMatrix> omegaEqn
531  (
532  fvm::ddt(alpha, rho, omega_)
533  + fvm::div(alphaRhoPhi, omega_)
534  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
535  ==
536  alpha()*rho()*gamma*GbyNu0
537  - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
538  - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
539  - fvm::SuSp
540  (
541  alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
542  omega_
543  )
544  + alpha()*rho()*beta*sqr(omegaInf_)
545  + Qsas(S2(), gamma, beta)
546  + omegaSource()
547  + fvOptions(alpha, rho, omega_)
548  );
549 
550  omegaEqn.ref().relax();
551  fvOptions.constrain(omegaEqn.ref());
552  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
553  solve(omegaEqn);
554  fvOptions.correct(omega_);
555  bound(omega_, this->omegaMin_);
556  }
557 
558  // Turbulent kinetic energy equation
560  (
561  fvm::ddt(alpha, rho, k_)
562  + fvm::div(alphaRhoPhi, k_)
563  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
564  ==
565  alpha()*rho()*Pk(G)
566  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
567  - fvm::Sp(alpha()*rho()*epsilonByk(F1, tgradU()), k_)
568  + alpha()*rho()*betaStar_*omegaInf_*kInf_
569  + kSource()
570  + fvOptions(alpha, rho, k_)
571  );
572 
573  tgradU.clear();
574 
575  kEqn.ref().relax();
576  fvOptions.constrain(kEqn.ref());
577  solve(kEqn);
578  fvOptions.correct(k_);
579  bound(k_, this->kMin_);
580 
581  correctNut(S2);
582 }
583 
584 
585 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
586 
587 } // End namespace Foam
588 
589 // ************************************************************************* //
wallDist.H
Foam::kOmegaSSTBase::Qsas
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
Definition: kOmegaSSTBase.C:204
Foam::kOmegaSSTBase::F2
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:72
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
F1
#define F1(B, C, D)
Definition: SHA1.C:150
Foam::fvMatrix::boundaryManipulate
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1348
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::kOmegaSSTBase::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmegaSSTBase.C:175
Foam::kOmegaSSTBase::F1
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:43
Foam::fv::optionList::constrain
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Definition: fvOptionListTemplates.C:314
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::kOmegaSSTBase::omegaSource
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmegaSSTBase.C:189
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::kOmegaSSTBase::epsilonByk
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k which for standard RAS is betaStar*omega.
Definition: kOmegaSSTBase.C:148
Foam::kOmegaSSTBase::setDecayControl
void setDecayControl(const dictionary &dict)
Definition: kOmegaSSTBase.C:433
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::kOmegaSSTBase::Pk
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
Definition: kOmegaSSTBase.C:137
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::kOmegaSSTBase::correctNut
virtual void correctNut()
Definition: kOmegaSSTBase.C:129
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::kOmegaSSTBase::GbyNu
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
Definition: kOmegaSSTBase.C:159
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
kOmegaSSTBase.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
F3
#define F3(B, C, D)
Definition: SHA1.C:152
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::kOmegaSSTBase::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTBase.C:456
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
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:123
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:98
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:97
Foam::kOmegaSSTBase::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:484
Foam::kOmegaSSTBase::F23
virtual tmp< volScalarField > F23() const
Definition: kOmegaSSTBase.C:102
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:99
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
gamma
const scalar gamma
Definition: EEqn.H:9
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
F2
#define F2(B, C, D)
Definition: SHA1.C:151
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::kOmegaSSTBase::F3
virtual tmp< volScalarField > F3() const
Definition: kOmegaSSTBase.C:89
Foam::kOmegaSSTBase
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
Definition: kOmegaSSTBase.H:128
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::Switch::getOrAddToDict
static Switch getOrAddToDict(const word &key, dictionary &dict, const Switch deflt=switchType::FALSE)
Definition: Switch.C:164
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106