Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
AvalancheMicroscopic.hh
Go to the documentation of this file.
1#ifndef G_AVALANCHE_MICROSCOPIC_H
2#define G_AVALANCHE_MICROSCOPIC_H
3
4#include <string>
5#include <vector>
6
7#include <TH1.h>
8
10#include "Sensor.hh"
11#include "ViewDrift.hh"
12
13namespace Garfield {
14
15/// Calculate electron drift lines and avalanches using microscopic tracking.
16
18 public:
19 /// Constructor
21 /// Destructor
23
24 /// Set the sensor.
25 void SetSensor(Sensor* sensor);
26
27 /// Switch on drift line plotting.
28 void EnablePlotting(ViewDrift* view);
29 /// Switch off drift line plotting.
30 void DisablePlotting() { m_viewer = nullptr; }
31 /// Draw a marker at every excitation or not.
32 void EnableExcitationMarkers(const bool on = true) {
33 m_plotExcitations = on;
34 }
35 /// Draw a marker at every ionising collision or not.
36 void EnableIonisationMarkers(const bool on = true) {
37 m_plotIonisations = on;
38 }
39 /// Draw a marker at every attachment or not.
40 void EnableAttachmentMarkers(const bool on = true) {
41 m_plotAttachments = on;
42 }
43
44 /// Switch on calculation of induced currents (default: off).
45 void EnableSignalCalculation(const bool on = true) { m_doSignal = on; }
46 /// Use the weighting potential (as opposed to the weighting field)
47 /// for calculating the induced current.
48 void UseWeightingPotential(const bool on = true) {
49 m_useWeightingPotential = on;
50 }
51 /// Integrate the weighting field over a drift line step when
52 /// calculating the induced current (default: off).
53 void EnableWeightingFieldIntegration(const bool on = true) {
54 m_integrateWeightingField = on;
55 }
56
57 /// Switch on calculation of the total induced charge (default: off).
58 void EnableInducedChargeCalculation(const bool on = true) {
59 m_doInducedCharge = on;
60 }
61
62 /// Fill a histogram with the electron energy distribution.
63 void EnableElectronEnergyHistogramming(TH1* histo);
64 /// Stop histogramming the electron energy distribution.
65 void DisableElectronEnergyHistogramming() { m_histElectronEnergy = nullptr; }
66 /// Fill a histogram with the hole energy distribution.
67 void EnableHoleEnergyHistogramming(TH1* histo);
68 /// Stop histogramming the hole energy distribution.
69 void DisableHoleEnergyHistogramming() { m_histHoleEnergy = nullptr; }
70
71 /** Fill histograms of the distance between successive collisions.
72 * \param histo
73 pointer to the histogram to be filled
74 * \param opt
75 direction ('x', 'y', 'z', 'r')
76 */
77 void SetDistanceHistogram(TH1* histo, const char opt = 'r');
78 /// Fill distance distribution histograms for a given collision type.
79 void EnableDistanceHistogramming(const int type);
80 /// Stop filling distance distribution histograms for a given collision type.
81 void DisableDistanceHistogramming(const int type);
82 /// Stop filling distance distribution histograms.
84 /// Fill histograms of the energy of electrons emitted in ionising collisions.
86 /// Stop histogramming the secondary electron energy distribution.
87 void DisableSecondaryEnergyHistogramming() { m_histSecondary = nullptr; }
88
89 /// Switch on storage of drift lines (default: off).
90 void EnableDriftLines(const bool on = true) { m_useDriftLines = on; }
91
92 /** Switch on photon transport.
93 * \remark This feature has not been tested thoroughly. */
94 void EnablePhotonTransport(const bool on = true) { m_usePhotons = on; }
95
96 /// Switch on stepping according to band structure E(k), for semiconductors.
97 void EnableBandStructure(const bool on = true) {
98 m_useBandStructureDefault = on;
99 }
100
101 /// Switch on update of coordinates for null-collision steps (default: off).
102 void EnableNullCollisionSteps(const bool on = true) {
103 m_useNullCollisionSteps = on;
104 }
105
106 /** Set a (lower) energy threshold for electron transport.
107 * This can be useful for simulating delta electrons. */
108 void SetElectronTransportCut(const double cut) { m_deltaCut = cut; }
109 /// Retrieve the value of the energy threshold.
110 double GetElectronTransportCut() const { return m_deltaCut; }
111
112 /// Set an energy threshold for photon transport.
113 void SetPhotonTransportCut(const double cut) { m_gammaCut = cut; }
114 /// Retrieve the energy threshold for transporting photons.
115 double GetPhotonTransportCut() const { return m_gammaCut; }
116
117 /** Set a max. avalanche size (i. e. ignore ionising collisions
118 once this size has been reached). */
119 void EnableAvalancheSizeLimit(const unsigned int size) { m_sizeCut = size; }
120 /// Do not apply a limit on the avalanche size.
121 void DisableAvalancheSizeLimit() { m_sizeCut = 0; }
122 /// Retrieve the currently set size limit.
123 int GetAvalancheSizeLimit() const { return m_sizeCut; }
124
125 /// Enable magnetic field in stepping algorithm (default: off).
126 void EnableMagneticField(const bool on = true) { m_useBfield = on; }
127
128 /// Set number of collisions to be skipped for plotting
129 void SetCollisionSteps(const unsigned int n) { m_nCollSkip = n; }
130
131 /// Define a time interval (only carriers inside the interval are simulated).
132 void SetTimeWindow(const double t0, const double t1);
133 /// Do not restrict the time interval within which carriers are simulated.
134 void UnsetTimeWindow() { m_hasTimeWindow = false; }
135
136 /// Return the number of electrons and ions in the avalanche.
137 void GetAvalancheSize(int& ne, int& ni) const {
138 ne = m_nElectrons;
139 ni = m_nIons;
140 }
141 void GetAvalancheSize(int& ne, int& nh, int& ni) const {
142 ne = m_nElectrons;
143 nh = m_nHoles;
144 ni = m_nIons;
145 }
146
147 /** Return the number of electron trajectories in the last
148 * simulated avalanche (including captured electrons). */
149 unsigned int GetNumberOfElectronEndpoints() const {
150 return m_endpointsElectrons.size();
151 }
152 /** Return the coordinates and time of start and end point of a given
153 * electron drift line.
154 * \param i index of the drift line
155 * \param x0,y0,z0,t0 coordinates and time of the starting point
156 * \param x1,y1,z1,t1 coordinates and time of the end point
157 * \param e0,e1 initial and final energy
158 * \param status status code (see GarfieldConstants.hh)
159 */
160 void GetElectronEndpoint(const unsigned int i, double& x0, double& y0,
161 double& z0, double& t0, double& e0, double& x1,
162 double& y1, double& z1, double& t1, double& e1,
163 int& status) const;
164 void GetElectronEndpoint(const unsigned int i, double& x0, double& y0,
165 double& z0, double& t0, double& e0, double& x1,
166 double& y1, double& z1, double& t1, double& e1,
167 double& dx1, double& dy1, double& dz1,
168 int& status) const;
170 const unsigned int i = 0) const;
171 unsigned int GetNumberOfHoleDriftLinePoints(const unsigned int i = 0) const;
172 void GetElectronDriftLinePoint(double& x, double& y, double& z, double& t,
173 const int ip,
174 const unsigned int iel = 0) const;
175 void GetHoleDriftLinePoint(double& x, double& y, double& z, double& t,
176 const int ip, const unsigned int iel = 0) const;
177
178 unsigned int GetNumberOfHoleEndpoints() const {
179 return m_endpointsHoles.size();
180 }
181 void GetHoleEndpoint(const unsigned int i, double& x0, double& y0, double& z0,
182 double& t0, double& e0, double& x1, double& y1,
183 double& z1, double& t1, double& e1, int& status) const;
184
185 unsigned int GetNumberOfPhotons() const { return m_photons.size(); }
186 // Status codes:
187 // -2: photon absorbed by gas molecule
188 void GetPhoton(const unsigned int i, double& e, double& x0, double& y0,
189 double& z0, double& t0, double& x1, double& y1, double& z1,
190 double& t1, int& status) const;
191
192 /** Calculate an electron drift line.
193 * \param x0,y0,z0,t0 starting point of the electron
194 * \param e0 initial energy of the electron
195 * \param dx0,dy0,dz0 initial direction of the electron
196 * If the initial direction is not specified, it is sampled randomly.
197 * Secondary electrons are not transported. */
198 bool DriftElectron(const double x0, const double y0, const double z0,
199 const double t0, const double e0, const double dx0 = 0.,
200 const double dy0 = 0., const double dz0 = 0.);
201
202 /// Calculate an avalanche initiated by a given electron.
203 bool AvalancheElectron(const double x0, const double y0, const double z0,
204 const double t0, const double e0,
205 const double dx0 = 0., const double dy0 = 0.,
206 const double dz0 = 0.);
207
208 /// Set a user handling procedure. This function is called at every step.
209 void SetUserHandleStep(void (*f)(double x, double y, double z, double t,
210 double e, double dx, double dy, double dz,
211 bool hole));
212 /// Deactivate the user handle called at every step.
213 void UnsetUserHandleStep() { m_userHandleStep = nullptr; }
214 /// Set a user handling procedure, to be called at every (real) collision.
215 void SetUserHandleCollision(void (*f)(double x, double y, double z, double t,
216 int type, int level, Medium* m,
217 double e0, double e1,
218 double dx0, double dy0, double dz0,
219 double dx1, double dy1, double dz1));
220 /// Deactivate the user handle called at every collision.
221 void UnsetUserHandleCollision() { m_userHandleCollision = nullptr; }
222 /// Set a user handling procedure, to be called at every attachment.
223 void SetUserHandleAttachment(void (*f)(double x, double y, double z, double t,
224 int type, int level, Medium* m));
225 /// Deactivate the user handle called at every attachment.
226 void UnsetUserHandleAttachment() { m_userHandleAttachment = nullptr; }
227 /// Set a user handling procedure, to be called at every inelastic collision.
228 void SetUserHandleInelastic(void (*f)(double x, double y, double z, double t,
229 int type, int level, Medium* m));
230 /// Deactivate the user handle called at every inelastic collision.
231 void UnsetUserHandleInelastic() { m_userHandleInelastic = nullptr; }
232 /// Set a user handling procedure, to be called at every ionising collision
233 /// or excitation followed by Penning transfer.
234 void SetUserHandleIonisation(void (*f)(double x, double y, double z, double t,
235 int type, int level, Medium* m));
236 /// Deactivate the user handle called at every ionisation.
237 void UnsetUserHandleIonisation() { m_userHandleIonisation = nullptr; }
238
239 /// Switch on debugging messages.
240 void EnableDebugging() { m_debug = true; }
241 void DisableDebugging() { m_debug = false; }
242
243 private:
244 std::string m_className = "AvalancheMicroscopic";
245
246 Sensor* m_sensor = nullptr;
247
248 struct point {
249 double x, y, z, t;
250 };
251
252 struct Electron {
253 int status; //< Status.
254 bool hole; //< Electron or hole.
255 double x0, y0, z0, t0; //< Starting point and time.
256 double e0; //< Initial kinetic energy.
257 int band; //< Band.
258 double x, y, z, t; //< Current position and time.
259 double kx, ky, kz; //< Current direction/wave vector.
260 double energy; //< Current kinetic energy.
261 std::vector<point> driftLine; //< Drift line.
262 double xLast, yLast, zLast; //< Previous position.
263 };
264 std::vector<Electron> m_endpointsElectrons;
265 std::vector<Electron> m_endpointsHoles;
266
267 struct photon {
268 int status; //< Status
269 double energy; //< Energy
270 double x0, y0, z0, t0; //< Starting point and time.
271 double x1, y1, z1, t1; //< End point and time.
272 };
273 std::vector<photon> m_photons;
274
275 /// Number of electrons produced
276 int m_nElectrons = 0;
277 /// Number of holes produced
278 int m_nHoles = 0;
279 /// Number of ions produced
280 int m_nIons = 0;
281
282 ViewDrift* m_viewer = nullptr;
283 bool m_plotExcitations = true;
284 bool m_plotIonisations = true;
285 bool m_plotAttachments = true;
286
287 TH1* m_histElectronEnergy = nullptr;
288 TH1* m_histHoleEnergy = nullptr;
289 TH1* m_histDistance = nullptr;
290 char m_distanceOption = 'r';
291 std::vector<int> m_distanceHistogramType;
292
293 TH1* m_histSecondary = nullptr;
294
295 bool m_doSignal = false;
296 bool m_useWeightingPotential = false;
297 bool m_integrateWeightingField = false;
298 bool m_doInducedCharge = false;
299 bool m_useDriftLines = false;
300 bool m_usePhotons = false;
301 bool m_useBandStructureDefault = true;
302 bool m_useNullCollisionSteps = false;
303 bool m_useBfield = false;
304
305 // Transport cuts
306 double m_deltaCut = 0.;
307 double m_gammaCut = 0.;
308
309 // Max. avalanche size
310 unsigned int m_sizeCut = 0;
311
312 unsigned int m_nCollSkip = 100;
313
314 bool m_hasTimeWindow = false;
315 double m_tMin = 0.;
316 double m_tMax = 0.;
317
318 // User procedures
319 void (*m_userHandleStep)(double x, double y, double z, double t, double e,
320 double dx, double dy, double dz,
321 bool hole) = nullptr;
322 void (*m_userHandleCollision)(double x, double y, double z, double t,
323 int type, int level, Medium* m,
324 double e0, double e1,
325 double dx0, double dy0, double dz0,
326 double dx1, double dy1, double dz1) = nullptr;
327 void (*m_userHandleAttachment)(double x, double y, double z, double t,
328 int type, int level, Medium* m) = nullptr;
329 void (*m_userHandleInelastic)(double x, double y, double z, double t,
330 int type, int level, Medium* m) = nullptr;
331 void (*m_userHandleIonisation)(double x, double y, double z, double t,
332 int type, int level, Medium* m) = nullptr;
333
334 // Switch on/off debugging messages
335 bool m_debug = false;
336
337 // Electron transport
338 bool TransportElectron(const double x0, const double y0, const double z0,
339 const double t0, double e0, const double dx0,
340 const double dy0, const double dz0, const bool aval,
341 bool hole);
342 // Photon transport
343 void TransportPhoton(const double x, const double y, const double z,
344 const double t, const double e,
345 std::vector<Electron>& stack);
346
347 static bool IsInactive(const Electron& item) {
348 return item.status == StatusLeftDriftMedium ||
349 item.status == StatusBelowTransportCut ||
350 item.status == StatusOutsideTimeWindow ||
351 item.status == StatusLeftDriftArea || item.status == StatusAttached;
352 }
353 void Update(std::vector<Electron>::iterator it, const double x,
354 const double y, const double z, const double t,
355 const double energy, const double kx, const double ky,
356 const double kz, const int band);
357 void AddToEndPoints(const Electron& item, const bool hole) {
358 if (hole) {
359 m_endpointsHoles.push_back(item);
360 } else {
361 m_endpointsElectrons.push_back(item);
362 }
363 }
364
365 /// Add a new electron/hole (with random direction) to a container.
366 void AddToStack(const double x, const double y, const double z,
367 const double t, const double energy, const bool hole,
368 std::vector<Electron>& container) const;
369 /// Add a new electron/hole to a container.
370 void AddToStack(const double x, const double y, const double z,
371 const double t, const double energy, const double dx,
372 const double dy, const double dz, const int band,
373 const bool hole, std::vector<Electron>& container) const;
374 void Terminate(double x0, double y0, double z0, double t0, double& x1,
375 double& y1, double& z1, double& t1);
376};
377}
378
379#endif
Calculate electron drift lines and avalanches using microscopic tracking.
void UnsetUserHandleCollision()
Deactivate the user handle called at every collision.
void EnableDistanceHistogramming(const int type)
Fill distance distribution histograms for a given collision type.
void UnsetTimeWindow()
Do not restrict the time interval within which carriers are simulated.
void EnableInducedChargeCalculation(const bool on=true)
Switch on calculation of the total induced charge (default: off).
void UseWeightingPotential(const bool on=true)
void SetDistanceHistogram(TH1 *histo, const char opt='r')
void EnableWeightingFieldIntegration(const bool on=true)
void EnableAvalancheSizeLimit(const unsigned int size)
void SetUserHandleStep(void(*f)(double x, double y, double z, double t, double e, double dx, double dy, double dz, bool hole))
Set a user handling procedure. This function is called at every step.
void EnableHoleEnergyHistogramming(TH1 *histo)
Fill a histogram with the hole energy distribution.
void EnableNullCollisionSteps(const bool on=true)
Switch on update of coordinates for null-collision steps (default: off).
void EnablePhotonTransport(const bool on=true)
void SetUserHandleCollision(void(*f)(double x, double y, double z, double t, int type, int level, Medium *m, double e0, double e1, double dx0, double dy0, double dz0, double dx1, double dy1, double dz1))
Set a user handling procedure, to be called at every (real) collision.
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
void EnableDriftLines(const bool on=true)
Switch on storage of drift lines (default: off).
void SetUserHandleIonisation(void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
void UnsetUserHandleAttachment()
Deactivate the user handle called at every attachment.
void DisableHoleEnergyHistogramming()
Stop histogramming the hole energy distribution.
bool DriftElectron(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0=0., const double dy0=0., const double dz0=0.)
void SetPhotonTransportCut(const double cut)
Set an energy threshold for photon transport.
void EnableIonisationMarkers(const bool on=true)
Draw a marker at every ionising collision or not.
void DisableSecondaryEnergyHistogramming()
Stop histogramming the secondary electron energy distribution.
unsigned int GetNumberOfElectronDriftLinePoints(const unsigned int i=0) const
unsigned int GetNumberOfHoleEndpoints() const
void DisableAvalancheSizeLimit()
Do not apply a limit on the avalanche size.
void EnableBandStructure(const bool on=true)
Switch on stepping according to band structure E(k), for semiconductors.
double GetPhotonTransportCut() const
Retrieve the energy threshold for transporting photons.
void SetCollisionSteps(const unsigned int n)
Set number of collisions to be skipped for plotting.
void GetAvalancheSize(int &ne, int &ni) const
Return the number of electrons and ions in the avalanche.
void SetUserHandleAttachment(void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
Set a user handling procedure, to be called at every attachment.
void SetSensor(Sensor *sensor)
Set the sensor.
void GetHoleDriftLinePoint(double &x, double &y, double &z, double &t, const int ip, const unsigned int iel=0) const
void GetElectronEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const
void EnableExcitationMarkers(const bool on=true)
Draw a marker at every excitation or not.
void EnableSignalCalculation(const bool on=true)
Switch on calculation of induced currents (default: off).
void GetHoleEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &e0, double &x1, double &y1, double &z1, double &t1, double &e1, int &status) const
void DisablePlotting()
Switch off drift line plotting.
void GetElectronDriftLinePoint(double &x, double &y, double &z, double &t, const int ip, const unsigned int iel=0) const
void UnsetUserHandleStep()
Deactivate the user handle called at every step.
void EnableMagneticField(const bool on=true)
Enable magnetic field in stepping algorithm (default: off).
double GetElectronTransportCut() const
Retrieve the value of the energy threshold.
void EnableElectronEnergyHistogramming(TH1 *histo)
Fill a histogram with the electron energy distribution.
void UnsetUserHandleIonisation()
Deactivate the user handle called at every ionisation.
void EnableAttachmentMarkers(const bool on=true)
Draw a marker at every attachment or not.
void DisableDistanceHistogramming()
Stop filling distance distribution histograms.
void GetPhoton(const unsigned int i, double &e, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
void GetAvalancheSize(int &ne, int &nh, int &ni) const
void EnableSecondaryEnergyHistogramming(TH1 *histo)
Fill histograms of the energy of electrons emitted in ionising collisions.
bool AvalancheElectron(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0=0., const double dy0=0., const double dz0=0.)
Calculate an avalanche initiated by a given electron.
void SetElectronTransportCut(const double cut)
int GetAvalancheSizeLimit() const
Retrieve the currently set size limit.
unsigned int GetNumberOfElectronEndpoints() const
void DisableElectronEnergyHistogramming()
Stop histogramming the electron energy distribution.
void EnableDebugging()
Switch on debugging messages.
unsigned int GetNumberOfHoleDriftLinePoints(const unsigned int i=0) const
void SetUserHandleInelastic(void(*f)(double x, double y, double z, double t, int type, int level, Medium *m))
Set a user handling procedure, to be called at every inelastic collision.
void UnsetUserHandleInelastic()
Deactivate the user handle called at every inelastic collision.
void SetTimeWindow(const double t0, const double t1)
Define a time interval (only carriers inside the interval are simulated).
Abstract base class for media.
Definition: Medium.hh:13
Visualize drift lines and tracks.
Definition: ViewDrift.hh:19