PrinceBlanch Class Reference

Model of Prince and Blanch (1990). The coalescence rate is calculated by. More...

Inheritance diagram for PrinceBlanch:
Collaboration diagram for PrinceBlanch:

Public Member Functions

 TypeName ("PrinceBlanch")
 Runtime type information. More...
 
 PrinceBlanch (const populationBalanceModel &popBal, const dictionary &dict)
 
virtual ~PrinceBlanch ()
 Destructor. More...
 
virtual void correct ()
 Correct diameter independent expressions. More...
 
virtual void addToCoalescenceRate (volScalarField &coalescenceRate, const label i, const label j)
 Add to coalescenceRate. More...
 
- Public Member Functions inherited from coalescenceModel
 TypeName ("coalescenceModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, coalescenceModel, dictionary,(const populationBalanceModel &popBal, const dictionary &dict),(popBal, dict))
 
 coalescenceModel (const populationBalanceModel &popBal, const dictionary &dict)
 
autoPtr< coalescenceModelclone () const
 
virtual ~coalescenceModel ()
 Destructor. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from coalescenceModel
static autoPtr< coalescenceModelNew (const word &type, const populationBalanceModel &popBal, const dictionary &dict)
 
- Protected Attributes inherited from coalescenceModel
const populationBalanceModelpopBal_
 Reference to the populationBalanceModel. More...
 

Detailed Description

Model of Prince and Blanch (1990). The coalescence rate is calculated by.

\[ \left( \theta_{ij}^{T} + \theta_{ij}^{B} + \theta_{ij}^{LS} \right) \lambda_{ij}\;, \]

with the coalescence efficiency

\[ \lambda_{ij} = \mathrm{exp} \left( - \sqrt{\frac{r_{ij}^3 \rho_c}{16 \sigma}} \mathrm{ln} \left(\frac{h_0}{h_f}\right) \epsilon_c^{1/3}/r_{ij}^{2/3} \right)\;, \]

the turbulent collision rate

\[ \theta_{ij}^{T} = C_1 \pi (d_i + d_j)^{2} \epsilon_c^{1/3} \sqrt{d_{i}^{2/3} + d_{j}^{2/3}}\;, \]

the buoyancy-driven collision rate

\[ \theta_{ij}^{B} = S_{ij} \left| u_{ri} - u_{rj} \right|\;, \]

and the laminar shear collision rate

\[ \theta_{ij}^{LS} = \frac{1}{6} \left(d_i + d_j\right)^{3} \gamma_c\;. \]

The rise velocity of bubble i is calculated by

\[ u_{ri} = \sqrt{2.14 \sigma / \left(\rho_c d_i \right) + 0.505 g d_i}\;, \]

the equivalent radius by

\[ r_{ij} = \left( \frac{1}{d_i} + \frac{1}{d_j} \right)^{-1}\;, \]

the collision cross sectional area by

\[ S_{ij} = \frac{\pi}{4} \left(d_i + d_j\right)^{2}\;, \]

and the shear strain rate by

\[ \dot{\gamma_{b}} = \mathrm{mag}(\mathrm{symm}(\mathrm{grad}(U_c)))\;. \]

Note that in equation 2, the bubble radius has been substituted by the bubble diameter. Also the expression for the equivalent radius r_ij (equation 19 in the paper of Prince and Blanch (1990)) was corrected.

$ \theta_{ij}^{T} $ = Turbulent collision rate [m3/s]
$ \theta_{ij}^{B} $ = Buoyancy-driven collision rate [m3/s]
$ \theta_{ij}^{LS} $ = Laminar shear collision rate [m3/s]
$ \lambda_{ij} $ = Coalescence efficiency [-]
$ r_{ij} $ = Equivalent radius [m]
$ \rho_c $ = Density of continuous phase [kg/m3]
$ \sigma $ = Surface tension [N/m]
$ h_0 $ = Initial film thickness [m]
$ h_f $ = Critical film thickness [m]
$ \epsilon_c $ = Continuous phase turbulent dissipation rate [m2/s3]
$ d_i $ = Diameter of bubble i [m]
$ d_j $ = Diameter of bubble j [m]
$ u_{ri} $ = Rise velocity of bubble i [m/s]
$ S_{ij} $ = Collision cross sectional area [m2]
$ g $ = Gravitational constant [m/s2]
$ \gamma_c $ = Continuous phase shear strain rate [1/s]
$ U_c $ = Continuous phase velocity field [m/s]

Reference:

        Prince, M. J., & Blanch, H. W. (1990).
        Bubble coalescence and break‐up in air‐sparged bubble columns.
        AIChE journal, 36(10), 1485-1499.
Usage
Property Description Required Default value
C1 coefficient C1 no 0.089
h0 initial film thickness no 1e-4m
hf critical film thickness no 1e-8m
turbulence Switch for collisions due to turbulence yes none
buoyancy Switch for collisions due to buoyancy yes none
laminarShear Switch for collisions due to laminar shear yes none
Source files

Definition at line 259 of file PrinceBlanch.H.

Constructor & Destructor Documentation

◆ PrinceBlanch()

PrinceBlanch ( const populationBalanceModel popBal,
const dictionary dict 
)

◆ ~PrinceBlanch()

virtual ~PrinceBlanch ( )
inlinevirtual

Destructor.

Definition at line 302 of file PrinceBlanch.H.

References PrinceBlanch::addToCoalescenceRate(), and PrinceBlanch::correct().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "PrinceBlanch"  )

Runtime type information.

◆ correct()

virtual void correct ( )
virtual

Correct diameter independent expressions.

Reimplemented from coalescenceModel.

Referenced by PrinceBlanch::~PrinceBlanch().

Here is the caller graph for this function:

◆ addToCoalescenceRate()

virtual void addToCoalescenceRate ( volScalarField coalescenceRate,
const label  i,
const label  j 
)
virtual

Add to coalescenceRate.

Implements coalescenceModel.

Referenced by PrinceBlanch::~PrinceBlanch().

Here is the caller graph for this function:

The documentation for this class was generated from the following file: