Loading...
Searching...
No Matches
massRosinRammler Class Reference

Particle-size distribution model wherein random samples are drawn from the two-parameter Rosin-Rammler (Weibull) probability density function corrected to take into account varying number of particles per parcels for fixed-mass parcels. More...

#include <massRosinRammler.H>

Inheritance diagram for massRosinRammler:
Collaboration diagram for massRosinRammler:

Public Member Functions

 TypeName ("massRosinRammler")
 Runtime type information.
 massRosinRammler (const dictionary &dict, Random &rndGen)
 Construct from components.
 massRosinRammler (const massRosinRammler &p)
 Copy construct.
virtual autoPtr< distributionModelclone () const
 Construct and return a clone.
void operator= (const massRosinRammler &)=delete
 No copy assignment.
virtual ~massRosinRammler ()=default
 Destructor.
virtual scalar sample () const
 Sample the distribution.
virtual scalar meanValue () const
 Return the theoretical mean of the distribution.
Public Member Functions inherited from distributionModel
 TypeName ("distributionModel")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, distributionModel, dictionary,(const dictionary &dict, Random &rndGen),(dict, rndGen))
 Declare runtime constructor selection table.
 distributionModel (const word &name, const dictionary &dict, Random &rndGen)
 Construct from dictionary.
 distributionModel (const distributionModel &p)
 Copy construct.
virtual ~distributionModel ()=default
 Destructor.
virtual scalar minValue () const
 Return the minimum of the distribution.
virtual scalar maxValue () const
 Return the maximum of the distribution.

Additional Inherited Members

Static Public Member Functions inherited from distributionModel
static autoPtr< distributionModelNew (const dictionary &dict, Random &rndGen)
 Selector.
Protected Member Functions inherited from distributionModel
virtual void check () const
 Check that the distribution model is valid.
Protected Attributes inherited from distributionModel
const dictionary distributionModelDict_
 Coefficients dictionary.
RandomrndGen_
 Reference to the random number generator.
scalar minValue_
 Minimum of the distribution.
scalar maxValue_
 Maximum of the distribution.

Detailed Description

Particle-size distribution model wherein random samples are drawn from the two-parameter Rosin-Rammler (Weibull) probability density function corrected to take into account varying number of particles per parcels for fixed-mass parcels.

"There is a factor of \_form#520 difference between the size distribution probability density function of individual particles and modeled parcels of particles, \_form#521, because the submodel parameter, \_form#522 (number of particles per parcel) is based on a fixed mass per parcel which weights the droplet distribution by a factor proportional to \_form#523 (i.e. \iline 34 \_form#524@_fakenl)." (YHD:p. 1374-1375).

massRosinRammler should be preferred over RosinRammler when parcelBasisType is based on mass:

    parcelBasisType mass;

The doubly-truncated mass-corrected Rosin-Rammler probability density function:

\‍[    f(x; \lambda, n, A, B) =
            x^3
            \frac{n}{\lambda}
            \left( \frac{x}{\lambda} \right)^{n-1}
            \frac{
                \exp\{ -(\frac{x}{\lambda})^n \}
            }
            {
                \exp\{ -(\frac{A}{\lambda})^n \}
              - \exp\{ -(\frac{B}{\lambda})^n \}
            }
\‍]

where

$        f(x; \lambda, n, A, B) $=Rosin-Rammler probability density function
$        \lambda   $=Scale parameter
$        n         $=Shape parameter
$        x         $=Sample
$        A         $=Minimum of the distribution
$        B         $=Maximum of the distribution

Constraints:

  • $ \lambda > 0 $
  • $ n > 0 $

Random samples are generated by the inverse transform sampling technique by using the quantile function of the doubly-truncated two-parameter mass-corrected Rosin-Rammler (Weibull) probability density function:

\‍[        d = \lambda \, Q(a, u)^{1/n}
\‍]

with

