Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
AvalancheMC.hh
Go to the documentation of this file.
1#ifndef G_AVALANCHE_MC_H
2#define G_AVALANCHE_MC_H
3
4#include <vector>
5#include <string>
6
7#include "Sensor.hh"
8#include "ViewDrift.hh"
9
10namespace Garfield {
11
13
14 public:
15 // Constructor
17 // Destructor
19
20 void SetSensor(Sensor* s);
21
22 // Switch on/off drift line plotting
23 void EnablePlotting(ViewDrift* view);
24 void DisablePlotting();
25
26 // Switch on/off calculation of induced currents (default: disabled)
27 void EnableSignalCalculation() { m_useSignal = true; }
28 void DisableSignalCalculation() { m_useSignal = false; }
29
30 // Switch on/off calculation of induced charge (default: disabled)
31 void EnableInducedChargeCalculation() { m_useInducedCharge = true; }
32 void DisableInducedChargeCalculation() { m_useInducedCharge = false; }
33
34 // Switch on/off equilibration of multiplication and attachment
35 // over the drift line (default: enabled)
36 void EnableProjectedPathIntegration() { m_useEquilibration = true; }
37 void DisableProjectedPathIntegration() { m_useEquilibration = false; }
38
39 // Switch on/off diffusion (default: enabled)
40 void EnableDiffusion() { m_useDiffusion = true; }
41 void DisableDiffusion() { m_useDiffusion = false; }
42
43 // Switch on/off attachment (and multiplication) for
44 // drift line calculation (default: enabled)
45 // For avalanches the flag is ignored
46 void EnableAttachment() { m_useAttachment = true; }
47 void DisableAttachment() { m_useAttachment = false; }
48
49 // Enable/disable magnetic field in stepping algorithm.
50 void EnableMagneticField() { m_useBfield = true; }
51 void DisableMagneticField() { m_useBfield = false; }
52
53 // Stepping model
54 // Fixed time step (default 20 ps)
55 void SetTimeSteps(const double d = 0.02);
56 // Fixed distance step (default 10 um)
57 void SetDistanceSteps(const double d = 0.001);
58 // Exponentially distributed time step with mean equal
59 // to the specified multiple of the collision time (default model)
60 void SetCollisionSteps(const int n = 100);
61
62 void SetTimeWindow(const double t0, const double t1);
63 void UnsetTimeWindow();
64
65 // Treat positive charge carriers as holes or ions (default: ions)
66 void SetHoles() { m_useIons = false; }
67 void SetIons() { m_useIons = true; }
68
69 void SetElectronSignalScalingFactor(const double scale) {
70 m_scaleElectronSignal = scale;
71 }
72 void SetHoleSignalScalingFactor(const double scale) {
73 m_scaleHoleSignal = scale;
74 }
75 void SetIonSignalScalingFactor(const double scale) {
76 m_scaleIonSignal = scale;
77 }
78
79 void GetAvalancheSize(int& ne, int& ni) const {
80 ne = m_nElectrons;
81 ni = m_nIons;
82 }
83
84 unsigned int GetNumberOfDriftLinePoints() const { return m_nDrift; }
85 void GetDriftLinePoint(const unsigned int i, double& x, double& y, double& z,
86 double& t);
87
88 unsigned int GetNumberOfElectronEndpoints() const {
89 return m_nEndpointsElectrons;
90 }
91 unsigned int GetNumberOfHoleEndpoints() const {
92 return m_nEndpointsHoles;
93 }
94 unsigned int GetNumberOfIonEndpoints() const {
95 return m_nEndpointsIons;
96 }
97
98 void GetElectronEndpoint(const unsigned int i,
99 double& x0, double& y0, double& z0,
100 double& t0, double& x1, double& y1, double& z1,
101 double& t1, int& status) const;
102 void GetHoleEndpoint(const unsigned int i,
103 double& x0, double& y0, double& z0,
104 double& t0, double& x1, double& y1, double& z1,
105 double& t1, int& status) const;
106 void GetIonEndpoint(const unsigned int i,
107 double& x0, double& y0, double& z0,
108 double& t0, double& x1, double& y1, double& z1,
109 double& t1, int& status) const;
110
111 bool DriftElectron(const double x0, const double y0, const double z0,
112 const double t0);
113 bool DriftHole(const double x0, const double y0, const double z0,
114 const double t0);
115 bool DriftIon(const double x0, const double y0, const double z0,
116 const double t0);
117 bool AvalancheElectron(const double x0, const double y0, const double z0,
118 const double t0, const bool hole = false);
119 bool AvalancheHole(const double x0, const double y0, const double z0,
120 const double t0, const bool electron = false);
121 bool AvalancheElectronHole(const double x0, const double y0, const double z0,
122 const double t0);
123
124 // Switch on/off debugging messages
125 void EnableDebugging() { m_debug = true; }
126 void DisableDebugging() { m_debug = false; }
127
128 private:
129 std::string m_className;
130
131 // Numerical prefactors
132 static double c1;
133
134 Sensor* m_sensor;
135
136 unsigned int m_nDrift;
137 struct driftPoint {
138 // Position
139 double x, y, z, t;
140 // Townsend and attachment coefficient
141 double alpha, eta;
142 // Number of secondaries produced at this point
143 int ne, nh, ni;
144 };
145 std::vector<driftPoint> m_drift;
146
147 struct avalPoint {
148 double x, y, z, t;
149 int ne, nh, ni;
150 };
151 std::vector<avalPoint> m_aval;
152
153 // Step size model
154 int m_stepModel;
155 // Fixed time step
156 double m_tMc;
157 // Fixed distance step
158 double m_dMc;
159 // Sample step size according to collision time
160 int m_nMc;
161
162 // Time window
163 bool m_hasTimeWindow;
164 double m_tMin, m_tMax;
165
166 // Number of electrons, holes and ions produced
167 unsigned int m_nElectrons;
168 unsigned int m_nHoles;
169 unsigned int m_nIons;
170
171 // Number of endpoints (including captured electrons)
172 unsigned int m_nEndpointsElectrons;
173 unsigned int m_nEndpointsHoles;
174 unsigned int m_nEndpointsIons;
175 struct endpoint {
176 double x0, y0, z0, t0;
177 double x1, y1, z1, t1;
178 int status;
179 };
180 std::vector<endpoint> m_endpointsElectrons;
181 std::vector<endpoint> m_endpointsHoles;
182 std::vector<endpoint> m_endpointsIons;
183
184 bool m_usePlotting;
185 ViewDrift* m_viewer;
186
187 bool m_useSignal;
188 bool m_useInducedCharge;
189 bool m_useEquilibration;
190 bool m_useDiffusion;
191 bool m_useAttachment;
192 bool m_useBfield;
193 bool m_useIons;
194 bool m_withElectrons;
195 bool m_withHoles;
196 double m_scaleElectronSignal;
197 double m_scaleHoleSignal;
198 double m_scaleIonSignal;
199
200 bool m_debug;
201
202 // Compute a drift line with starting point (x0, y0, z0)
203 bool DriftLine(const double x0, const double y0, const double z0,
204 const double t0, const int type, const bool aval = false);
205 bool Avalanche();
206 // Compute effective multiplication and ionisation
207 // for the current drift line
208 bool ComputeAlphaEta(const int q);
209 // Compute the induced signal for the current drift line
210 void ComputeSignal(const double q);
211 void ComputeInducedCharge(const double q);
212
213};
214}
215
216#endif
void EnableInducedChargeCalculation()
Definition: AvalancheMC.hh:31
void EnablePlotting(ViewDrift *view)
Definition: AvalancheMC.cc:63
void SetTimeSteps(const double d=0.02)
Definition: AvalancheMC.cc:81
bool AvalancheElectronHole(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:910
bool DriftIon(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:279
unsigned int GetNumberOfIonEndpoints() const
Definition: AvalancheMC.hh:94
void SetIonSignalScalingFactor(const double scale)
Definition: AvalancheMC.hh:75
void DisableProjectedPathIntegration()
Definition: AvalancheMC.hh:37
void SetTimeWindow(const double t0, const double t1)
Definition: AvalancheMC.cc:133
void GetIonEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
Definition: AvalancheMC.cc:185
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:229
void GetElectronEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
Definition: AvalancheMC.cc:206
void SetCollisionSteps(const int n=100)
Definition: AvalancheMC.cc:115
void DisableSignalCalculation()
Definition: AvalancheMC.hh:28
void GetHoleEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
Definition: AvalancheMC.cc:163
void GetDriftLinePoint(const unsigned int i, double &x, double &y, double &z, double &t)
Definition: AvalancheMC.cc:148
void SetSensor(Sensor *s)
Definition: AvalancheMC.cc:52
void GetAvalancheSize(int &ne, int &ni) const
Definition: AvalancheMC.hh:79
unsigned int GetNumberOfHoleEndpoints() const
Definition: AvalancheMC.hh:91
void SetElectronSignalScalingFactor(const double scale)
Definition: AvalancheMC.hh:69
void SetHoleSignalScalingFactor(const double scale)
Definition: AvalancheMC.hh:72
unsigned int GetNumberOfDriftLinePoints() const
Definition: AvalancheMC.hh:84
unsigned int GetNumberOfElectronEndpoints() const
Definition: AvalancheMC.hh:88
bool AvalancheElectron(const double x0, const double y0, const double z0, const double t0, const bool hole=false)
Definition: AvalancheMC.cc:848
void DisableInducedChargeCalculation()
Definition: AvalancheMC.hh:32
void EnableProjectedPathIntegration()
Definition: AvalancheMC.hh:36
void SetDistanceSteps(const double d=0.001)
Definition: AvalancheMC.cc:98
void EnableSignalCalculation()
Definition: AvalancheMC.hh:27
bool AvalancheHole(const double x0, const double y0, const double z0, const double t0, const bool electron=false)
Definition: AvalancheMC.cc:879
bool DriftHole(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:254