velocityGroup Class Reference

This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispersed bubbly or particulate flows. It can hold any number of sizeGroups from which the Sauter mean diameter is calculated. It can also be used as a diameterModel without a populationBalance and would then behave like a constantDiameter model. In this case, some arbitrary name must be entered for the populationBalance keyword. 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 components. More...
 
virtual ~velocityGroup ()
 Destructor. More...
 
const wordpopBalName () const
 Return name of populationBalance this velocityGroup belongs to. More...
 
const volScalarFieldf () const
 Return reference field for sizeGroup's. More...
 
const dimensionedScalarformFactor () const
 Return the form factor. More...
 
const PtrList< sizeGroup > & sizeGroups () const
 Return sizeGroups belonging to this velocityGroup. More...
 
const tmp< fv::convectionScheme< scalar > > & mvConvection () const
 Return const-reference to multivariate convectionScheme. More...
 
const volScalarFielddmdt () const
 Return const-reference to the mass transfer rate. More...
 
volScalarFielddmdtRef ()
 Return reference to the mass transfer rate. More...
 
void preSolve ()
 Corrections before populationBalanceModel::solve() More...
 
void postSolve ()
 Corrections after populationBalanceModel::solve() More...
 
virtual bool read (const dictionary &diameterProperties)
 Read diameterProperties dictionary. More...
 
virtual tmp< volScalarFieldd () const
 Return the Sauter-mean diameter. More...
 
- Public Member Functions inherited from diameterModel
 TypeName ("diameterModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, diameterModel, dictionary,(const dictionary &dict, const phaseModel &phase),(dict, phase))
 
 diameterModel (const dictionary &dict, const phaseModel &phase)
 
virtual ~diameterModel ()
 Destructor. More...
 
 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)
 
virtual ~diameterModel ()
 Destructor. More...
 
const dictionarydiameterProperties () const
 Return the phase diameter properties dictionary. More...
 
const phaseModelphase () const
 Return the phase. More...
 
virtual void correct ()
 Correct the diameter field. More...
 
 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)
 
virtual ~diameterModel ()
 Destructor. More...
 
const dictionarydiameterProperties () const
 Return the phase diameter properties dictionary. More...
 
const phaseModelphase () const
 Return the phase. More...
 
virtual void correct ()
 Correct the diameter field. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from diameterModel
static autoPtr< diameterModelNew (const dictionary &dict, const phaseModel &phase)
 
static autoPtr< diameterModelNew (const dictionary &diameterProperties, const phaseModel &phase)
 
static autoPtr< diameterModelNew (const dictionary &diameterProperties, const phaseModel &phase)
 
- Protected Attributes inherited from diameterModel
const dictionarydict_
 
const phaseModelphase_
 
dictionary diameterProperties_
 

Detailed Description

This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispersed bubbly or particulate flows. It can hold any number of sizeGroups from which the Sauter mean diameter is calculated. It can also be used as a diameterModel without a populationBalance and would then behave like a constantDiameter model. In this case, some arbitrary name must be entered for the populationBalance keyword.

Usage
Property Description
populationBalance Name of the corresponding populationBalance
formFactor Form factor for converting diameter into volume
sizeGroups List of sizeGroups

Example

    diameterModel   velocityGroup;
    velocityGroupCoeffs
    {
        populationBalance    bubbles;

        formFactor      0.5235987756;

        sizeGroups
        (
            f0{d  1.00e-3; value 0;}
            f1{d  1.08e-3; value 0;}
            f2{d  1.16e-3; value 0.25;}
            f3{d  1.25e-3; value 0.5;}
            f4{d  1.36e-3; value 0.25;}
            f5{d  1.46e-3; value 0;}
            ...
        );
    }
See also
Foam::diameterModels::sizeGroup Foam::diameterModels::populationBalanceModel
Source files

Definition at line 105 of file velocityGroup.H.

Constructor & Destructor Documentation

◆ velocityGroup()

velocityGroup ( const dictionary diameterProperties,
const phaseModel phase 
)

Construct from components.

◆ ~velocityGroup()

virtual ~velocityGroup ( )
virtual

Destructor.

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 29 of file velocityGroupI.H.

◆ f()

const Foam::volScalarField & f ( ) const
inline

Return reference field for sizeGroup's.

Definition at line 36 of file velocityGroupI.H.

◆ formFactor()

const Foam::dimensionedScalar & formFactor ( ) const
inline

Return the form factor.

Definition at line 43 of file velocityGroupI.H.

◆ sizeGroups()

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

Return sizeGroups belonging to this velocityGroup.

Definition at line 50 of file velocityGroupI.H.

◆ mvConvection()

const Foam::tmp< Foam::fv::convectionScheme< Foam::scalar > > & mvConvection ( ) const
inline

Return const-reference to multivariate convectionScheme.

Definition at line 57 of file velocityGroupI.H.

References velocityGroup::dmdt().

Here is the call graph for this function:

◆ dmdt()

const Foam::volScalarField & dmdt ( ) const
inline

Return const-reference to the mass transfer rate.

Definition at line 64 of file velocityGroupI.H.

Referenced by velocityGroup::mvConvection().

Here is the caller graph for this function:

◆ dmdtRef()

Foam::volScalarField & dmdtRef ( )
inline

Return reference to the mass transfer rate.

Definition at line 70 of file velocityGroupI.H.

◆ preSolve()

void preSolve ( )

Corrections before populationBalanceModel::solve()

◆ postSolve()

void postSolve ( )

◆ read()

virtual bool read ( const dictionary diameterProperties)
virtual

Read diameterProperties dictionary.

Implements diameterModel.

◆ d()

virtual tmp<volScalarField> d ( ) const
virtual

Return the Sauter-mean diameter.

Implements diameterModel.


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