Skip to content
Snippets Groups Projects
Commit 3885c493 authored by Henry's avatar Henry Committed by Andrew Heather
Browse files

driftFluxFoam: add support for alphaMax and d of the dispersed-phase

parent 78eaf506
No related branches found
No related tags found
No related merge requests found
...@@ -43,9 +43,23 @@ ...@@ -43,9 +43,23 @@
{ {
Info<< "Applying the previous iteration correction flux" << endl; Info<< "Applying the previous iteration correction flux" << endl;
#ifdef LTSSOLVE #ifdef LTSSOLVE
MULES::LTScorrect(alpha1, phiAlpha, tphiAlphaCorr0(), 1, 0); MULES::LTScorrect
(
alpha1,
phiAlpha,
tphiAlphaCorr0(),
mixture.alphaMax(),
0
);
#else #else
MULES::correct(alpha1, phiAlpha, tphiAlphaCorr0(), 1, 0); MULES::correct
(
alpha1,
phiAlpha,
tphiAlphaCorr0(),
mixture.alphaMax(),
0
);
#endif #endif
phiAlpha += tphiAlphaCorr0(); phiAlpha += tphiAlphaCorr0();
...@@ -79,9 +93,23 @@ ...@@ -79,9 +93,23 @@
volScalarField alpha10(alpha1); volScalarField alpha10(alpha1);
#ifdef LTSSOLVE #ifdef LTSSOLVE
MULES::LTScorrect(alpha1, tphiAlphaUn(), tphiAlphaCorr(), 1, 0); MULES::LTScorrect
(
alpha1,
tphiAlphaUn(),
tphiAlphaCorr(),
mixture.alphaMax(),
0
);
#else #else
MULES::correct(alpha1, tphiAlphaUn(), tphiAlphaCorr(), 1, 0); MULES::correct
(
alpha1,
tphiAlphaUn(),
tphiAlphaCorr(),
mixture.alphaMax(),
0
);
#endif #endif
// Under-relax the correction for all but the 1st corrector // Under-relax the correction for all but the 1st corrector
...@@ -100,9 +128,23 @@ ...@@ -100,9 +128,23 @@
phiAlpha = tphiAlphaUn; phiAlpha = tphiAlphaUn;
#ifdef LTSSOLVE #ifdef LTSSOLVE
MULES::explicitLTSSolve(alpha1, phi, phiAlpha, 1, 0); MULES::explicitLTSSolve
(
alpha1,
phi,
phiAlpha,
mixture.alphaMax(),
0
);
#else #else
MULES::explicitSolve(alpha1, phi, phiAlpha, 1, 0); MULES::explicitSolve
(
alpha1,
phi,
phiAlpha,
mixture.alphaMax(),
0
);
#endif #endif
} }
} }
......
...@@ -82,6 +82,13 @@ incompressibleTwoPhaseInteractingMixture ...@@ -82,6 +82,13 @@ incompressibleTwoPhaseInteractingMixture
rhod_("rho", dimDensity, muModel_->viscosityProperties().lookup("rho")), rhod_("rho", dimDensity, muModel_->viscosityProperties().lookup("rho")),
rhoc_("rho", dimDensity, nucModel_->viscosityProperties().lookup("rho")), rhoc_("rho", dimDensity, nucModel_->viscosityProperties().lookup("rho")),
dd_
(
"d",
dimLength,
muModel_->viscosityProperties().lookupOrDefault("d", 0.0)
),
alphaMax_(muModel_->viscosityProperties().lookupOrDefault("alphaMax", 1.0)),
U_(U), U_(U),
phi_(phi), phi_(phi),
...@@ -118,6 +125,20 @@ bool Foam::incompressibleTwoPhaseInteractingMixture::read() ...@@ -118,6 +125,20 @@ bool Foam::incompressibleTwoPhaseInteractingMixture::read()
muModel_->viscosityProperties().lookup("rho") >> rhod_; muModel_->viscosityProperties().lookup("rho") >> rhod_;
nucModel_->viscosityProperties().lookup("rho") >> rhoc_; nucModel_->viscosityProperties().lookup("rho") >> rhoc_;
dd_ = dimensionedScalar
(
"d",
dimLength,
muModel_->viscosityProperties().lookupOrDefault("d", 0)
);
alphaMax_ =
muModel_->viscosityProperties().lookupOrDefault
(
"alphaMax",
1.0
);
return true; return true;
} }
else else
......
...@@ -69,6 +69,12 @@ protected: ...@@ -69,6 +69,12 @@ protected:
dimensionedScalar rhod_; dimensionedScalar rhod_;
dimensionedScalar rhoc_; dimensionedScalar rhoc_;
//- Optional diameter of the dispersed phase particles
dimensionedScalar dd_;
//- Optional maximum dispersed phase-fraction (e.g. packing limit)
scalar alphaMax_;
const volVectorField& U_; const volVectorField& U_;
const surfaceScalarField& phi_; const surfaceScalarField& phi_;
...@@ -121,6 +127,19 @@ public: ...@@ -121,6 +127,19 @@ public:
return rhoc_; return rhoc_;
}; };
//- Return the diameter of the dispersed-phase particles
const dimensionedScalar& dd() const
{
return dd_;
}
//- Optional maximum phase-fraction (e.g. packing limit)
// Defaults to 1
scalar alphaMax() const
{
return alphaMax_;
}
//- Return const-access to the mixture velocity //- Return const-access to the mixture velocity
const volVectorField& U() const const volVectorField& U() const
{ {
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment