39 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
43 omega_ =
max(omega_, k_/(this->nutMaxCoeff_*this->nu()));
47 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
66 (scalar(1)/betaStar_)*
sqrt(k_)/(omega_*this->
y()),
67 scalar(500)*this->nu()/(
sqr(this->
y())*omega_)
69 (4*alphaOmega2_)*k_/(CDkOmegaPlus*
sqr(this->
y()))
77 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
85 (scalar(2)/betaStar_)*
sqrt(k_)/(omega_*this->
y()),
86 scalar(500)*this->nu()/(
sqr(this->
y())*omega_)
94 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
100 150*this->nu()/(omega_*
sqr(this->
y())),
107 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
122 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
129 this->nut_ = a1_*k_/
max(a1_*omega_, b1_*
F2*
sqrt(S2));
130 this->nut_.correctBoundaryConditions();
135 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
143 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
150 return min(
G, (c1_*betaStar_)*this->k_()*this->omega_());
154 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
162 return betaStar_*omega_();
166 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
181 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
197 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
219 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
231 MomentumTransportModel
364 this->groupName(
"k"),
365 this->runTime_.
name(),
376 this->groupName(
"omega"),
377 this->runTime_.
name(),
385 bound(k_, this->kMin_);
392 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
397 alphaK1_.readIfPresent(this->coeffDict());
398 alphaK2_.readIfPresent(this->coeffDict());
399 alphaOmega1_.readIfPresent(this->coeffDict());
400 alphaOmega2_.readIfPresent(this->coeffDict());
401 gamma1_.readIfPresent(this->coeffDict());
402 gamma2_.readIfPresent(this->coeffDict());
403 beta1_.readIfPresent(this->coeffDict());
404 beta2_.readIfPresent(this->coeffDict());
405 betaStar_.readIfPresent(this->coeffDict());
406 a1_.readIfPresent(this->coeffDict());
407 b1_.readIfPresent(this->coeffDict());
408 c1_.readIfPresent(this->coeffDict());
409 F3_.readIfPresent(
"F3", this->coeffDict());
420 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
423 if (!this->turbulence_)
454 omega_.boundaryFieldRef().updateCoeffs();
479 (c1_/a1_)*betaStar_*omega_()
480 *
max(a1_*omega_(), b1_*F23()*
sqrt(S2()))
486 alpha()*
rho()*(
F1() - scalar(1))*CDkOmega()/omega_(),
489 + Qsas(S2(), gamma, beta)
494 omegaEqn.
ref().relax();
496 omegaEqn.
ref().boundaryManipulate(omega_.boundaryFieldRef());
520 bound(k_, this->kMin_);
Bound the given scalar field where it is below the specified minimum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Generic dimensioned Type class.
Finite volume constraints.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
BasicMomentumTransportModel::alphaField alphaField
void boundOmega()
Bound omega.
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
virtual tmp< volScalarField > F23() const
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
kOmegaSST(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
virtual bool read()
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
virtual tmp< volScalarField > F3() const
A class for managing temporary objects.
void clear() const
If object pointer points to valid object:
T & ref() const
Return non-const reference or generate a fatal error.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
dimensionedScalar tanh(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.