Go to the documentation of this file.
33 namespace phaseChangeTwoPhaseMixtures
48 phaseChangeTwoPhaseMixture(typeName,
U,
phi),
50 UInf_(phaseChangeTwoPhaseMixtureCoeffs_.
lookup(
"UInf")),
51 tInf_(phaseChangeTwoPhaseMixtureCoeffs_.
lookup(
"tInf")),
52 Cc_(phaseChangeTwoPhaseMixtureCoeffs_.
lookup(
"Cc")),
53 Cv_(phaseChangeTwoPhaseMixtureCoeffs_.
lookup(
"Cv")),
55 p0_(
"0",
pSat().dimensions(), 0.0),
57 mcCoeff_(Cc_/(0.5*
sqr(UInf_)*tInf_)),
71 return Pair<tmp<volScalarField> >
84 return Pair<tmp<volScalarField> >
86 mcCoeff_*(1.0 - limitedAlpha1)*
pos(
p -
pSat()),
87 (-mvCoeff_)*limitedAlpha1*
neg(
p -
pSat())
100 phaseChangeTwoPhaseMixtureCoeffs_ = subDict(
type() +
"Coeffs");
102 phaseChangeTwoPhaseMixtureCoeffs_.lookup(
"UInf") >> UInf_;
103 phaseChangeTwoPhaseMixtureCoeffs_.lookup(
"tInf") >> tInf_;
104 phaseChangeTwoPhaseMixtureCoeffs_.lookup(
"Cc") >> Cc_;
105 phaseChangeTwoPhaseMixtureCoeffs_.lookup(
"Cv") >> Cv_;
107 mcCoeff_ = Cc_/(0.5*
sqr(UInf_)*tInf_);
108 mvCoeff_ = Cv_*
rho1()/(0.5*
sqr(UInf_)*tInf_*
rho2());
virtual void correct()
Correct the Merkle phaseChange model.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Merkle(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual bool read()=0
Read the transportProperties dictionary and update.
const dimensionedScalar & pSat
Macros for easy insertion into run-time selection tables.
virtual void correct()
Correct the Kunz phaseChange model.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An ordered pair of two objects of type <T> with first() and second() elements.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
virtual bool read()
Read the transportProperties dictionary and update.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar pos(const dimensionedScalar &ds)