velocityGroup Class Reference

Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of size class fractions. Intended for use with a population balance model to account for the evolution of a size distribution by means of coalescence, breakup, drift and nucleation. More...

Inheritance diagram for velocityGroup:
Collaboration diagram for velocityGroup:

Public Member Functions

 TypeName ("velocityGroup")
 Runtime type information. More...
 
 velocityGroup (const dictionary &diameterProperties, const phaseModel &phase)
 Construct from dictionary and phase. More...
 
virtual ~velocityGroup ()
 Destructor. More...
 
const wordpopBalName () const
 Return name of populationBalance this velocityGroup belongs to. More...
 
const populationBalanceModelpopBal () const
 Return the populationBalance this velocityGroup belongs to. More...
 
const PtrList< sizeGroup > & sizeGroups () const
 Return sizeGroups belonging to this velocityGroup. More...
 
PtrList< sizeGroup > & sizeGroups ()
 Return sizeGroups belonging to this velocityGroup. More...
 
virtual tmp< volScalarFieldd () const
 Get the diameter field. More...
 
virtual tmp< volScalarFieldAv () const
 Get the surface area per unit volume field. More...
 
virtual void correct ()
 Correct the model. More...
 
virtual bool read (const dictionary &diameterProperties)
 Read diameterProperties dictionary. More...
 
- Public Member Functions inherited from diameterModel
 TypeName ("diameterModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, diameterModel, dictionary,(const dictionary &diameterProperties, const phaseModel &phase),(diameterProperties, phase))
 
 diameterModel (const dictionary &diameterProperties, const phaseModel &phase)
 Construct from dictionary and phase. More...
 
virtual ~diameterModel ()
 Destructor. More...
 
const dictionarydiameterProperties () const
 Return the phase diameter properties dictionary. More...
 
const phaseModelphase () const
 Return the phase. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from diameterModel
static autoPtr< diameterModelNew (const dictionary &diameterProperties, const phaseModel &phase)
 Select from dictionary and phase. More...
 

Detailed Description

Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of size class fractions. Intended for use with a population balance model to account for the evolution of a size distribution by means of coalescence, breakup, drift and nucleation.

Usage
Excerpt from an exemplary phaseProperties dictionary:
diameterModel   velocityGroup;

velocityGroupCoeffs
{
    populationBalance    bubbles;

    shapeModel           spherical;

    sizeGroups
    (
        f1 {dSph  1e-3; value 1.0;}
        f2 {dSph  2e-3; value 0.0;}
        f3 {dSph  3e-3; value 0.0;}
        f4 {dSph  4e-3; value 0.0;}
        f5 {dSph  5e-3; value 0.0;}
        ...
    );
}
See also
Foam::diameterModels::sizeGroup Foam::diameterModels::populationBalanceModel
Source files

Definition at line 84 of file velocityGroup.H.

Constructor & Destructor Documentation

◆ velocityGroup()

velocityGroup ( const dictionary diameterProperties,
const phaseModel phase 
)

Construct from dictionary and phase.

Definition at line 138 of file velocityGroup.C.

References populationBalanceModel::groups::insert(), DimensionedField< Type, GeoMesh >::mesh(), populationBalanceModel::groups::New(), and diameterModel::phase().

Here is the call graph for this function:

◆ ~velocityGroup()

~velocityGroup ( )
virtual

Destructor.

Definition at line 172 of file velocityGroup.C.

Member Function Documentation

◆ TypeName()

TypeName ( "velocityGroup"  )

Runtime type information.

◆ popBalName()

const Foam::word & popBalName ( ) const
inline

Return name of populationBalance this velocityGroup belongs to.

Definition at line 31 of file velocityGroupI.H.

Referenced by fractal::fractal(), sizeGroup::iNew::operator()(), and sizeGroup::sizeGroup().

Here is the caller graph for this function:

◆ popBal()

Return the populationBalance this velocityGroup belongs to.

Definition at line 179 of file velocityGroup.C.

References objectRegistry::lookupObject(), and populationBalanceModel::mesh().

Referenced by nucleationSizeGroupFvScalarFieldSource::sourceValue().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sizeGroups() [1/2]

const Foam::PtrList< Foam::diameterModels::sizeGroup > & sizeGroups ( ) const
inline

Return sizeGroups belonging to this velocityGroup.

Definition at line 38 of file velocityGroupI.H.

Referenced by phaseChange::addToDriftRate(), phaseChange::precompute(), and reactionDriven::reactionDriven().

Here is the caller graph for this function:

◆ sizeGroups() [2/2]

Foam::PtrList< Foam::diameterModels::sizeGroup > & sizeGroups ( )
inline

Return sizeGroups belonging to this velocityGroup.

Definition at line 45 of file velocityGroupI.H.

◆ d()

Foam::tmp< Foam::volScalarField > d ( ) const
virtual

Get the diameter field.

Implements diameterModel.

Definition at line 191 of file velocityGroup.C.

◆ Av()

Foam::tmp< Foam::volScalarField > Av ( ) const
virtual

Get the surface area per unit volume field.

Implements diameterModel.

Definition at line 197 of file velocityGroup.C.

References sizeGroup::a(), Foam::dimLength, forAll, Foam::inv(), GeometricField< Type, PatchField, GeoMesh >::New(), tmp< T >::ref(), sizeGroup::x(), and Foam::Zero.

Here is the call graph for this function:

◆ correct()

◆ read()

bool read ( const dictionary diameterProperties)
virtual

Read diameterProperties dictionary.

Reimplemented from diameterModel.

Definition at line 265 of file velocityGroup.C.

References diameterModel::read().

Here is the call graph for this function:

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