Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Random.hh
Go to the documentation of this file.
1#ifndef G_RANDOM_H
2#define G_RANDOM_H
3
4#include <cmath>
5#include <functional>
6
8
9namespace Garfield {
10
11class Random {
12 public:
13 Random() = default;
14 template <typename T>
15 static inline void SetEngine(T engine) {
16 draw = std::bind(&T::Draw, engine);
17 }
18 inline double static Draw() noexcept { return draw(); };
19
20 private:
21 static std::function<double(void)> draw;
22};
23
25inline double RndmUniform() { return Random::Draw(); }
26
28inline double RndmUniformPos() {
29 double r = RndmUniform();
30 while (r <= 0.) r = RndmUniform();
31 return r;
32}
33
35inline double RndmGaussian() {
36 static bool cached = false;
37 static double u = 0.;
38 if (cached) {
39 cached = false;
40 return u;
41 }
42 // Box-Muller algorithm
43 u = 2. * RndmUniform() - 1.;
44 double v = 2. * RndmUniform() - 1.;
45 double r2 = u * u + v * v;
46 while (r2 > 1.) {
47 u = 2. * RndmUniform() - 1.;
48 v = 2. * RndmUniform() - 1.;
49 r2 = u * u + v * v;
50 }
51 const double p = sqrt(-2. * log(r2) / r2);
52 u *= p;
53 cached = true;
54 return v * p;
55}
56
58inline double RndmGaussian(const double mu, const double sigma) {
59 return mu + sigma * RndmGaussian();
60}
61
64inline double RndmLorentzian(const double mu, const double gamma) {
65 return mu + gamma * tan(Pi * (RndmUniform() - 0.5));
66}
67
72inline double RndmVoigt(const double mu, const double sigma,
73 const double gamma) {
74 if (sigma <= 0.) return RndmLorentzian(mu, gamma);
75 const double a = gamma / (Sqrt2 * sigma);
76 const double x = RndmLorentzian(0., a) + RndmGaussian(0., 1. / Sqrt2);
77 return mu + x * Sqrt2 * sigma;
78}
79
81inline unsigned int RndmYuleFurry(const double mean) {
82 if (mean <= 0.) return 0;
83 return 1 + static_cast<unsigned int>(std::log(RndmUniformPos()) /
84 std::log1p(-1. / mean));
85 /*
86 const double u = RndmUniform();
87 double p = 1. / mean;
88 const double q = 1. - p;
89 double sum = p;
90 unsigned int k = 1;
91 while (sum < u) {
92 ++k;
93 p *= q;
94 sum += p;
95 }
96 return k;
97 */
98}
99
101inline double RndmPolya(const double theta) {
102 // Algorithm from Review of Particle Physics
103 // C. Amsler et al, Phys. Lett. B 667 (2008)
104 if (theta <= 0.) return -log(RndmUniformPos());
105 const double c = 3 * (theta + 1.) - 0.75;
106 double u1, u2, v1, v2, v3;
107 double x;
108 while (1) {
109 u1 = RndmUniformPos();
110 v1 = u1 * (1. - u1);
111 if (v1 == 0.) continue;
112 v2 = (u1 - 0.5) * sqrt(c / v1);
113 x = theta + v2;
114 if (x <= 0.) continue;
115 u2 = RndmUniformPos();
116 v3 = 64 * u2 * u2 * pow(v1, 3);
117 if (v3 <= 1. - 2 * v2 * v2 / x ||
118 log(v3) <= 2 * (theta * log(x / theta) - v2)) {
119 return x / (theta + 1.);
120 }
121 }
122}
123
125double RndmLandau();
127double RndmVavilov(const double rkappa, const double beta2);
128
130int RndmPoisson(const double mean);
131
135double RndmHeedWF(const double w, const double f);
136
138inline void RndmDirection(double& dx, double& dy, double& dz,
139 const double length = 1.) {
140 const double phi = TwoPi * RndmUniform();
141 const double ctheta = 2 * RndmUniform() - 1.;
142 const double stheta = sqrt(1. - ctheta * ctheta);
143 dx = length * cos(phi) * stheta;
144 dy = length * sin(phi) * stheta;
145 dz = length * ctheta;
146}
147} // namespace Garfield
148
149#endif
Random()=default
static std::function< double(void)> draw
Definition Random.hh:21
static void SetEngine(T engine)
Definition Random.hh:15
static double Draw() noexcept
Definition Random.hh:18
double RndmLorentzian(const double mu, const double gamma)
Draw a Lorentzian random variate with mean mu and half-width at half maximum gamma.
Definition Random.hh:64
unsigned int RndmYuleFurry(const double mean)
Draw a random number from a geometric distribution.
Definition Random.hh:81
double RndmHeedWF(const double w, const double f)
Draw a random energy needed to create a single electron in a material asymptotic work function W and ...
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition Random.hh:25
double RndmVavilov(const double rkappa, const double beta2)
Draw a random number from a Vavilov distribution.
double RndmPolya(const double theta)
Draw a Polya distributed random number.
Definition Random.hh:101
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition Random.hh:138
double RndmVoigt(const double mu, const double sigma, const double gamma)
Draw a random number according to a Voigt function with mean mu.
Definition Random.hh:72
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition Random.hh:28
double RndmGaussian()
Draw a Gaussian random variate with mean zero and standard deviation one.
Definition Random.hh:35
double RndmLandau()
Draw a random number from a Landau distribution.
int RndmPoisson(const double mean)
Draw a random number from a Poisson distribution.