\‍[        a = \frac{3}{n} + 1
\‍]

where $ Q(z_1, z_2) $ is the inverse of regularised lower incomplete gamma function, and $ u $ is a sample drawn from the uniform probability density function on the interval $ (x, y) $:

\‍[        x = \gamma \left( a, \frac{A}{\lambda}^n \right)
\‍]

\‍[        y = \gamma \left( a, \frac{B}{\lambda}^n \right)
\‍]

where $ \gamma(z_1, z_2) $ is the lower incomplete gamma function.

Reference:

        Standard model (tag:YHD):
            Yoon, S. S., Hewson, J. C., DesJardin, P. E.,
            Glaze, D. J., Black, A. R., & Skaggs, R. R. (2004).
            Numerical modeling and experimental measurements of a high speed
            solid-cone water spray for use in fire suppression applications.
            International Journal of Multiphase Flow, 30(11), 1369-1388.
            DOI:10.1016/j.ijmultiphaseflow.2004.07.006
Usage
Minimal example by using constant/<CloudProperties>:
subModels
{
    injectionModels
    {
        <name>
        {
            ...
            parcelBasisType  mass;
            ...

            sizeDistribution
            {
                type        massRosinRammler;
                massRosinRammlerDistribution
                {
                    lambda      <scaleParameterValue>;
                    n           <shapeParameterValue>;
                    minValue    <minValue>;
                    maxValue    <maxValue>;
                }
            }
        }
    }
}

where the entries mean:

Property Description Type Reqd Deflt
type Type name: massRosinRammler word yes -
massRosinRammlerDistribution Distribution settings dict yes -
lambda Scale parameter scalar yes -
n Shape parameter scalar yes -
minValue Minimum of the distribution scalar yes -
maxValue Maximum of the distribution scalar yes -
Note
  • In the current context, lambda (or d) is called "characteristic droplet size", and n "empirical dimensionless constant to specify the distribution width, sometimes referred to as the dispersion coefficient." (YHD:p. 1374).
Source files

Definition at line 245 of file massRosinRammler.H.

Constructor & Destructor Documentation

◆ massRosinRammler() [1/2]

massRosinRammler ( const dictionary & dict,
Random & rndGen )

Construct from components.

Definition at line 39 of file massRosinRammler.C.

References dict, distributionModel::distributionModel(), distributionModel::distributionModelDict_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, and rndGen.

Referenced by clone(), massRosinRammler(), operator=(), and TypeName().

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

◆ massRosinRammler() [2/2]

massRosinRammler ( const massRosinRammler & p)

Copy construct.

Definition at line 62 of file massRosinRammler.C.

References distributionModel::distributionModel(), massRosinRammler(), and p.

Here is the call graph for this function:

◆ ~massRosinRammler()

virtual ~massRosinRammler ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "massRosinRammler" )

Runtime type information.

References dict, massRosinRammler(), p, and rndGen.

Here is the call graph for this function:

◆ clone()

virtual autoPtr< distributionModel > clone ( ) const
inlinevirtual

Construct and return a clone.

Implements distributionModel.

Definition at line 285 of file massRosinRammler.H.

References massRosinRammler().

Here is the call graph for this function:

◆ operator=()

void operator= ( const massRosinRammler & )
delete

No copy assignment.

References massRosinRammler().

Here is the call graph for this function:

◆ sample()

Foam::scalar sample ( ) const
virtual

Sample the distribution.

Implements distributionModel.

Definition at line 75 of file massRosinRammler.C.

References Foam::Math::incGamma_P(), Foam::Math::invIncGamma(), distributionModel::maxValue_, distributionModel::minValue_, Foam::pow(), distributionModel::rndGen_, and x.

Here is the call graph for this function:

◆ meanValue()

Foam::scalar meanValue ( ) const
virtual

Return the theoretical mean of the distribution.

Implements distributionModel.

Definition at line 95 of file massRosinRammler.C.


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