BrownianCollisions Class Reference

Model describing coagulation due to Brownian motion. Utilises collisional diameters and the Cunningham slip correction. The slip correction coefficient is implemented in the following form: More...

Inheritance diagram for BrownianCollisions:
Collaboration diagram for BrownianCollisions:

Public Member Functions

 TypeName ("BrownianCollisions")
 Runtime type information. More...
 
 BrownianCollisions (const populationBalanceModel &popBal, const dictionary &dict)
 
virtual ~BrownianCollisions ()
 Destructor. More...
 
virtual void precompute ()
 Precompute 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 describing coagulation due to Brownian motion. Utilises collisional diameters and the Cunningham slip correction. The slip correction coefficient is implemented in the following form:

\[ C_{c_i} = 1 + \lambda [A_1 + A_2 \exp(-A_3 d_i/\lambda)]/d_i\,. \]

The coefficients default to the values proposed by Davis (1945). The mean free path is computed by

\[ \lambda = \frac{kT}{\sqrt{2} \pi p \sigma^{2}}\,. \]

$ A_1 $ = Coefficient [-]
$ A_2 $ = Coefficient [-]
$ A_3 $ = Coefficient [-]
$ \sigma $ = Lennard-Jones parameter [m]

Reference:

        Davies, C. N. (1945).
        Definitive equations for the fluid resistance of spheres.
        Proceedings of the Physical Society, 57(4), 259.
Usage
Property Description Required Default value
A1 Coefficient A1 no 2.514
A2 Coefficient A2 no 0.8
A3 Coefficient A2 no 0.55
sigma Lennard-Jones parameter yes none
Source files

Definition at line 129 of file BrownianCollisions.H.

Constructor & Destructor Documentation

◆ BrownianCollisions()

BrownianCollisions ( const populationBalanceModel popBal,
const dictionary dict 
)

Definition at line 55 of file BrownianCollisions.C.

◆ ~BrownianCollisions()

virtual ~BrownianCollisions ( )
inlinevirtual

Destructor.

Definition at line 166 of file BrownianCollisions.H.

Member Function Documentation

◆ TypeName()

TypeName ( "BrownianCollisions"  )

Runtime type information.

◆ precompute()

void precompute ( )
virtual

Precompute diameter independent expressions.

Reimplemented from coalescenceModel.

Definition at line 88 of file BrownianCollisions.C.

References k, p, Foam::constant::mathematical::pi(), Foam::sqr(), Foam::sqrt(), and Foam::T().

Here is the call graph for this function:

◆ addToCoalescenceRate()

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

Add to coalescenceRate.

Implements coalescenceModel.

Definition at line 98 of file BrownianCollisions.C.

References sizeGroup::d(), Foam::exp(), k, Foam::constant::physicoChemical::mu, and Foam::T().

Here is the call graph for this function:

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