Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield Namespace Reference

Namespaces

namespace  HeedInterface
 
namespace  Magboltz
 
namespace  Numerics
 

Classes

class  AvalancheMC
 
class  AvalancheMicroscopic
 
class  ComponentAnalyticField
 
class  ComponentAnsys121
 
class  ComponentAnsys123
 
class  ComponentBase
 
class  ComponentConstant
 
class  ComponentCST
 
class  ComponentElmer
 
class  ComponentFieldMap
 
class  ComponentNeBem2d
 
class  ComponentTcad2d
 
class  ComponentTcad3d
 
class  ComponentUser
 
class  ComponentUserMapBase
 
class  ComponentVoxel
 
class  DriftLineRKF
 
class  GeometryBase
 
class  GeometryRoot
 
class  GeometrySimple
 
class  HeedChamber
 
class  Medium
 
class  MediumCdTe
 
class  MediumConductor
 
class  MediumGaAs
 
class  MediumGas
 
class  MediumMagboltz
 
class  MediumPlastic
 
class  MediumSilicon
 
class  OpticalData
 
class  PlottingEngine
 
class  PlottingEngineRoot
 
struct  PolygonInfo
 
class  RandomEngine
 
class  RandomEngineRoot
 
class  Sensor
 
class  Solid
 
class  SolidBox
 
class  SolidSphere
 
class  SolidTube
 
class  Track
 
class  TrackBichsel
 
class  TrackElectron
 
class  TrackHeed
 
class  TrackPAI
 
class  TrackSimple
 
class  ViewCell
 
class  ViewDrift
 
class  ViewFEMesh
 
class  ViewField
 
class  ViewGeometry
 
class  ViewMedium
 
class  ViewSignal
 

Functions

double RndmUniform ()
 
double RndmUniformPos ()
 
double RndmGaussian ()
 
double RndmGaussian (const double mu, const double sigma)
 
double RndmLorentzian (const double mu, const double gamma)
 
double RndmVoigt (const double mu, const double sigma, const double gamma)
 
double RndmPolya (const double theta)
 

Variables

PlottingEngineRoot plottingEngine
 
RandomEngineRoot randomEngine
 

Function Documentation

◆ RndmGaussian() [1/2]

double Garfield::RndmGaussian ( )
inline

Definition at line 27 of file Random.hh.

27 {
28
29 static bool cached = false;
30 static double u = 0.;
31 if (cached) {
32 cached = false;
33 return u;
34 }
35 // Box-Muller algorithm
36 u = 2. * RndmUniform() - 1.;
37 double v = 2. * RndmUniform() - 1.;
38 double r2 = u * u + v * v;
39 while (r2 > 1.) {
40 u = 2. * RndmUniform() - 1.;
41 v = 2. * RndmUniform() - 1.;
42 r2 = u * u + v * v;
43 }
44 const double p = sqrt(-2. * log(r2) / r2);
45 u *= p;
46 cached = true;
47 return v * p;
48}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313

Referenced by RndmGaussian(), and RndmVoigt().

◆ RndmGaussian() [2/2]

double Garfield::RndmGaussian ( const double  mu,
const double  sigma 
)
inline

Definition at line 51 of file Random.hh.

51 {
52
53 return mu + sigma * RndmGaussian();
54}
double RndmGaussian()
Definition: Random.hh:27

◆ RndmLorentzian()

double Garfield::RndmLorentzian ( const double  mu,
const double  gamma 
)
inline

Definition at line 58 of file Random.hh.

58 {
59
60 return mu + gamma * tan(Pi * (RndmUniform() - 0.5));
61}
double RndmUniform()
Definition: Random.hh:16

Referenced by RndmVoigt().

◆ RndmPolya()

double Garfield::RndmPolya ( const double  theta)
inline

Definition at line 77 of file Random.hh.

77 {
78
79 // Algorithm from Review of Particle Physics
80 // C. Amsler et al, Phys. Lett. B 667 (2008)
81 if (theta <= 0.) return -log(RndmUniformPos());
82 const double c = 3 * (theta + 1.) - 0.75;
83 double u1, u2, v1, v2, v3;
84 double x;
85 while (1) {
86 u1 = RndmUniformPos();
87 v1 = u1 * (1. - u1);
88 if (v1 == 0.) continue;
89 v2 = (u1 - 0.5) * sqrt(c / v1);
90 x = theta + v2;
91 if (x <= 0.) continue;
92 u2 = RndmUniformPos();
93 v3 = 64 * u2 * u2 * pow(v1, 3);
94 if (v3 <= 1. - 2 * v2 * v2 / x ||
95 log(v3) <= 2 * (theta * log(x / theta) - v2)) {
96 return x / (theta + 1.);
97 }
98 }
99}
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:336
double RndmUniformPos()
Definition: Random.hh:19

◆ RndmUniform()

◆ RndmUniformPos()

◆ RndmVoigt()

double Garfield::RndmVoigt ( const double  mu,
const double  sigma,
const double  gamma 
)
inline

Definition at line 67 of file Random.hh.

68 {
69
70 if (sigma <= 0.) return RndmLorentzian(mu, gamma);
71 const double a = gamma / (Sqrt2 * sigma);
72 const double x = RndmLorentzian(0., a) + RndmGaussian(0., 1. / Sqrt2);
73 return mu + x * Sqrt2 * sigma;
74}
double RndmLorentzian(const double mu, const double gamma)
Definition: Random.hh:58

Variable Documentation

◆ plottingEngine

◆ randomEngine

RandomEngineRoot Garfield::randomEngine

Definition at line 6 of file RandomEngineRoot.cc.

Referenced by main(), and RndmUniform().