normal.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) 2011-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 "normal.H"
30 #include "mathematicalConstants.H"
31 #include "MathFunctions.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace distributionModels
39 {
40  defineTypeNameAndDebug(normal, 0);
41  addToRunTimeSelectionTable(distributionModel, normal, dictionary);
42 }
43 }
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const dictionary& dict,
50  Random& rndGen
51 )
52 :
53  distributionModel(typeName, dict, rndGen),
54  mu_
55  (
56  distributionModelDict_.getCompat<scalar>
57  (
58  "mu",
59  {{"expectation", 2112}}
60  )
61  ),
62  sigma_
63  (
64  distributionModelDict_.getCompat<scalar>
65  (
66  "sigma",
67  {{"variance", 2112}}
68  )
69  )
70 {
71  if (mag(sigma_) == 0)
72  {
74  << "Standard deviation cannot be zero." << nl
75  << " sigma = " << sigma_ << nl
76  << exit(FatalError);
77  }
78 
79  check();
80 }
81 
82 
84 :
86  mu_(p.mu_),
87  sigma_(p.sigma_)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
94 {
95  const scalar a = (minValue_ - mu_)/sigma_;
96  const scalar b = (maxValue_ - mu_)/sigma_;
97 
98  const scalar aPhi = 0.5*(scalar(1) + erf(a/Foam::sqrt(scalar(2))));
99  const scalar bPhi = 0.5*(scalar(1) + erf(b/Foam::sqrt(scalar(2))));
100 
101  const scalar u = rndGen_.sample01<scalar>();
102  const scalar p = u*(bPhi - aPhi) + aPhi;
103 
104  // (B:p. 20-24)
105  const scalar x =
106  mu_ + sigma_*Foam::sqrt(scalar(2))*Math::erfInv(scalar(2)*p - scalar(1));
107 
108  // Note: numerical approximation of the inverse function yields slight
109  // inaccuracies
110 
111  return min(max(x, minValue_), maxValue_);
112 }
113 
114 
116 {
117  const scalar a = (minValue_ - mu_)/sigma_;
118  const scalar b = (maxValue_ - mu_)/sigma_;
119 
120  // (B:p. 2)
121  const scalar aphi =
122  scalar(1)/Foam::sqrt(scalar(2)*constant::mathematical::pi)
123  *exp(-0.5*sqr(a));
124  const scalar bphi =
125  scalar(1)/Foam::sqrt(scalar(2)*constant::mathematical::pi)
126  *exp(-0.5*sqr(b));
127 
128  // (B:p. 4)
129  const scalar aPhi = 0.5*(scalar(1) + erf(a/Foam::sqrt(scalar(2))));
130  const scalar bPhi = 0.5*(scalar(1) + erf(b/Foam::sqrt(scalar(2))));
131 
132  // (B:p. 25)
133  return mu_ - sigma_*(bphi - aphi)/(bPhi - aPhi + VSMALL);
134 }
135 
136 
137 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
mathematicalConstants.H
normal.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::distributionModels::normal::normal
normal(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: normal.C:48
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::distributionModel
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Definition: distributionModel.H:72
Foam::erf
dimensionedScalar erf(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:276
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
MathFunctions.H
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::distributionModels::normal
Particle-size distribution model wherein random samples are drawn from the doubly-truncated univariat...
Definition: normal.H:248
Foam::Math::erfInv
scalar erfInv(const scalar y)
Inverse error function of a real-number argument.
Definition: erfInv.C:34
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::distributionModels::normal::meanValue
virtual scalar meanValue() const
Return the theoretical mean of the distribution.
Definition: normal.C:115
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::normal::sample
virtual scalar sample() const
Sample the distribution.
Definition: normal.C:93
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
x
x
Definition: LISASMDCalcMethod2.H:52
rndGen
Random rndGen
Definition: createFields.H:23