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
242 alphaK1_(
"alphaK1", this->coeffDict(), 0.85),
243 alphaK2_(
"alphaK2", this->coeffDict(), 1.0),
244 alphaOmega1_(
"alphaOmega1", this->coeffDict(), 0.5),
245 alphaOmega2_(
"alphaOmega2", this->coeffDict(), 0.856),
246 gamma1_(
"gamma1", this->coeffDict(), 5.0/9.0),
247 gamma2_(
"gamma2", this->coeffDict(), 0.44),
248 beta1_(
"beta1", this->coeffDict(), 0.075),
249 beta2_(
"beta2", this->coeffDict(), 0.0828),
250 betaStar_(
"betaStar", this->coeffDict(), 0.09),
251 a1_(
"a1", this->coeffDict(), 0.31),
252 b1_(
"b1", this->coeffDict(), 1.0),
253 c1_(
"c1", this->coeffDict(), 10.0),
254 F3_(this->coeffDict().template lookupOrDefault<
Switch>(
"F3", false)),
260 this->groupName(
"k"),
261 this->runTime_.
name(),
272 this->groupName(
"omega"),
273 this->runTime_.
name(),
281 bound(k_, this->kMin_);
288 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
293 alphaK1_.readIfPresent(this->coeffDict());
294 alphaK2_.readIfPresent(this->coeffDict());
295 alphaOmega1_.readIfPresent(this->coeffDict());
296 alphaOmega2_.readIfPresent(this->coeffDict());
297 gamma1_.readIfPresent(this->coeffDict());
298 gamma2_.readIfPresent(this->coeffDict());
299 beta1_.readIfPresent(this->coeffDict());
300 beta2_.readIfPresent(this->coeffDict());
301 betaStar_.readIfPresent(this->coeffDict());
302 a1_.readIfPresent(this->coeffDict());
303 b1_.readIfPresent(this->coeffDict());
304 c1_.readIfPresent(this->coeffDict());
305 F3_.readIfPresent(
"F3", this->coeffDict());
316 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
319 if (!this->turbulence_)
350 omega_.boundaryFieldRef().updateCoeffs();
375 (c1_/a1_)*betaStar_*omega_()
376 *
max(a1_*omega_(), b1_*F23()*
sqrt(S2()))
382 alpha()*
rho()*(
F1() - scalar(1))*CDkOmega()/omega_(),
385 + Qsas(S2(), gamma, beta)
390 omegaEqn.
ref().relax();
392 omegaEqn.
ref().boundaryManipulate(omega_.boundaryFieldRef());
416 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.
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 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)
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimless
dimensionedScalar tanh(const dimensionedScalar &ds)
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimVolume
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.
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.