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-2023 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"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace distributions
35 {
38 }
39 }
40 
41 
42 const Foam::scalar Foam::distributions::normal::a_ = 0.147;
43 
44 
45 using namespace Foam::constant::mathematical;
46 
47 
48 // * * * * * * * * * * * Private Static Member Functions * * * * * * * * * * //
49 
50 Foam::tmp<Foam::scalarField> Foam::distributions::normal::approxErf
51 (
52  const scalarField& x
53 )
54 {
55  return sign(x)*sqrt(1 - exp(-sqr(x)*(4/pi + a_*sqr(x))/(1 + a_*sqr(x))));
56 }
57 
58 
59 Foam::scalar Foam::distributions::normal::approxErfInv(const scalar y)
60 {
61  const scalar l = log(1 - y*y), b = 2/(pi*a_) + l/2;
62  return sign(y)*sqrt(-b + sqrt(b*b - l/a_));
63 }
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
68 Foam::tmp<Foam::scalarField> Foam::distributions::normal::phi
69 (
70  const label q,
71  const scalarField& x
72 ) const
73 {
74  static const scalar sqrt2Pi = sqrt(2*pi);
75  return integerPow(x, q)/(sigma_*sqrt2Pi)*exp(- sqr((x - mu_)/sigma_)/2);
76 }
77 
78 
79 Foam::tmp<Foam::scalarField> Foam::distributions::normal::Phi
80 (
81  const label q,
82  const scalarField& x
83 ) const
84 {
85  if (q == 0)
86  {
87  static const scalar sqrt2 = sqrt(scalar(2));
88  return (1 + approxErf((x - mu_)/(sigma_*sqrt2)))/2;
89  }
90  else
91  {
92  return unintegrableForNonZeroQ::Phi(q, x);
93  }
94 }
95 
96 
97 Foam::scalar Foam::distributions::normal::sampleForZeroQ(const scalar s) const
98 {
99  static const scalar sqrt2 = sqrt(scalar(2));
100  const Pair<scalar>& Phi01 = this->Phi01();
101  const scalar PhiS = (1 - s)*Phi01[0] + s*Phi01[1];
102  return approxErfInv(2*PhiS - 1)*sigma_*sqrt2 + mu_;
103 }
104 
105 
106 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107 
109 (
110  const dictionary& dict,
111  Random& rndGen,
112  const label sampleQ
113 )
114 :
116  (
117  typeName,
118  dict,
119  rndGen,
120  sampleQ
121  ),
122  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"})),
123  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"})),
124  mu_(dict.lookupBackwardsCompatible<scalar>({"mu", "expectation"})),
125  sigma_(dict.lookup<scalar>("sigma"))
126 {
127  validateBounds(dict);
128  if (q() != 0) validatePositive(dict);
129  mean();
130  report();
131 }
132 
133 
135 (
136  Random& rndGen,
137  const label Q,
138  const label sampleQ,
139  const label n,
140  const scalar min,
141  const scalar max,
142  const scalar mu,
143  const scalar sigma
144 )
145 :
147  min_(min),
148  max_(max),
149  mu_(mu),
150  sigma_(sigma)
151 {}
152 
153 
155 :
157  min_(d.min_),
158  max_(d.max_),
159  mu_(d.mu_),
160  sigma_(d.sigma_)
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
165 
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
173 {
174  if (q() == 0)
175  {
176  return sampleForZeroQ(rndGen_.sample01<scalar>());
177  }
178  else
179  {
181  }
182 }
183 
184 
186 {
187  return min_;
188 }
189 
190 
192 {
193  return max_;
194 }
195 
196 
198 {
199  if (q() == 0)
200  {
201  // !!! This isn't technically correct. The min/max clipping affects the
202  // mean. To calculate it properly we have to integrate the first moment
203  // of the PDF, and that can't be done analytically; we have to hand
204  // over to the unintegrable handling. This probably isn't worth it for
205  // what should be a small deviation. Also, it looks a bit odd in the
206  // log for the distribution to be reporting a different mean to that
207  // which was specified. So, just use the mean value of the un-clipped
208  // distribution for now.
209  return mu_;
210  }
211  else
212  {
214  }
215 }
216 
217 
220 {
222  tx.ref()[0] = Foam::max(tx.ref()[0], q() < 0 ? min_/2 : -vGreat);
223  return tx;
224 }
225 
226 
227 // ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
Random number generator.
Definition: Random.H:58
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:140
Normal distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: normal.H:74
virtual scalar min() const
Return the minimum value.
Definition: normal.C:185
virtual scalar sample() const
Sample the distribution.
Definition: normal.C:172
normal(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: normal.C:109
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: normal.C:219
virtual ~normal()
Destructor.
Definition: normal.C:166
virtual scalar max() const
Return the maximum value.
Definition: normal.C:191
virtual scalar mean() const
Return the mean value.
Definition: normal.C:197
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:349
virtual scalar mean() const
Return the mean value.
Definition: unintegrable.C:366
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
mathematical constants.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar integerPow(const scalar x, const label e)
Compute the power of the number x to the integer e.
Definition: scalarI.H:30
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
Random rndGen(label(0))