Garfield++ v2r0
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
12/// Calculate drift lines and avalanches based on macroscopic transport
13/// coefficients, using Monte Carlo integration.
14
16
17 public:
18 /// Constructor
20 /// Destructor
22
23 /// Set the sensor.
24 void SetSensor(Sensor* s);
25
26 /// Switch on drift line plotting.
27 void EnablePlotting(ViewDrift* view);
28 void DisablePlotting() { m_viewer = NULL; }
29
30 /// Switch on calculation of induced currents (default: disabled).
31 void EnableSignalCalculation() { m_useSignal = true; }
32 void DisableSignalCalculation() { m_useSignal = false; }
33
34 /// Switch on calculation of induced charge (default: disabled).
35 void EnableInducedChargeCalculation() { m_useInducedCharge = true; }
36 void DisableInducedChargeCalculation() { m_useInducedCharge = false; }
37
38 /** Switch on equilibration of multiplication and attachment
39 * over the drift line (default: enabled) */
40 void EnableProjectedPathIntegration() { m_useEquilibration = true; }
41 void DisableProjectedPathIntegration() { m_useEquilibration = false; }
42
43 /// Switch on diffusion (default: enabled)
44 void EnableDiffusion() { m_useDiffusion = true; }
45 void DisableDiffusion() { m_useDiffusion = false; }
46
47 /** Switch on attachment (and multiplication) for drift line calculation
48 * (default: enabled). For avalanches the flag is ignored. */
49 void EnableAttachment() { m_useAttachment = true; }
50 void DisableAttachment() { m_useAttachment = false; }
51
52 /// Switch on calculating trapping with TCAD traps.
53 void EnableTcadTraps() { m_useTcadTrapping = true; }
54 void DisableTcadTraps() { m_useTcadTrapping = false; }
55
56 /// Switch on TCAD velocity maps
57 void EnableTcadVelocity() { m_useTcadVelocity = true; }
58 void DisableTcadVelocity() { m_useTcadVelocity = false; }
59
60 /// Enable use of magnetic field in stepping algorithm.
61 void EnableMagneticField() { m_useBfield = true; }
62 void DisableMagneticField() { m_useBfield = false; }
63
64 /// Use fixed-time steps (default 20 ps).
65 void SetTimeSteps(const double d = 0.02);
66 /// Use fixed distance steps (default 10 um).
67 void SetDistanceSteps(const double d = 0.001);
68 /** Use exponentially distributed time steps with mean equal
69 * to the specified multiple of the collision time (default model).*/
70 void SetCollisionSteps(const int n = 100);
71
72 /// Define a time interval (only carriers inside the interval are drifted).
73 void SetTimeWindow(const double t0, const double t1);
74 void UnsetTimeWindow() { m_hasTimeWindow = false; }
75
76 /// Treat positive charge carriers as holes (default: ions).
77 void SetHoles() { m_useIons = false; }
78 void SetIons() { m_useIons = true; }
79
80 /// Set multiplication factor for the signal induced by electrons.
81 void SetElectronSignalScalingFactor(const double scale) {
82 m_scaleElectronSignal = scale;
83 }
84 /// Set multiplication factor for the signal induced by holes.
85 void SetHoleSignalScalingFactor(const double scale) {
86 m_scaleHoleSignal = scale;
87 }
88 /// Set multiplication factor for the signal induced by ions.
89 void SetIonSignalScalingFactor(const double scale) {
90 m_scaleIonSignal = scale;
91 }
92
93 /// Return the number of electrons and ions/holes in the avalanche.
94 void GetAvalancheSize(int& ne, int& ni) const {
95 ne = m_nElectrons;
96 ni = m_nIons;
97 }
98
99 /// Return the number of points along the last simulated drift line.
100 unsigned int GetNumberOfDriftLinePoints() const { return m_drift.size(); }
101 /// Return the coordinates and time of a point along the last drift line.
102 void GetDriftLinePoint(const unsigned int i, double& x, double& y, double& z,
103 double& t) const;
104
105 /** Return the number of electron trajectories in the last
106 * simulated avalanche (including captured electrons). */
107 unsigned int GetNumberOfElectronEndpoints() const {
108 return m_endpointsElectrons.size();
109 }
110 unsigned int GetNumberOfHoleEndpoints() const {
111 return m_endpointsHoles.size();
112 }
113 unsigned int GetNumberOfIonEndpoints() const {
114 return m_endpointsIons.size();
115 }
116
117 /** Return the coordinates and time of start and end point of a given
118 * electron drift line.
119 * \param i index of the drift line
120 * \param x0,y0,z0,t0 coordinates and time of the starting point
121 * \param x1,y1,z1,t1 coordinates and time of the end point
122 * \param status status code (see GarfieldConstants.hh)
123 */
124 void GetElectronEndpoint(const unsigned int i, double& x0, double& y0,
125 double& z0, double& t0, double& x1, double& y1,
126 double& z1, double& t1, int& status) const;
127 void GetHoleEndpoint(const unsigned int i, double& x0, double& y0, double& z0,
128 double& t0, double& x1, double& y1, double& z1,
129 double& t1, int& status) const;
130 void GetIonEndpoint(const unsigned int i, double& x0, double& y0, double& z0,
131 double& t0, double& x1, double& y1, double& z1,
132 double& t1, int& status) const;
133
134 /// Simulate the drift line of an electron with a given starting point.
135 bool DriftElectron(const double x0, const double y0, const double z0,
136 const double t0);
137 bool DriftHole(const double x0, const double y0, const double z0,
138 const double t0);
139 bool DriftIon(const double x0, const double y0, const double z0,
140 const double t0);
141 /// Simulate an avalanche initiated by an electron with given starting point.
142 bool AvalancheElectron(const double x0, const double y0, const double z0,
143 const double t0, const bool hole = false);
144 bool AvalancheHole(const double x0, const double y0, const double z0,
145 const double t0, const bool electron = false);
146 bool AvalancheElectronHole(const double x0, const double y0, const double z0,
147 const double t0);
148
149 /// Switch on debugging messages (default: off).
150 void EnableDebugging() { m_debug = true; }
151 void DisableDebugging() { m_debug = false; }
152
153 private:
154 std::string m_className;
155
156 /// Numerical prefactor
157 static double c1;
158
159 Sensor* m_sensor;
160
161 struct DriftPoint {
162 double x, y, z, t; //< Position.
163 int ne, nh, ni; //< Number of secondaries produced at this point.
164 };
165 /// Current drift line
166 std::vector<DriftPoint> m_drift;
167
168 struct AvalPoint {
169 double x, y, z, t;
170 int ne, nh, ni;
171 };
172
173 /// Step size model
174 int m_stepModel;
175 /// Fixed time step
176 double m_tMc;
177 /// Fixed distance step
178 double m_dMc;
179 /// Sample step size according to collision time
180 int m_nMc;
181
182 /// Flag whether a time window should be used.
183 bool m_hasTimeWindow;
184 /// Time window
185 double m_tMin, m_tMax;
186
187 /// Number of electrons produced
188 unsigned int m_nElectrons;
189 /// Number of holes produced
190 unsigned int m_nHoles;
191 /// Number of ions produced
192 unsigned int m_nIons;
193
194 struct EndPoint {
195 double x0, y0, z0, t0; //< Starting point.
196 double x1, y1, z1, t1; //< End point.
197 int status; //< Status flag at the end point.
198 };
199 /// Endpoints of all electrons in the avalanche (including captured ones)
200 std::vector<EndPoint> m_endpointsElectrons;
201 /// Endpoints of all holes in the avalanche (including captured ones)
202 std::vector<EndPoint> m_endpointsHoles;
203 /// Endpoints of all ions in the avalanche
204 std::vector<EndPoint> m_endpointsIons;
205
206 ViewDrift* m_viewer;
207
208 bool m_useSignal;
209 bool m_useInducedCharge;
210 bool m_useEquilibration;
211 bool m_useDiffusion;
212 bool m_useAttachment;
213 bool m_useBfield;
214 bool m_useIons;
215 bool m_withElectrons;
216 bool m_withHoles;
217 double m_scaleElectronSignal;
218 double m_scaleHoleSignal;
219 double m_scaleIonSignal;
220
221 /// Use traps from the field component, ComponentTcad2d;
222 bool m_useTcadTrapping;
223 /// Take the velocity from the field component, ComponentTcad2d;
224 bool m_useTcadVelocity;
225
226 bool m_debug;
227
228 /// Compute a drift line with starting point (x0, y0, z0).
229 bool DriftLine(const double x0, const double y0, const double z0,
230 const double t0, const int type, const bool aval = false);
231 /// Compute an avalanche with starting point (x0, y0, z0).
232 bool Avalanche(const double x0, const double y0, const double z0,
233 const double t0, const unsigned int ne, const unsigned int nh,
234 const unsigned int ni);
235
236 /// Compute electric and magnetic field at a given position.
237 int GetField(const double x, const double y, const double z,
238 double& ex, double& ey, double& ez,
239 double& bx, double& by, double& bz, Medium*& medium);
240 /// Compute the drift velocity.
241 bool GetVelocity(const int type, Medium* medium,
242 const double x,const double y, const double z,
243 const double ex,const double ey, const double ez,
244 const double bx,const double by, const double bz,
245 double& vx, double& vy, double& vz);
246 /// Add a diffusion step.
247 bool AddDiffusion(const int type, Medium* medium, const double step,
248 double& x, double& y, double& z,
249 const double vx, const double vy, const double vz,
250 const double ex, const double ey, const double ez,
251 const double bx, const double by, const double bz);
252 /// Terminate a drift line close to the boundary.
253 void TerminateLine(double x0, double y0, double z0, double t0,
254 double& x, double& y, double& z, double& t) const;
255 /// Compute multiplication and losses along the current drift line.
256 bool ComputeGainLoss(const int type, int& status);
257 /// Compute Townsend and attachment coefficients along the current drift line.
258 bool ComputeAlphaEta(const int q, std::vector<double>& alphas,
259 std::vector<double>& etas) const;
260 bool Equilibrate(std::vector<double>& alphas) const;
261 /// Compute the induced signal for the current drift line.
262 void ComputeSignal(const double q) const;
263 /// Compute the induced charge for the current drift line.
264 void ComputeInducedCharge(const double q) const;
265};
266}
267
268#endif
void EnableInducedChargeCalculation()
Switch on calculation of induced charge (default: disabled).
Definition: AvalancheMC.hh:35
void GetDriftLinePoint(const unsigned int i, double &x, double &y, double &z, double &t) const
Return the coordinates and time of a point along the last drift line.
Definition: AvalancheMC.cc:130
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
Definition: AvalancheMC.cc:59
void EnableDebugging()
Switch on debugging messages (default: off).
Definition: AvalancheMC.hh:150
void SetTimeSteps(const double d=0.02)
Use fixed-time steps (default 20 ps).
Definition: AvalancheMC.cc:69
bool AvalancheElectronHole(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:527
bool DriftIon(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:249
unsigned int GetNumberOfIonEndpoints() const
Definition: AvalancheMC.hh:113
void SetIonSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by ions.
Definition: AvalancheMC.hh:89
void DisableProjectedPathIntegration()
Definition: AvalancheMC.hh:41
void SetTimeWindow(const double t0, const double t1)
Define a time interval (only carriers inside the interval are drifted).
Definition: AvalancheMC.cc:117
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:167
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an electron with a given starting point.
Definition: AvalancheMC.cc:210
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:188
void SetCollisionSteps(const int n=100)
Definition: AvalancheMC.cc:101
void DisableSignalCalculation()
Definition: AvalancheMC.hh:32
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:145
~AvalancheMC()
Destructor.
Definition: AvalancheMC.hh:21
AvalancheMC()
Constructor.
Definition: AvalancheMC.cc:15
void SetSensor(Sensor *s)
Set the sensor.
Definition: AvalancheMC.cc:49
void EnableTcadVelocity()
Switch on TCAD velocity maps.
Definition: AvalancheMC.hh:57
void EnableDiffusion()
Switch on diffusion (default: enabled)
Definition: AvalancheMC.hh:44
void EnableMagneticField()
Enable use of magnetic field in stepping algorithm.
Definition: AvalancheMC.hh:61
void GetAvalancheSize(int &ne, int &ni) const
Return the number of electrons and ions/holes in the avalanche.
Definition: AvalancheMC.hh:94
unsigned int GetNumberOfHoleEndpoints() const
Definition: AvalancheMC.hh:110
void EnableTcadTraps()
Switch on calculating trapping with TCAD traps.
Definition: AvalancheMC.hh:53
void SetElectronSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by electrons.
Definition: AvalancheMC.hh:81
void SetHoleSignalScalingFactor(const double scale)
Set multiplication factor for the signal induced by holes.
Definition: AvalancheMC.hh:85
void SetHoles()
Treat positive charge carriers as holes (default: ions).
Definition: AvalancheMC.hh:77
unsigned int GetNumberOfDriftLinePoints() const
Return the number of points along the last simulated drift line.
Definition: AvalancheMC.hh:100
unsigned int GetNumberOfElectronEndpoints() const
Definition: AvalancheMC.hh:107
bool AvalancheElectron(const double x0, const double y0, const double z0, const double t0, const bool hole=false)
Simulate an avalanche initiated by an electron with given starting point.
Definition: AvalancheMC.cc:511
void DisableInducedChargeCalculation()
Definition: AvalancheMC.hh:36
void EnableProjectedPathIntegration()
Definition: AvalancheMC.hh:40
void SetDistanceSteps(const double d=0.001)
Use fixed distance steps (default 10 um).
Definition: AvalancheMC.cc:85
void EnableSignalCalculation()
Switch on calculation of induced currents (default: disabled).
Definition: AvalancheMC.hh:31
bool AvalancheHole(const double x0, const double y0, const double z0, const double t0, const bool electron=false)
Definition: AvalancheMC.cc:519
bool DriftHole(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:230
Visualize drift lines and tracks.
Definition: ViewDrift.hh:15