KochFriedlanderSintering Class Reference

Sintering model of Koch and Friedlander (1990). The characteristic time for sintering is given by. More...

Inheritance diagram for KochFriedlanderSintering:
Collaboration diagram for KochFriedlanderSintering:

Public Member Functions

 TypeName ("KochFriedlanderSintering")
 Runtime type information. More...
 
 KochFriedlanderSintering (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
virtual bool addsSupToField (const word &fieldName) const
 Return true if the fvModel adds a source term to the given. More...
 
virtual tmp< volScalarField::Internaltau (const volScalarField::Internal &kappa) const
 Return the characteristic time for sintering. More...
 
void addSup (const volScalarField &alphaFi, const volScalarField &rho, const volScalarField &kappa, fvMatrix< scalar > &eqn) const
 Add a source term to the surface-area-volume-ratio equation. More...
 
virtual bool movePoints ()
 Update for mesh motion. More...
 
virtual void topoChange (const polyTopoChangeMap &)
 Update topology using the given map. More...
 
virtual void mapMesh (const polyMeshMap &)
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)
 Redistribute or update using the given distribution map. More...
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 
- Public Member Functions inherited from fvModel
 TypeName ("fvModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
 
 fvModel (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
autoPtr< fvModelclone () const
 Return clone. More...
 
virtual ~fvModel ()
 Destructor. More...
 
const wordname () const
 Return const access to the source name. More...
 
const wordkeyword () const
 Return name as the keyword. More...
 
const fvMeshmesh () const
 Return const access to the mesh database. More...
 
const dictionarycoeffs (const dictionary &) const
 Return the coefficients sub-dictionary. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the fvModel adds source term. More...
 
virtual scalar maxDeltaT () const
 Return the maximum time-step for stable operation. More...
 
virtual void addSup (fvMatrix< scalar > &eqn) const
 Add a source term to a field-less proxy equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const VolField< Type > &eqnField) const
 Add a source term to an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const VolField< Type > &field) const
 Return source for an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 Return source for a compressible equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for a compressible equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 Return source for an equation with a second time derivative. More...
 
virtual void preUpdateMesh ()
 Prepare for mesh update. More...
 
virtual void correct ()
 Correct the fvModel. More...
 
virtual bool write (const bool write=true) const
 Write fvModel data. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
Foam::dimensionSet sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
const Foam::wordfieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 
template<class AlphaRhoFieldType >
const Foam::wordfieldName (const AlphaRhoFieldType &alphaRhoField)
 
template<class Type , class ... AlphaRhoFieldTypes>
Foam::tmp< Foam::fvMatrix< Type > > sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 

Additional Inherited Members

- Static Public Member Functions inherited from fvModel
static const dictionarycoeffs (const word &modelType, const dictionary &)
 Return the coefficients sub-dictionary for a given model type. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static dimensionSet sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 Return the dimensions of the matrix of a source term. More...
 
static const dimensionSetsourceDims (const dimensionSet &ds)
 Return the dimensions of the matrix of a source term (base. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static const wordfieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 Return the name of the field associated with a source term. More...
 
template<class AlphaRhoFieldType >
static const wordfieldName (const AlphaRhoFieldType &alphaRhoField)
 Return the name of the field associated with a source term (base. More...
 
static const wordfieldName ()
 Return the name of the field associated with a source term. Special. More...
 
static autoPtr< fvModelNew (const word &name, const fvMesh &mesh, const dictionary &dict)
 Return a reference to the selected fvModel. More...
 
- Static Public Attributes inherited from fvModel
static const wordHashSet keywords
 The keywords read by this class. More...
 
- Protected Member Functions inherited from fvModel
template<class Type >
void addSupType (const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to an equation. More...
 
template<class Type >
void addSupType (const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to a compressible equation. More...
 
template<class Type >
void addSupType (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to a phase equation. More...
 
template<class Type , class ... AlphaRhoFieldTypes>
tmp< fvMatrix< Type > > sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
 Return a source for an equation. More...
 

Detailed Description

Sintering model of Koch and Friedlander (1990). The characteristic time for sintering is given by.

\[ \tau = c_s d_{p_i}^n T^m \exp(T_a/T \cdot [1 - d_{p,min}/d_{p_i}])\;. \]

Note that the correction factor in the exponential function can be eliminated by setting $d_{p,min}$ to zero which is done by default.

Reference:

    Koch, W., & Friedlander, S. K. (1990).
    The effect of particle coalescence on the surface area of a coagulating
    aerosol.
    Journal of Colloid and Interface Science, 140(2), 419-427.
Usage
Property Description Required Default value
Cs Sintering time coefficient yes none
n Particle diameter exponent yes none
m Temperature exponent yes none
Ta Activation temperature yes none
dpMin Minimum primary particle diameter no 0

Example usage:

    sintering
    {
        type            KochFriedlanderSintering;
        libs            ("libmultiphaseEulerFvModels.so");

        populationBalance aggregates;

        Cs              8.3e24;
        n               4.0;
        m               1.0;
        Ta              3700.0;
    }
Source files

Definition at line 122 of file KochFriedlanderSintering.H.

Constructor & Destructor Documentation

◆ KochFriedlanderSintering()

KochFriedlanderSintering ( const word name,
const word modelType,
const fvMesh mesh,
const dictionary dict 
)

Construct from explicit source name and mesh.

Definition at line 71 of file KochFriedlanderSintering.C.

References fvModel::coeffs(), dict, fractal::fld(), forAll, HashTable< T, Key, Hash >::insert(), IOobject::name(), and populationBalanceModel::sizeGroups().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "KochFriedlanderSintering"  )

Runtime type information.

◆ addsSupToField()

bool addsSupToField ( const word fieldName) const
virtual

Return true if the fvModel adds a source term to the given.

field's transport equation

Reimplemented from fvModel.

Definition at line 112 of file KochFriedlanderSintering.C.

◆ tau()

Foam::tmp< Foam::volScalarField::Internal > tau ( const volScalarField::Internal kappa) const
virtual

Return the characteristic time for sintering.

Definition at line 122 of file KochFriedlanderSintering.C.

References sizeGroup::dSph(), Foam::exp(), Foam::constant::electromagnetic::kappa, Foam::max(), dimensioned< Type >::name(), sizeGroup::phase(), Foam::pow(), basicThermo::T(), Foam::T(), and phaseModel::thermo().

Here is the call graph for this function:

◆ addSup()

void addSup ( const volScalarField alphaFi,
const volScalarField rho,
const volScalarField kappa,
fvMatrix< scalar > &  eqn 
) const

Add a source term to the surface-area-volume-ratio equation.

Definition at line 138 of file KochFriedlanderSintering.C.

References sizeGroup::dSph(), Foam::constant::electromagnetic::kappa, dimensioned< Type >::name(), Foam::R(), rho, Foam::fvm::Sp(), and tau.

Here is the call graph for this function:

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 155 of file KochFriedlanderSintering.C.

◆ topoChange()

void topoChange ( const polyTopoChangeMap )
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 161 of file KochFriedlanderSintering.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap )
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 165 of file KochFriedlanderSintering.C.

◆ distribute()

void distribute ( const polyDistributionMap )
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 169 of file KochFriedlanderSintering.C.

◆ read()

bool read ( const dictionary dict)
virtual

Read source dictionary.

Reimplemented from fvModel.

Definition at line 173 of file KochFriedlanderSintering.C.

References dict, and fvModel::read().

Here is the call graph for this function:

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