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-2017, 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 "bound.H"
31 #include "wallDist.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
39 
40 template<class BasicEddyViscosityModel>
42 (
43  const volScalarField& CDkOmega
44 ) const
45 {
46  tmp<volScalarField> CDkOmegaPlus = max
47  (
48  CDkOmega,
49  dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
50  );
51 
53  (
54  min
55  (
56  max
57  (
58  (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
59  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
60  ),
61  (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
62  ),
63  scalar(10)
64  );
65 
66  return tanh(pow4(arg1));
67 }
68 
69 
70 template<class BasicEddyViscosityModel>
72 {
74  (
75  max
76  (
77  (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
78  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
79  ),
80  scalar(100)
81  );
82 
83  return tanh(sqr(arg2));
84 }
85 
86 
87 template<class BasicEddyViscosityModel>
89 {
91  (
92  150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
93  scalar(10)
94  );
95 
96  return 1 - tanh(pow4(arg3));
97 }
98 
99 
100 template<class BasicEddyViscosityModel>
102 {
103  tmp<volScalarField> f23(F2());
104 
105  if (F3_)
106  {
107  f23.ref() *= F3();
108  }
109 
110  return f23;
111 }
112 
113 
114 template<class BasicEddyViscosityModel>
116 (
117  const volScalarField& S2
118 )
119 {
120  // Correct the turbulence viscosity
121  this->nut_ = a1_*k_/max(a1_*omega_, b1_*F23()*sqrt(S2));
122  this->nut_.correctBoundaryConditions();
123  fv::options::New(this->mesh_).correct(this->nut_);
124 }
125 
126 
127 template<class BasicEddyViscosityModel>
129 {
130  correctNut(2*magSqr(symm(fvc::grad(this->U_))));
131 }
132 
133 
134 template<class BasicEddyViscosityModel>
136 (
138 ) const
139 {
140  return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
141 }
142 
143 
144 template<class BasicEddyViscosityModel>
147 (
148  const volScalarField& F1,
149  const volTensorField& gradU
150 ) const
151 {
152  return betaStar_*omega_();
153 }
154 
155 
156 template<class BasicEddyViscosityModel>
158 (
159  const volScalarField::Internal& GbyNu0,
161  const volScalarField::Internal& S2
162 ) const
163 {
164  return min
165  (
166  GbyNu0,
167  (c1_/a1_)*betaStar_*omega_()
168  *max(a1_*omega_(), b1_*F2*sqrt(S2))
169  );
170 }
171 
172 
173 template<class BasicEddyViscosityModel>
175 {
176  return tmp<fvScalarMatrix>
177  (
178  new fvScalarMatrix
179  (
180  k_,
181  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
182  )
183  );
184 }
185 
186 
187 template<class BasicEddyViscosityModel>
189 {
190  return tmp<fvScalarMatrix>
191  (
192  new fvScalarMatrix
193  (
194  omega_,
195  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
196  )
197  );
198 }
199 
200 
201 template<class BasicEddyViscosityModel>
203 (
204  const volScalarField::Internal& S2,
207 ) const
208 {
209  return tmp<fvScalarMatrix>
210  (
211  new fvScalarMatrix
212  (
213  omega_,
214  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
215  )
216  );
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
221 
222 template<class BasicEddyViscosityModel>
224 (
225  const word& type,
226  const alphaField& alpha,
227  const rhoField& rho,
228  const volVectorField& U,
229  const surfaceScalarField& alphaRhoPhi,
230  const surfaceScalarField& phi,
231  const transportModel& transport,
232  const word& propertiesName
233 )
234 :
235  BasicEddyViscosityModel
236  (
237  type,
238  alpha,
239  rho,
240  U,
241  alphaRhoPhi,
242  phi,
243  transport,
244  propertiesName
245  ),
246 
247  alphaK1_
248  (
250  (
251  "alphaK1",
252  this->coeffDict_,
253  0.85
254  )
255  ),
256  alphaK2_
257  (
259  (
260  "alphaK2",
261  this->coeffDict_,
262  1.0
263  )
264  ),
265  alphaOmega1_
266  (
268  (
269  "alphaOmega1",
270  this->coeffDict_,
271  0.5
272  )
273  ),
274  alphaOmega2_
275  (
277  (
278  "alphaOmega2",
279  this->coeffDict_,
280  0.856
281  )
282  ),
283  gamma1_
284  (
286  (
287  "gamma1",
288  this->coeffDict_,
289  5.0/9.0
290  )
291  ),
292  gamma2_
293  (
295  (
296  "gamma2",
297  this->coeffDict_,
298  0.44
299  )
300  ),
301  beta1_
302  (
304  (
305  "beta1",
306  this->coeffDict_,
307  0.075
308  )
309  ),
310  beta2_
311  (
313  (
314  "beta2",
315  this->coeffDict_,
316  0.0828
317  )
318  ),
319  betaStar_
320  (
322  (
323  "betaStar",
324  this->coeffDict_,
325  0.09
326  )
327  ),
328  a1_
329  (
331  (
332  "a1",
333  this->coeffDict_,
334  0.31
335  )
336  ),
337  b1_
338  (
340  (
341  "b1",
342  this->coeffDict_,
343  1.0
344  )
345  ),
346  c1_
347  (
349  (
350  "c1",
351  this->coeffDict_,
352  10.0
353  )
354  ),
355  F3_
356  (
358  (
359  "F3",
360  this->coeffDict_,
361  false
362  )
363  ),
364 
365  y_(wallDist::New(this->mesh_).y()),
366 
367  k_
368  (
369  IOobject
370  (
371  IOobject::groupName("k", alphaRhoPhi.group()),
372  this->runTime_.timeName(),
373  this->mesh_,
376  ),
377  this->mesh_
378  ),
379  omega_
380  (
381  IOobject
382  (
383  IOobject::groupName("omega", alphaRhoPhi.group()),
384  this->runTime_.timeName(),
385  this->mesh_,
388  ),
389  this->mesh_
390  ),
391  decayControl_
392  (
394  (
395  "decayControl",
396  this->coeffDict_,
397  false
398  )
399  ),
400  kInf_
401  (
403  (
404  "kInf",
405  this->coeffDict_,
406  k_.dimensions(),
407  0
408  )
409  ),
410  omegaInf_
411  (
413  (
414  "omegaInf",
415  this->coeffDict_,
416  omega_.dimensions(),
417  0
418  )
419  )
420 {
421  bound(k_, this->kMin_);
422  bound(omega_, this->omegaMin_);
423 
424  setDecayControl(this->coeffDict_);
425 }
426 
427 
428 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
429 
430 template<class BasicEddyViscosityModel>
432 (
433  const dictionary& dict
434 )
435 {
436  decayControl_.readIfPresent("decayControl", dict);
437 
438  if (decayControl_)
439  {
440  kInf_.read(dict);
441  omegaInf_.read(dict);
442 
443  Info<< " Employing decay control with kInf:" << kInf_
444  << " and omegaInf:" << omegaInf_ << endl;
445  }
446  else
447  {
448  kInf_.value() = 0;
449  omegaInf_.value() = 0;
450  }
451 }
452 
453 
454 template<class BasicEddyViscosityModel>
456 {
458  {
459  alphaK1_.readIfPresent(this->coeffDict());
460  alphaK2_.readIfPresent(this->coeffDict());
461  alphaOmega1_.readIfPresent(this->coeffDict());
462  alphaOmega2_.readIfPresent(this->coeffDict());
463  gamma1_.readIfPresent(this->coeffDict());
464  gamma2_.readIfPresent(this->coeffDict());
465  beta1_.readIfPresent(this->coeffDict());
466  beta2_.readIfPresent(this->coeffDict());
467  betaStar_.readIfPresent(this->coeffDict());
468  a1_.readIfPresent(this->coeffDict());
469  b1_.readIfPresent(this->coeffDict());
470  c1_.readIfPresent(this->coeffDict());
471  F3_.readIfPresent("F3", this->coeffDict());
472 
473  setDecayControl(this->coeffDict());
474 
475  return true;
476  }
477 
478  return false;
479 }
480 
481 
482 template<class BasicEddyViscosityModel>
484 {
485  if (!this->turbulence_)
486  {
487  return;
488  }
489 
490  // Local references
491  const alphaField& alpha = this->alpha_;
492  const rhoField& rho = this->rho_;
493  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
494  const volVectorField& U = this->U_;
495  volScalarField& nut = this->nut_;
496  fv::options& fvOptions(fv::options::New(this->mesh_));
497 
499 
501 
502  tmp<volTensorField> tgradU = fvc::grad(U);
503  volScalarField S2(2*magSqr(symm(tgradU())));
504  volScalarField::Internal GbyNu0((tgradU() && dev(twoSymm(tgradU()))));
505  volScalarField::Internal G(this->GName(), nut*GbyNu0);
506 
507  // Update omega and G at the wall
508  omega_.boundaryFieldRef().updateCoeffs();
509 
510  volScalarField CDkOmega
511  (
512  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
513  );
514 
515  volScalarField F1(this->F1(CDkOmega));
516  volScalarField F23(this->F23());
517 
518  {
521 
522  // Turbulent frequency equation
523  tmp<fvScalarMatrix> omegaEqn
524  (
525  fvm::ddt(alpha, rho, omega_)
526  + fvm::div(alphaRhoPhi, omega_)
527  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
528  ==
529  alpha()*rho()*gamma*GbyNu(GbyNu0, F23(), S2())
530  - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
531  - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
532  - fvm::SuSp
533  (
534  alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
535  omega_
536  )
537  + alpha()*rho()*beta*sqr(omegaInf_)
538  + Qsas(S2(), gamma, beta)
539  + omegaSource()
540  + fvOptions(alpha, rho, omega_)
541  );
542 
543  omegaEqn.ref().relax();
544  fvOptions.constrain(omegaEqn.ref());
545  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
546  solve(omegaEqn);
547  fvOptions.correct(omega_);
548  bound(omega_, this->omegaMin_);
549  }
550 
551  // Turbulent kinetic energy equation
553  (
554  fvm::ddt(alpha, rho, k_)
555  + fvm::div(alphaRhoPhi, k_)
556  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
557  ==
558  alpha()*rho()*Pk(G)
559  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
560  - fvm::Sp(alpha()*rho()*epsilonByk(F1, tgradU()), k_)
561  + alpha()*rho()*betaStar_*omegaInf_*kInf_
562  + kSource()
563  + fvOptions(alpha, rho, k_)
564  );
565 
566  tgradU.clear();
567 
568  kEqn.ref().relax();
569  fvOptions.constrain(kEqn.ref());
570  solve(kEqn);
571  fvOptions.correct(k_);
572  bound(k_, this->kMin_);
573 
574  correctNut(S2);
575 }
576 
577 
578 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
579 
580 } // End namespace Foam
581 
582 // ************************************************************************* //
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:203
Foam::kOmegaSSTBase::F2
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:71
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
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:742
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Switch::lookupOrAddToDict
static Switch lookupOrAddToDict(const word &name, dictionary &dict, const Switch deflt=switchType::FALSE)
Definition: Switch.H:255
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:321
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:325
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< volScalarField >
Foam::kOmegaSSTBase::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmegaSSTBase.C:174
Foam::kOmegaSSTBase::F1
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:42
Foam::fv::optionList::constrain
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
Definition: fvOptionListTemplates.C:288
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
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:188
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:147
Foam::kOmegaSSTBase::setDecayControl
void setDecayControl(const dictionary &dict)
Definition: kOmegaSSTBase.C:432
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:337
Foam::kOmegaSSTBase::Pk
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
Definition: kOmegaSSTBase.C:136
rho
rho
Definition: readInitialConditions.H:96
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:128
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
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:258
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:158
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 (uses stdout - output is on the master only)
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:574
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:455
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:43
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
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
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:121
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:97
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:96
Foam::kOmegaSSTBase::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:483
Foam::kOmegaSSTBase::F23
virtual tmp< volScalarField > F23() const
Definition: kOmegaSSTBase.C:101
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:98
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::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
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:76
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
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:88
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
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106