41 template<
class BasicMomentumTransportModel>
66 hasDispersedPhaseNames_(this->coeffDict().
found(
"dispersedPhases")),
70 this->coeffDict().lookupOrDefault(
"dispersedPhases",
hashedWordList())
73 Cmub_(
"Cmub", this->coeffDict(), 0.6)
79 template<
class BasicMomentumTransportModel>
84 Cmub_.readIfPresent(this->coeffDict());
95 template<
class BasicMomentumTransportModel>
103 if (hasDispersedPhaseNames_)
105 dispersedPhases.
resize(dispersedPhaseNames_.size());
107 forAll(dispersedPhaseNames_, dispersedPhasei)
112 &fluid.
phases()[dispersedPhaseNames_[dispersedPhasei]]
120 label dispersedPhasei = 0;
124 const phaseModel& otherPhase = fluid.
movingPhases()[movingPhasei];
126 if (&otherPhase != &phase_)
137 return dispersedPhases;
141 template<
class BasicMomentumTransportModel>
150 pow(this->betaStar_, 0.25)*this->
y()*
sqrt(this->k_)/this->nu()
154 this->a1_*this->k_/
max(this->a1_*this->omega_, this->b1_*
F2*
sqrt(S2));
161 *Cmub_*phaseIter().d()*phaseIter()
162 *(
mag(this->U_ - phaseIter().
U()));
165 this->nut_.correctBoundaryConditions();
170 template<
class BasicMomentumTransportModel>
#define forAll(list, i)
Loop across all elements in list.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Generic GeometricField class.
Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble...
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
kOmegaSSTSato(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual bool read()
Read model coefficients if they have changed.
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information,...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
bool set(const label) const
Is element set.
label size() const
Return the number of elements in the UPtrList.
void resize(const label)
Reset size of UPtrList. This can only be used to set the size.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
A wordList with hashed indices for faster lookup by name.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual void correctNut()
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Class to represent a system of phases.
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
const phaseModelList & phases() const
Return the phase models.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
To & refCast(From &r)
Reference type cast template function.
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.