constant
Constant dispersed-phase particle diameter model.
1 2 3 4 5 6 7 8 9 10 |
water { diameterModel constant; constantCoeffs { d 1e-4; } residualAlpha 1e-6; } |
68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 |
Foam::tmp<Foam::volScalarField> Foam::diameterModels::constant::d() const { return tmp<Foam::volScalarField> ( new volScalarField ( IOobject ( "d", phase_.U().time().timeName(), phase_.U().mesh(), IOobject::NO_READ, IOobject::NO_WRITE, false ), phase_.U().mesh(), d_ ) ); } |
isothermal
Isothermal dispersed-phase particle diameter model.
\begin{align}
p \times \frac{4}{3} \pi \left( \frac{d}{2} \right)^3 &= p_0 \times \frac{4}{3} \pi \left( \frac{d_0}{2} \right)^3 \\
d &= d_0 \left( \frac{p_0}{p} \right)^{1/3} \tag{1} \label{eq:levmdiffusion}
\end{align}
1 2 3 4 5 6 7 8 9 10 11 |
air { diameterModel isothermal; isothermalCoeffs { d0 3e-3; p0 1e5; } residualAlpha 1e-6; } |
69 70 71 72 73 74 75 76 77 |
Foam::tmp<Foam::volScalarField> Foam::diameterModels::isothermal::d() const { const volScalarField& p = phase_.U().db().lookupObject<volScalarField> ( "p" ); return d0_*cbrt(p0_/p); } |
IATE
IATE (Interfacial Area Transport Equation) bubble diameter model.
Solves for the interfacial curvature per unit volume of the phase rather than interfacial area per unit volume to avoid stability issues relating to the consistency requirements between the phase fraction and interfacial area per unit volume.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 |
air { diameterModel IATE; IATECoeffs { dMax 1e-2; dMin 1e-4; residualAlpha 1e-6; sources ( wakeEntrainmentCoalescence { Cwe 0.002; } randomCoalescence { Crc 0.04; C 3; alphaMax 0.75; } turbulentBreakUp { Cti 0.085; WeCr 6; } ); } residualAlpha 1e-6; } |
The following bubble interaction mechanisms are considered in this model.
-
- wakeEntrainmentCoalescence
Bubble coalescence by wake entrainment - randomCoalescence
Bubble coalescence by random collision - turbulentBreakUp
Bubble break-up due to the impact of turbulent eddies
\begin{align}
\phi_{TI} = \frac{C_{ti}}{3} \frac{U_t}{d} \sqrt{1 – \frac{We_{Cr}}{We}} e^{-\frac{We_{Cr}}{We}}
\end{align}- \(We\): Weber number [-]
- \(U_t\): Turbulent fluctuation velocity [m/s]
115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const{return max(6/max(kappai_, 6/dMax_), dMin_);}void Foam::diameterModels::IATE::correct(){// Initialise the accumulated source term to the dilatation effectvolScalarField R(((1.0/3.0)/max(fvc::average(phase_ + phase_.oldTime()),residualAlpha_))*(fvc::ddt(phase_) + fvc::div(phase_.alphaPhi())));// Accumulate the run-time selectable sourcesforAll(sources_, j){R -= sources_[j].R();}fv::options& fvOptions(fv::options::New(phase_.mesh()));// Construct the interfacial curvature equationfvScalarMatrix kappaiEqn(fvm::ddt(kappai_) + fvm::div(phase_.phi(), kappai_)- fvm::Sp(fvc::div(phase_.phi()), kappai_)==- fvm::SuSp(R, kappai_)//+ Rph() // Omit the nucleation/condensation term+ fvOptions(kappai_));kappaiEqn.relax();fvOptions.constrain(kappaiEqn);kappaiEqn.solve();// Update the Sauter-mean diameterd_ = dsm();} - wakeEntrainmentCoalescence