Go to the documentation of this file.
57 "transportProperties",
60 IOobject::MUST_READ_IF_MODIFIED,
73 IOobject::groupName(
"alpha", phase1Name_),
86 IOobject::groupName(
"alpha", phase2Name_),
99 IOobject::groupName(
"alpha", phase3Name_),
121 calculatedFvPatchScalarField::typeName
129 subDict(phase1Name_),
139 subDict(phase2Name_),
149 subDict(phase3Name_),
155 rho1_(
"rho",
dimDensity, nuModel1_->viscosityProperties()),
156 rho2_(
"rho",
dimDensity, nuModel2_->viscosityProperties()),
157 rho3_(
"rho",
dimDensity, nuModel3_->viscosityProperties())
169 return tmp<volScalarField>
174 alpha1_*rho1_*nuModel1_->nu()
175 + alpha2_*rho2_*nuModel2_->nu()
176 + alpha3_*rho3_*nuModel3_->nu()
189 return tmp<surfaceScalarField>
209 return tmp<surfaceScalarField>
218 )/(alpha1f*rho1_ + alpha2f*rho2_ + alpha3f*rho3_)
230 nuModel1_().
read(*
this)
231 && nuModel2_().
read(*
this)
232 && nuModel3_().
read(*
this)
235 nuModel1_->viscosityProperties().lookup(
"rho") >> rho1_;
236 nuModel2_->viscosityProperties().lookup(
"rho") >> rho2_;
237 nuModel3_->viscosityProperties().lookup(
"rho") >> rho3_;
incompressibleThreePhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
A class for managing temporary objects.
const dimensionSet dimDensity
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
bool read(const char *, int32_t &)
autoPtr< viscosityModel > nuModel1_
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
autoPtr< viscosityModel > nuModel2_
Constant dispersed-phase particle diameter model.
List< word > wordList
A List of words.
virtual bool read()=0
Read transportProperties dictionary.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Macros for easy insertion into run-time selection tables.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read()
Read base transportProperties dictionary.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void calcNu()
Calculate and return the laminar viscosity.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
autoPtr< viscosityModel > nuModel3_
stressControl lookup("compactNormalStress") >> compactNormalStress