standardNormal.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) 2024 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 "standardNormal.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace distributions
35 {
37 }
38 }
39 
40 
41 const Foam::scalar Foam::distributions::standardNormal::a_ = 0.147;
42 
43 
44 using namespace Foam::constant::mathematical;
45 
46 
47 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
48 
50 (
51  const scalarField& x
52 )
53 {
54  return sign(x)*sqrt(1 - exp(-sqr(x)*(4/pi + a_*sqr(x))/(1 + a_*sqr(x))));
55 }
56 
57 
59 {
60  const scalar l = log(Foam::max(1 - y*y, small/2)), b = 2/(pi*a_) + l/2;
61  return sign(y)*sqrt(-b + sqrt(b*b - l/a_));
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
68 :
70 {}
71 
72 
74 (
75  const randomGenerator::seed& s,
76  const bool global
77 )
78 :
80 {}
81 
82 
84 :
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  static const scalar sqrt2 = sqrt(scalar(2));
100  const scalar s = rndGen_.sample01<scalar>();
101  return approxErfInv(2*s - 1)*sqrt2;
102 }
103 
104 
106 {
107  return -vGreat;
108 }
109 
110 
112 {
113  return vGreat;
114 }
115 
116 
118 {
119  return 0;
120 }
121 
122 
125 (
126  const scalarField& x,
127  const label e,
128  const bool
129 ) const
130 {
131  if (e != 0) NotImplemented;
132 
133  static const scalar sqrt2 = sqrt(scalar(2));
134  return (1 + approxErf(x/sqrt2))/2;
135 }
136 
137 
140 {
141  const scalar x0 = approxErfInv(1 - rootSmall);
142  const scalar x1 = approxErfInv(rootSmall - 1);
143 
144  tmp<scalarField> tResult(new scalarField(n));
145  scalarField& result = tResult.ref();
146 
147  forAll(result, i)
148  {
149  const scalar f = scalar(i)/(n - 1);
150  result[i] = (1 - f)*x0 + f*x1;
151  }
152 
153  return tResult;
154 }
155 
156 
158 (
159  const scalarField& x
160 ) const
161 {
162  static const scalar sqrt2Pi = sqrt(2*pi);
163  return exp(- sqr(x)/2)/sqrt2Pi;
164 }
165 
166 
167 // ************************************************************************* //
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Base class for statistical distributions.
Definition: distribution.H:76
Standard normal distribution. Not selectable.
static tmp< scalarField > approxErf(const scalarField &x)
Approximate error function.
virtual tmp< scalarField > integralPDFxPow(const scalarField &x, const label e, const bool consistent=false) const
Return the integral of the PDF multiplied by an integer power of x.
virtual scalar min() const
Return the minimum value.
virtual ~standardNormal()
Destructor.
virtual scalar sample() const
Sample the distribution.
static scalar approxErfInv(const scalar y)
Approximate error function inverse.
virtual tmp< scalarField > plotPDF(const scalarField &x) const
Return values to plot the probability density function.
standardNormal(randomGenerator &&rndGen)
Construct from a random generator.
virtual scalar max() const
Return the maximum value.
virtual scalar mean() const
Return the mean value.
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
Random number generator.
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:197
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
mathematical constants.
defineTypeNameAndDebug(exponential, 0)
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
labelList f(nPoints)
randomGenerator rndGen(653213)