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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "normal.H"
30 #include "mathematicalConstants.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace distributionModels
37 {
38  defineTypeNameAndDebug(normal, 0);
39  addToRunTimeSelectionTable(distributionModel, normal, dictionary);
40 }
41 }
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
48  Random& rndGen
49 )
50 :
51  distributionModel(typeName, dict, rndGen),
52  minValue_(distributionModelDict_.get<scalar>("minValue")),
53  maxValue_(distributionModelDict_.get<scalar>("maxValue")),
54  expectation_(distributionModelDict_.get<scalar>("expectation")),
55  variance_(distributionModelDict_.get<scalar>("variance")),
56  a_(0.147)
57 {
58  if (minValue_ < 0)
59  {
61  << "Minimum value must be greater than zero. "
62  << "Supplied minValue = " << minValue_
63  << abort(FatalError);
64  }
65 
66  if (maxValue_ < minValue_)
67  {
69  << "Maximum value is smaller than the minimum value:"
70  << " maxValue = " << maxValue_ << ", minValue = " << minValue_
71  << abort(FatalError);
72  }
73 }
74 
75 
77 :
79  minValue_(p.minValue_),
80  maxValue_(p.maxValue_),
81  expectation_(p.expectation_),
82  variance_(p.variance_),
83  a_(p.a_)
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97 
98  scalar a = erf((minValue_ - expectation_)/variance_);
99  scalar b = erf((maxValue_ - expectation_)/variance_);
100 
101  scalar y = rndGen_.sample01<scalar>();
102  scalar x = erfInv(y*(b - a) + a)*variance_ + expectation_;
103 
104  // Note: numerical approximation of the inverse function yields slight
105  // inaccuracies
106 
107  x = min(max(x, minValue_), maxValue_);
108 
109  return x;
110 }
111 
112 
114 {
115  return minValue_;
116 }
117 
118 
120 {
121  return maxValue_;
122 }
123 
124 
126 {
127  return expectation_;
128 }
129 
130 
131 Foam::scalar Foam::distributionModels::normal::erfInv(const scalar y) const
132 {
133  scalar k = 2.0/(constant::mathematical::pi*a_) + 0.5*log(1.0 - y*y);
134  scalar h = log(1.0 - y*y)/a_;
135  scalar x = sqrt(-k + sqrt(k*k - h));
136  if (y < 0.0)
137  {
138  x *= -1.0;
139  }
140  return x;
141 }
142 
143 
144 // ************************************************************************* //
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:46
Foam::distributionModels::normal::erfInv
virtual scalar erfInv(const scalar y) const
Definition: normal.C:131
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
Foam::distributionModel
A library of runtime-selectable distribution models.
Definition: distributionModel.H:70
Foam::erf
dimensionedScalar erf(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:276
Foam::distributionModels::normal::~normal
virtual ~normal()
Destructor.
Definition: normal.C:89
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::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::distributionModels::normal
A normal distribution model.
Definition: normal.H:58
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
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 mean value.
Definition: normal.C:125
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:121
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::distributionModels::normal::sample
virtual scalar sample() const
Sample the distributionModel.
Definition: normal.C:95
Foam::distributionModels::normal::minValue
virtual scalar minValue() const
Return the minimum value.
Definition: normal.C:113
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
x
x
Definition: LISASMDCalcMethod2.H:52
rndGen
Random rndGen
Definition: createFields.H:23
Foam::distributionModels::normal::maxValue
virtual scalar maxValue() const
Return the maximum value.
Definition: normal.C:119
y
scalar y
Definition: LISASMDCalcMethod1.H:14