normal.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "normal.H"
28 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace distributionModels
35 {
36  defineTypeNameAndDebug(normal, 0);
37  addToRunTimeSelectionTable(distributionModel, normal, dictionary);
38 }
39 }
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  Random& rndGen
47 )
48 :
49  distributionModel(typeName, dict, rndGen),
50  minValue_(distributionModelDict_.template lookup<scalar>("minValue")),
51  maxValue_(distributionModelDict_.template lookup<scalar>("maxValue")),
52  expectation_(distributionModelDict_.template lookup<scalar>("expectation")),
53  variance_(distributionModelDict_.template lookup<scalar>("variance")),
54  a_(0.147)
55 {
56  check();
57  info();
58 }
59 
60 
62 :
64  minValue_(p.minValue_),
65  maxValue_(p.maxValue_),
66  expectation_(p.expectation_),
67  variance_(p.variance_),
68  a_(p.a_)
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82 
83  scalar a = erf((minValue_ - expectation_)/variance_);
84  scalar b = erf((maxValue_ - expectation_)/variance_);
85 
86  scalar y = rndGen_.sample01<scalar>();
87  scalar x = erfInv(y*(b - a) + a)*variance_ + expectation_;
88 
89  // Note: numerical approximation of the inverse function yields slight
90  // inaccuracies
91 
92  x = min(max(x, minValue_), maxValue_);
93 
94  return x;
95 }
96 
97 
99 {
100  return minValue_;
101 }
102 
103 
105 {
106  return maxValue_;
107 }
108 
109 
111 {
112  return expectation_;
113 }
114 
115 
116 Foam::scalar Foam::distributionModels::normal::erfInv(const scalar y) const
117 {
118  scalar k = 2.0/(constant::mathematical::pi*a_) + 0.5*log(1.0 - y*y);
119  scalar h = log(1.0 - y*y)/a_;
120  scalar x = sqrt(-k + sqrt(k*k - h));
121  if (y < 0.0)
122  {
123  x *= -1.0;
124  }
125  return x;
126 }
127 
128 
129 // ************************************************************************* //
A normal distribution model.
Definition: normal.H:56
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual scalar sample() const
Sample the distributionModel.
Definition: normal.C:80
dimensionedScalar log(const dimensionedScalar &ds)
Random & rndGen_
Reference to the random number generator.
virtual scalar minValue() const
Return the minimum value.
Definition: normal.C:98
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual scalar erfInv(const scalar y) const
Definition: normal.C:116
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
label k
Boltzmann constant.
const dimensionedScalar h
Planck constant.
Macros for easy insertion into run-time selection tables.
virtual scalar maxValue() const
Return the maximum value.
Definition: normal.C:104
scalar y
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:70
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
virtual scalar meanValue() const
Return the mean value.
Definition: normal.C:110
virtual ~normal()
Destructor.
Definition: normal.C:74
Random number generator.
Definition: Random.H:57
normal(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: normal.C:44
defineTypeNameAndDebug(exponential, 0)
dimensionedScalar erf(const dimensionedScalar &ds)
A library of runtime-selectable distribution models.
addToRunTimeSelectionTable(distributionModel, exponential, dictionary)
volScalarField & p
Namespace for OpenFOAM.