multiNormal.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-------------------------------------------------------------------------------
11License
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 "multiNormal.H"
31#include "MathFunctions.H"
32#include "ListOps.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace distributionModels
40{
43}
44}
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const dictionary& dict,
52)
53:
54 distributionModel(typeName, dict, rndGen),
55 mu_
56 (
57 distributionModelDict_.lookupCompat
58 (
59 "mu",
60 {{"expectation", 2112}}
61 )
62 ),
63 sigma_
64 (
65 distributionModelDict_.lookupCompat
66 (
67 "sigma",
68 {{"variance", 2112}}
69 )
70 ),
71 weight_
72 (
73 distributionModelDict_.lookupCompat
74 (
75 "weight",
76 {{"strength", 2112}}
77 )
78 )
79{
80 check();
81
82 scalar sum = 0;
83 for (label i = 0; i < weight_.size(); ++i)
84 {
85 if (i > 0 && weight_[i] < weight_[i-1])
86 {
88 << type() << "distribution: "
89 << "Weights must be specified in a monotonic order." << nl
90 << "Please see the row i = " << i << nl
91 << "weight[i-1] = " << weight_[i-1] << nl
92 << "weight[i] = " << weight_[i]
93 << exit(FatalError);
94 }
95
96 sum += weight_[i];
97 }
98
99 if (sum < VSMALL)
100 {
102 << type() << "distribution: "
103 << "The sum of weights cannot be zero." << nl
104 << "weight = " << weight_
105 << exit(FatalError);
106 }
107
108 for (label i = 1; i < weight_.size(); ++i)
109 {
110 weight_[i] += weight_[i-1];
111 }
112
113 for (auto& w : weight_)
114 {
115 w /= sum;
116 }
117}
118
119
121:
123 mu_(p.mu_),
124 sigma_(p.sigma_),
125 weight_(p.weight_)
126{}
127
128
129// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130
132{
133 const scalar u = rndGen_.sample01<scalar>();
134
135 for (label i = 0; i < weight_.size(); ++i)
136 {
137 if (weight_[i] > u)
138 {
139 return sample(mu_[i], sigma_[i]);
140 }
141 }
142
143 const label last = weight_.size() - 1;
144
145 return sample(mu_[last], sigma_[last]);
146}
147
148
150(
151 const scalar mu,
152 const scalar sigma
153) const
154{
155 const scalar a = (minValue_ - mu)/sigma;
156 const scalar b = (maxValue_ - mu)/sigma;
157
158 const scalar aPhi = 0.5*(1.0 + erf(a/Foam::sqrt(2.0)));
159 const scalar bPhi = 0.5*(1.0 + erf(b/Foam::sqrt(2.0)));
160
161 const scalar u = rndGen_.sample01<scalar>();
162 const scalar p = u*(bPhi - aPhi) + aPhi;
163
164 // (B:p. 20-24)
165 const scalar x =
166 mu + sigma*Foam::sqrt(2.0)*Math::erfInv(2.0*p - 1.0);
167
168 // Note: numerical approximation of the inverse function yields slight
169 // inaccuracies
170
171 return min(max(x, minValue_), maxValue_);
172}
173
174
176{
177 scalar mean = 0;
178 forAll(weight_, i)
179 {
180 mean += weight_[i]*mu_[i];
181 }
182
183 return mean;
184}
185
186
187// ************************************************************************* //
Various functions to operate on Lists.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Minimal example by using system/controlDict.functions:
Random number generator.
Definition: Random.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Particle-size distribution model wherein random samples are drawn from a mixture of a finite set of d...
Definition: multiNormal.H:288
virtual scalar meanValue() const
Return the theoretical mean of the distribution.
Definition: multiNormal.C:175
virtual scalar sample() const
Sample the distribution.
Definition: multiNormal.C:131
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & mu
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
scalar erfInv(const scalar y)
Inverse error function of a real-number argument.
Definition: erfInv.C:36
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar erf(const dimensionedScalar &ds)
static void check(const int retVal, const char *what)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23