massRosinRammler.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "massRosinRammler.H"
30 #include "MathFunctions.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace distributionModels
38 {
39  defineTypeNameAndDebug(massRosinRammler, 0);
40  addToRunTimeSelectionTable(distributionModel, massRosinRammler, dictionary);
41 }
42 }
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const dictionary& dict,
49  Random& rndGen
50 )
51 :
52  distributionModel(typeName, dict, rndGen),
53  lambda_(distributionModelDict_.getCompat<scalar>("lambda", {{"d", 2112}})),
54  n_(distributionModelDict_.get<scalar>("n"))
55 {
56  if (lambda_ < VSMALL || n_ < VSMALL)
57  {
59  << "Scale/Shape parameter cannot be equal to or less than zero:"
60  << " lambda = " << lambda_
61  << " n = " << n_
62  << exit(FatalError);
63  }
64 
65  check();
66 }
67 
68 
70 (
71  const massRosinRammler& p
72 )
73 :
75  lambda_(p.lambda_),
76  n_(p.n_)
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
83 {
84  scalar d = 0;
85  do
86  {
87  // (YHD:Inverse of Eq. 10)
88  const scalar a = scalar(3)/n_ + scalar(1);
89  const scalar cdfA = Math::incGamma_P(a, pow(minValue_/lambda_, n_) );
90  const scalar cdfB = Math::incGamma_P(a, pow(maxValue_/lambda_, n_) );
91 
92  const scalar u = rndGen_.position<scalar>(cdfA, cdfB);
93  const scalar x = Math::invIncGamma(a, u);
94  d = lambda_*pow(x, scalar(1)/n_);
95 
96  } while (std::isnan(d));
97 
98  return d;
99 }
100 
101 
103 {
104  // (YHD:Eqs. 11-12)
105  return lambda_*tgamma(scalar(1)/n_ + scalar(1));
106 }
107 
108 
109 // ************************************************************************* //
Foam::distributionModel::rndGen_
Random & rndGen_
Reference to the random number generator.
Definition: distributionModel.H:82
Foam::Random
Random number generator.
Definition: Random.H:59
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
Foam::distributionModels::massRosinRammler::meanValue
virtual scalar meanValue() const
Return the theoretical mean of the distribution.
Definition: massRosinRammler.C:102
Foam::distributionModel
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Definition: distributionModel.H:72
Foam::distributionModels::massRosinRammler
Particle-size distribution model wherein random samples are drawn from the two-parameter Rosin-Rammle...
Definition: massRosinRammler.H:249
Foam::distributionModel::maxValue_
scalar maxValue_
Maximum of the distribution.
Definition: distributionModel.H:88
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
MathFunctions.H
Foam::distributionModel::minValue_
scalar minValue_
Minimum of the distribution.
Definition: distributionModel.H:85
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distributionModels::massRosinRammler::sample
virtual scalar sample() const
Sample the distribution.
Definition: massRosinRammler.C:82
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::Math::invIncGamma
scalar invIncGamma(const scalar a, const scalar P)
Inverse of regularised lower incomplete gamma function.
Definition: invIncGamma.C:91
Foam::Math::incGamma_P
scalar incGamma_P(const scalar a, const scalar x)
Lower incomplete gamma function.
Definition: incGamma.C:449
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Random::position
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Definition: RandomTemplates.C:62
Foam::distributionModels::massRosinRammler::massRosinRammler
massRosinRammler(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: massRosinRammler.C:47
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
massRosinRammler.H
x
x
Definition: LISASMDCalcMethod2.H:52
rndGen
Random rndGen
Definition: createFields.H:23