Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
AvalancheMicroscopic.cc
Go to the documentation of this file.
1#include <algorithm>
2#include <cmath>
3#include <functional>
4#include <iostream>
5#include <string>
6
9#include "Garfield/Random.hh"
10
11namespace {
12
13double Mag(const double x, const double y, const double z) {
14 return sqrt(x * x + y * y + z * z);
15}
16
17void ToLocal(const std::array<std::array<double, 3>, 3>& rot, const double xg,
18 const double yg, const double zg, double& xl, double& yl,
19 double& zl) {
20 xl = rot[0][0] * xg + rot[0][1] * yg + rot[0][2] * zg;
21 yl = rot[1][0] * xg + rot[1][1] * yg + rot[1][2] * zg;
22 zl = rot[2][0] * xg + rot[2][1] * yg + rot[2][2] * zg;
23}
24
25void ToGlobal(const std::array<std::array<double, 3>, 3>& rot, const double xl,
26 const double yl, const double zl, double& xg, double& yg,
27 double& zg) {
28 xg = rot[0][0] * xl + rot[1][0] * yl + rot[2][0] * zl;
29 yg = rot[0][1] * xl + rot[1][1] * yl + rot[2][1] * zl;
30 zg = rot[0][2] * xl + rot[1][2] * yl + rot[2][2] * zl;
31}
32
33void RotationMatrix(double bx, double by, double bz, const double bmag,
34 const double ex, const double ey, const double ez,
35 std::array<std::array<double, 3>, 3>& rot) {
36 // Adopting the Magboltz convention, the stepping is performed
37 // in a coordinate system with the B field along the x axis
38 // and the electric field at an angle btheta in the x-z plane.
39
40 // Calculate the first rotation matrix (to align B with x axis).
41 std::array<std::array<double, 3>, 3> rB = {{{1, 0, 0}, {0, 1, 0}, {0, 0, 1}}};
42 if (bmag > Garfield::Small) {
43 bx /= bmag;
44 by /= bmag;
45 bz /= bmag;
46 const double bt = by * by + bz * bz;
47 if (bt > Garfield::Small) {
48 const double btInv = 1. / bt;
49 rB[0][0] = bx;
50 rB[0][1] = by;
51 rB[0][2] = bz;
52 rB[1][0] = -by;
53 rB[2][0] = -bz;
54 rB[1][1] = (bx * by * by + bz * bz) * btInv;
55 rB[2][2] = (bx * bz * bz + by * by) * btInv;
56 rB[1][2] = rB[2][1] = (bx - 1.) * by * bz * btInv;
57 }
58 }
59 // Calculate the second rotation matrix (rotation around x axis).
60 const double fy = rB[1][0] * ex + rB[1][1] * ey + rB[1][2] * ez;
61 const double fz = rB[2][0] * ex + rB[2][1] * ey + rB[2][2] * ez;
62 const double ft = sqrt(fy * fy + fz * fz);
63 std::array<std::array<double, 3>, 3> rX = {{{1, 0, 0}, {0, 1, 0}, {0, 0, 1}}};
64 if (ft > Garfield::Small) {
65 // E and B field are not parallel.
66 rX[1][1] = rX[2][2] = fz / ft;
67 rX[1][2] = -fy / ft;
68 rX[2][1] = -rX[1][2];
69 }
70 for (unsigned int i = 0; i < 3; ++i) {
71 for (unsigned int j = 0; j < 3; ++j) {
72 rot[i][j] = 0.;
73 for (unsigned int k = 0; k < 3; ++k) {
74 rot[i][j] += rX[i][k] * rB[k][j];
75 }
76 }
77 }
78}
79
80void PrintStatus(const std::string& hdr, const std::string& status,
81 const double x, const double y, const double z,
82 const bool hole) {
83 const std::string eh = hole ? "Hole " : "Electron ";
84 std::cout << hdr << eh << status << " at " << x << ", " << y << ", " << z
85 << "\n";
86}
87} // namespace
88
89namespace Garfield {
90
92 m_endpointsElectrons.reserve(10000);
93 m_endpointsHoles.reserve(10000);
94 m_photons.reserve(1000);
95}
96
98 if (!s) {
99 std::cerr << m_className << "::SetSensor: Null pointer.\n";
100 return;
101 }
102 m_sensor = s;
103}
104
106 if (!view) {
107 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
108 return;
109 }
110
111 m_viewer = view;
112 if (!m_storeDriftLines) {
113 std::cout << m_className << "::EnablePlotting:\n"
114 << " Enabling storage of drift line.\n";
116 }
117}
118
120 if (!histo) {
121 std::cerr << m_className << "::EnableElectronEnergyHistogramming:\n"
122 << " Null pointer.\n";
123 return;
124 }
125
126 m_histElectronEnergy = histo;
127}
128
130 if (!histo) {
131 std::cerr << m_className << "::EnableHoleEnergyHistogramming:\n"
132 << " Null pointer.\n";
133 return;
134 }
135
136 m_histHoleEnergy = histo;
137}
138
139void AvalancheMicroscopic::SetDistanceHistogram(TH1* histo, const char opt) {
140 if (!histo) {
141 std::cerr << m_className << "::SetDistanceHistogram: Null pointer.\n";
142 return;
143 }
144
145 m_histDistance = histo;
146
147 if (opt == 'x' || opt == 'y' || opt == 'z' || opt == 'r') {
148 m_distanceOption = opt;
149 } else {
150 std::cerr << m_className << "::SetDistanceHistogram:";
151 std::cerr << " Unknown option " << opt << ".\n";
152 std::cerr << " Valid options are x, y, z, r.\n";
153 std::cerr << " Using default value (r).\n";
154 m_distanceOption = 'r';
155 }
156
157 if (m_distanceHistogramType.empty()) {
158 std::cout << m_className << "::SetDistanceHistogram:\n";
159 std::cout << " Don't forget to call EnableDistanceHistogramming.\n";
160 }
161}
162
164 // Check if this type of collision is already registered
165 // for histogramming.
166 const unsigned int nDistanceHistogramTypes = m_distanceHistogramType.size();
167 if (nDistanceHistogramTypes > 0) {
168 for (unsigned int i = 0; i < nDistanceHistogramTypes; ++i) {
169 if (m_distanceHistogramType[i] != type) continue;
170 std::cout << m_className << "::EnableDistanceHistogramming:\n";
171 std::cout << " Collision type " << type
172 << " is already being histogrammed.\n";
173 return;
174 }
175 }
176
177 m_distanceHistogramType.push_back(type);
178 std::cout << m_className << "::EnableDistanceHistogramming:\n";
179 std::cout << " Histogramming of collision type " << type << " enabled.\n";
180 if (!m_histDistance) {
181 std::cout << " Don't forget to set the histogram.\n";
182 }
183}
184
186 if (std::find(m_distanceHistogramType.begin(), m_distanceHistogramType.end(),
187 type) == m_distanceHistogramType.end()) {
188 std::cerr << m_className << "::DisableDistanceHistogramming:\n"
189 << " Collision type " << type << " is not histogrammed.\n";
190 return;
191 }
192
193 m_distanceHistogramType.erase(
194 std::remove(m_distanceHistogramType.begin(),
195 m_distanceHistogramType.end(), type),
196 m_distanceHistogramType.end());
197}
198
200 m_histDistance = nullptr;
201 m_distanceHistogramType.clear();
202}
203
205 if (!histo) {
206 std::cerr << m_className << "::EnableSecondaryEnergyHistogramming:\n"
207 << " Null pointer.\n";
208 return;
209 }
210
211 m_histSecondary = histo;
212}
213
214void AvalancheMicroscopic::SetTimeWindow(const double t0, const double t1) {
215 if (fabs(t1 - t0) < Small) {
216 std::cerr << m_className << "::SetTimeWindow:\n";
217 std::cerr << " Time interval must be greater than zero.\n";
218 return;
219 }
220
221 m_tMin = std::min(t0, t1);
222 m_tMax = std::max(t0, t1);
223 m_hasTimeWindow = true;
224}
225
226void AvalancheMicroscopic::GetElectronEndpoint(const unsigned int i, double& x0,
227 double& y0, double& z0,
228 double& t0, double& e0,
229 double& x1, double& y1,
230 double& z1, double& t1,
231 double& e1, int& status) const {
232 if (i >= m_endpointsElectrons.size()) {
233 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
234 x0 = y0 = z0 = t0 = e0 = 0.;
235 x1 = y1 = z1 = t1 = e1 = 0.;
236 status = 0;
237 return;
238 }
239
240 x0 = m_endpointsElectrons[i].x0;
241 y0 = m_endpointsElectrons[i].y0;
242 z0 = m_endpointsElectrons[i].z0;
243 t0 = m_endpointsElectrons[i].t0;
244 e0 = m_endpointsElectrons[i].e0;
245 x1 = m_endpointsElectrons[i].x;
246 y1 = m_endpointsElectrons[i].y;
247 z1 = m_endpointsElectrons[i].z;
248 t1 = m_endpointsElectrons[i].t;
249 e1 = m_endpointsElectrons[i].energy;
250 status = m_endpointsElectrons[i].status;
251}
252
254 const unsigned int i, double& x0, double& y0, double& z0, double& t0,
255 double& e0, double& x1, double& y1, double& z1, double& t1, double& e1,
256 double& dx1, double& dy1, double& dz1, int& status) const {
257 if (i >= m_endpointsElectrons.size()) {
258 std::cerr << m_className << "::GetElectronEndpoint: Index out of range.\n";
259 x0 = y0 = z0 = t0 = e0 = 0.;
260 x1 = y1 = z1 = t1 = e1 = 0.;
261 dx1 = dy1 = dz1 = 0.;
262 status = 0;
263 return;
264 }
265
266 x0 = m_endpointsElectrons[i].x0;
267 y0 = m_endpointsElectrons[i].y0;
268 z0 = m_endpointsElectrons[i].z0;
269 t0 = m_endpointsElectrons[i].t0;
270 e0 = m_endpointsElectrons[i].e0;
271 x1 = m_endpointsElectrons[i].x;
272 y1 = m_endpointsElectrons[i].y;
273 z1 = m_endpointsElectrons[i].z;
274 t1 = m_endpointsElectrons[i].t;
275 e1 = m_endpointsElectrons[i].energy;
276 dx1 = m_endpointsElectrons[i].kx;
277 dy1 = m_endpointsElectrons[i].ky;
278 dz1 = m_endpointsElectrons[i].kz;
279 status = m_endpointsElectrons[i].status;
280}
281
282void AvalancheMicroscopic::GetHoleEndpoint(const unsigned int i, double& x0,
283 double& y0, double& z0, double& t0,
284 double& e0, double& x1, double& y1,
285 double& z1, double& t1, double& e1,
286 int& status) const {
287 if (i >= m_endpointsHoles.size()) {
288 std::cerr << m_className << "::GetHoleEndpoint: Index out of range.\n";
289 x0 = y0 = z0 = t0 = e0 = 0.;
290 x1 = y1 = z1 = t1 = e1 = 0.;
291 status = 0;
292 return;
293 }
294
295 x0 = m_endpointsHoles[i].x0;
296 y0 = m_endpointsHoles[i].y0;
297 z0 = m_endpointsHoles[i].z0;
298 t0 = m_endpointsHoles[i].t0;
299 e0 = m_endpointsHoles[i].e0;
300 x1 = m_endpointsHoles[i].x;
301 y1 = m_endpointsHoles[i].y;
302 z1 = m_endpointsHoles[i].z;
303 t1 = m_endpointsHoles[i].t;
304 e1 = m_endpointsHoles[i].energy;
305 status = m_endpointsHoles[i].status;
306}
307
309 const unsigned int i) const {
310 if (i >= m_endpointsElectrons.size()) {
311 std::cerr << m_className << "::GetNumberOfElectronDriftLinePoints:\n";
312 std::cerr << " Endpoint index (" << i << ") out of range.\n";
313 return 0;
314 }
315
316 if (!m_storeDriftLines) return 2;
317
318 return m_endpointsElectrons[i].driftLine.size() + 2;
319}
320
322 const unsigned int i) const {
323 if (i >= m_endpointsHoles.size()) {
324 std::cerr << m_className << "::GetNumberOfHoleDriftLinePoints:\n";
325 std::cerr << " Endpoint index (" << i << ") out of range.\n";
326 return 0;
327 }
328
329 if (!m_storeDriftLines) return 2;
330
331 return m_endpointsHoles[i].driftLine.size() + 2;
332}
333
335 double& x, double& y, double& z, double& t, const int ip,
336 const unsigned int iel) const {
337 if (iel >= m_endpointsElectrons.size()) {
338 std::cerr << m_className << "::GetElectronDriftLinePoint:\n";
339 std::cerr << " Endpoint index (" << iel << ") out of range.\n";
340 return;
341 }
342
343 if (ip <= 0) {
344 x = m_endpointsElectrons[iel].x0;
345 y = m_endpointsElectrons[iel].y0;
346 z = m_endpointsElectrons[iel].z0;
347 t = m_endpointsElectrons[iel].t0;
348 return;
349 }
350
351 const int np = m_endpointsElectrons[iel].driftLine.size();
352 if (ip > np) {
353 x = m_endpointsElectrons[iel].x;
354 y = m_endpointsElectrons[iel].y;
355 z = m_endpointsElectrons[iel].z;
356 t = m_endpointsElectrons[iel].t;
357 return;
358 }
359
360 x = m_endpointsElectrons[iel].driftLine[ip - 1].x;
361 y = m_endpointsElectrons[iel].driftLine[ip - 1].y;
362 z = m_endpointsElectrons[iel].driftLine[ip - 1].z;
363 t = m_endpointsElectrons[iel].driftLine[ip - 1].t;
364}
365
367 double& z, double& t,
368 const int ip,
369 const unsigned int ih) const {
370 if (ih >= m_endpointsHoles.size()) {
371 std::cerr << m_className << "::GetHoleDriftLinePoint:\n";
372 std::cerr << " Endpoint index (" << ih << ") out of range.\n";
373 return;
374 }
375
376 if (ip <= 0) {
377 x = m_endpointsHoles[ih].x0;
378 y = m_endpointsHoles[ih].y0;
379 z = m_endpointsHoles[ih].z0;
380 t = m_endpointsHoles[ih].t0;
381 return;
382 }
383
384 const int np = m_endpointsHoles[ih].driftLine.size();
385 if (ip > np) {
386 x = m_endpointsHoles[ih].x;
387 y = m_endpointsHoles[ih].y;
388 z = m_endpointsHoles[ih].z;
389 t = m_endpointsHoles[ih].t;
390 return;
391 }
392
393 x = m_endpointsHoles[ih].driftLine[ip - 1].x;
394 y = m_endpointsHoles[ih].driftLine[ip - 1].y;
395 z = m_endpointsHoles[ih].driftLine[ip - 1].z;
396 t = m_endpointsHoles[ih].driftLine[ip - 1].t;
397}
398
399void AvalancheMicroscopic::GetPhoton(const unsigned int i, double& e,
400 double& x0, double& y0, double& z0,
401 double& t0, double& x1, double& y1,
402 double& z1, double& t1,
403 int& status) const {
404 if (i >= m_photons.size()) {
405 std::cerr << m_className << "::GetPhoton: Index out of range.\n";
406 return;
407 }
408
409 x0 = m_photons[i].x0;
410 x1 = m_photons[i].x1;
411 y0 = m_photons[i].y0;
412 y1 = m_photons[i].y1;
413 z0 = m_photons[i].z0;
414 z1 = m_photons[i].z1;
415 t0 = m_photons[i].t0;
416 t1 = m_photons[i].t1;
417 status = m_photons[i].status;
418 e = m_photons[i].energy;
419}
420
422 void (*f)(double x, double y, double z, double t, double e, double dx,
423 double dy, double dz, bool hole)) {
424 if (!f) {
425 std::cerr << m_className << "::SetUserHandleStep: Null pointer.\n";
426 return;
427 }
428 m_userHandleStep = f;
429}
430
432 void (*f)(double x, double y, double z, double t, int type, int level,
433 Medium* m, double e0, double e1, double dx0, double dy0,
434 double dz0, double dx1, double dy1, double dz1)) {
435 m_userHandleCollision = f;
436}
437
439 double x, double y, double z, double t, int type, int level, Medium* m)) {
440 m_userHandleAttachment = f;
441}
442
444 double x, double y, double z, double t, int type, int level, Medium* m)) {
445 m_userHandleInelastic = f;
446}
447
449 double x, double y, double z, double t, int type, int level, Medium* m)) {
450 m_userHandleIonisation = f;
451}
452
453bool AvalancheMicroscopic::DriftElectron(const double x0, const double y0,
454 const double z0, const double t0,
455 const double e0, const double dx0,
456 const double dy0, const double dz0) {
457 std::vector<Electron> stack;
458 AddToStack(x0, y0, z0, t0, e0, dx0, dy0, dz0, 0, false, stack);
459 return TransportElectrons(stack, false);
460}
461
462bool AvalancheMicroscopic::AvalancheElectron(const double x0, const double y0,
463 const double z0, const double t0,
464 const double e0, const double dx0,
465 const double dy0,
466 const double dz0) {
467 std::vector<Electron> stack;
468 AddToStack(x0, y0, z0, t0, e0, dx0, dy0, dz0, 0, false, stack);
469 return TransportElectrons(stack, true);
470}
471
473 std::vector<Electron> stack;
474 for (const auto& p : m_endpointsElectrons) {
475 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
476 AddToStack(p.x, p.y, p.z, p.t, p.energy,
477 p.kx, p.ky, p.kz, p.band, false, stack);
478 }
479 }
480 for (const auto& p : m_endpointsHoles) {
481 if (p.status == StatusAlive || p.status == StatusOutsideTimeWindow) {
482 AddToStack(p.x, p.y, p.z, p.t, p.energy,
483 p.kx, p.ky, p.kz, p.band, true, stack);
484 }
485 }
486 return TransportElectrons(stack, true);
487}
488
489bool AvalancheMicroscopic::TransportElectrons(std::vector<Electron>& stack,
490 const bool aval) {
491
492 // Clear the list of electrons, holes and photons.
493 m_endpointsElectrons.clear();
494 m_endpointsHoles.clear();
495 m_photons.clear();
496
497 // Reset the particle counters.
498 m_nElectrons = m_nHoles = m_nIons = 0;
499
500 const std::string hdr = m_className + "::TransportElectron: ";
501 // Make sure that the sensor is defined.
502 if (!m_sensor) {
503 std::cerr << hdr << "Sensor is not defined.\n";
504 return false;
505 }
506
507 // Loop over the initial set of electrons/holes.
508 for (auto& p : stack) {
509 // Make sure that the starting point is inside a medium.
510 Medium* medium = nullptr;
511 if (!m_sensor->GetMedium(p.x0, p.y0, p.z0, medium) || !medium) {
512 std::cerr << hdr << "No medium at initial position.\n";
513 return false;
514 }
515 // Make sure that the medium is "driftable" and microscopic.
516 if (!medium->IsDriftable() || !medium->IsMicroscopic()) {
517 std::cerr << hdr << "Medium does not have cross-section data.\n";
518 return false;
519 }
520 // Make sure the initial energy is positive.
521 p.e0 = std::max(p.e0, Small);
522 if (medium->IsSemiconductor() && m_useBandStructure) {
523 if (p.band < 0) {
524 // Sample the initial momentum and band.
525 medium->GetElectronMomentum(p.e0, p.kx, p.ky, p.kz, p.band);
526 }
527 } else {
528 p.band = 0;
529 const double kmag = Mag(p.kx, p.ky, p.kz);
530 if (fabs(kmag) < Small) {
531 // Direction has zero norm, draw a random direction.
532 RndmDirection(p.kx, p.ky, p.kz);
533 } else {
534 // Normalise the direction to 1.
535 p.kx /= kmag;
536 p.ky /= kmag;
537 p.kz /= kmag;
538 }
539 }
540 if (p.hole) {
541 ++m_nHoles;
542 } else {
543 ++m_nElectrons;
544 }
545 }
546
547 // Numerical prefactors in equation of motion
548 const double c1 = SpeedOfLight * sqrt(2. / ElectronMass);
549 const double c2 = c1 * c1 / 4.;
550
551 std::vector<Electron> stackNew;
552 stackNew.reserve(1000);
553 std::vector<std::pair<double, double> > stackPhotons;
554 std::vector<std::pair<int, double> > secondaries;
555
556 while (true) {
557 // Remove all inactive items from the stack.
558 stack.erase(std::remove_if(stack.begin(), stack.end(), IsInactive),
559 stack.end());
560 // Add the secondaries produced in the last iteration.
561 if (aval && m_sizeCut > 0) {
562 // If needed, reduce the number of electrons to add.
563 if (stack.size() > m_sizeCut) {
564 stackNew.clear();
565 } else if (stack.size() + stackNew.size() > m_sizeCut) {
566 stackNew.resize(m_sizeCut - stack.size());
567 }
568 }
569 stack.insert(stack.end(), std::make_move_iterator(stackNew.begin()),
570 std::make_move_iterator(stackNew.end()));
571 stackNew.clear();
572 // If the list of electrons/holes is exhausted, we're done.
573 if (stack.empty()) break;
574 // Loop over all electrons/holes in the avalanche.
575 for (auto it = stack.begin(), end = stack.end(); it != end; ++it) {
576 // Get an electron/hole from the stack.
577 double x = (*it).x;
578 double y = (*it).y;
579 double z = (*it).z;
580 double t = (*it).t;
581 double en = (*it).energy;
582 int band = (*it).band;
583 double kx = (*it).kx;
584 double ky = (*it).ky;
585 double kz = (*it).kz;
586 bool hole = (*it).hole;
587
588 bool ok = true;
589
590 // Count number of collisions between updates.
591 unsigned int nCollTemp = 0;
592
593 // Get the local electric field and medium.
594 double ex = 0., ey = 0., ez = 0.;
595 Medium* medium = nullptr;
596 int status = 0;
597 m_sensor->ElectricField(x, y, z, ex, ey, ez, medium, status);
598 // Sign change for electrons.
599 if (!hole) {
600 ex = -ex;
601 ey = -ey;
602 ez = -ez;
603 }
604 if (m_debug) {
605 const std::string eh = hole ? "hole " : "electron ";
606 std::cout << hdr + "Drifting " + eh << it - stack.begin() << ".\n"
607 << " Field [V/cm] at (" << x << ", " << y << ", " << z
608 << "): " << ex << ", " << ey << ", " << ez
609 << "\n Status: " << status << "\n";
610 if (medium) std::cout << " Medium: " << medium->GetName() << "\n";
611 }
612
613 if (status != 0 || !medium || !medium->IsDriftable() ||
614 !medium->IsMicroscopic()) {
615 // Electron is not inside a microscopic drift medium.
616 Update(it, x, y, z, t, en, kx, ky, kz, band);
617 (*it).status = StatusLeftDriftMedium;
618 AddToEndPoints(*it, hole);
619 if (m_debug) PrintStatus(hdr, "left the drift medium", x, y, z, hole);
620 continue;
621 }
622 // Get the id number of the drift medium.
623 auto id = medium->GetId();
624 bool sc = (medium->IsSemiconductor() && m_useBandStructure);
625 // Get the null-collision rate.
626 double fLim = medium->GetElectronNullCollisionRate(band);
627 if (fLim <= 0.) {
628 std::cerr << hdr << "Got null-collision rate <= 0.\n";
629 return false;
630 }
631 double fInv = 1. / fLim;
632
633 // If switched on, get the local magnetic field.
634 double bx = 0., by = 0., bz = 0.;
635 // Cyclotron frequency.
636 double omega = 0.;
637 // Ratio of transverse electric field component and magnetic field.
638 double ezovb = 0.;
639 std::array<std::array<double, 3>, 3> rot;
640 if (m_useBfield) {
641 m_sensor->MagneticField(x, y, z, bx, by, bz, status);
642 const double scale = hole ? Tesla2Internal : -Tesla2Internal;
643 bx *= scale;
644 by *= scale;
645 bz *= scale;
646 const double bmag = Mag(bx, by, bz);
647 // Calculate the rotation matrix to a local coordinate system
648 // with B along x and E in the x-z plane.
649 RotationMatrix(bx, by, bz, bmag, ex, ey, ez, rot);
650 // Calculate the cyclotron frequency.
651 omega = OmegaCyclotronOverB * bmag;
652 // Calculate the electric field in the local frame.
653 ToLocal(rot, ex, ey, ez, ex, ey, ez);
654 ezovb = bmag > Small ? ez / bmag : 0.;
655 }
656
657 // Trace the electron/hole.
658 while (1) {
659 bool isNullCollision = false;
660
661 // Make sure the electron energy exceeds the transport cut.
662 if (en < m_deltaCut) {
663 Update(it, x, y, z, t, en, kx, ky, kz, band);
664 (*it).status = StatusBelowTransportCut;
665 AddToEndPoints(*it, hole);
666 if (m_debug) {
667 std::cout << hdr << "Kinetic energy (" << en
668 << ") below transport cut.\n";
669 }
670 ok = false;
671 break;
672 }
673
674 // Fill the energy distribution histogram.
675 if (hole && m_histHoleEnergy) {
676 m_histHoleEnergy->Fill(en);
677 } else if (!hole && m_histElectronEnergy) {
678 m_histElectronEnergy->Fill(en);
679 }
680
681 // Check if the electron is within the specified time window.
682 if (m_hasTimeWindow && (t < m_tMin || t > m_tMax)) {
683 Update(it, x, y, z, t, en, kx, ky, kz, band);
684 (*it).status = StatusOutsideTimeWindow;
685 AddToEndPoints(*it, hole);
686 if (m_debug) PrintStatus(hdr, "left the time window", x, y, z, hole);
687 ok = false;
688 break;
689 }
690
691 if (medium->GetId() != id) {
692 // Medium has changed.
693 if (!medium->IsMicroscopic()) {
694 // Electron/hole has left the microscopic drift medium.
695 Update(it, x, y, z, t, en, kx, ky, kz, band);
696 (*it).status = StatusLeftDriftMedium;
697 AddToEndPoints(*it, hole);
698 ok = false;
699 if (m_debug) {
700 std::cout << hdr << "\n Medium at " << x << ", " << y << ", "
701 << z << " does not have cross-section data.\n";
702 }
703 break;
704 }
705 id = medium->GetId();
706 sc = (medium->IsSemiconductor() && m_useBandStructure);
707 // Update the null-collision rate.
708 fLim = medium->GetElectronNullCollisionRate(band);
709 if (fLim <= 0.) {
710 std::cerr << hdr << "Got null-collision rate <= 0.\n";
711 return false;
712 }
713 fInv = 1. / fLim;
714 }
715
716 double a1 = 0., a2 = 0.;
717 // Initial velocity.
718 double vx = 0., vy = 0., vz = 0.;
719 if (m_useBfield) {
720 // Calculate the velocity vector in the local frame.
721 const double vmag = c1 * sqrt(en);
722 ToLocal(rot, vmag * kx, vmag * ky, vmag * kz, vx, vy, vz);
723 a1 = vx * ex;
724 a2 = c2 * ex * ex;
725 if (omega > Small) {
726 vy -= ezovb;
727 } else {
728 a1 += vz * ez;
729 a2 += c2 * ez * ez;
730 }
731 } else if (sc) {
732 en = medium->GetElectronEnergy(kx, ky, kz, vx, vy, vz, band);
733 } else {
734 // No band structure, no magnetic field.
735 // Calculate the velocity vector.
736 const double vmag = c1 * sqrt(en);
737 vx = vmag * kx;
738 vy = vmag * ky;
739 vz = vmag * kz;
740 a1 = vx * ex + vy * ey + vz * ez;
741 a2 = c2 * (ex * ex + ey * ey + ez * ez);
742 }
743
744 if (m_userHandleStep) {
745 m_userHandleStep(x, y, z, t, en, kx, ky, kz, hole);
746 }
747
748 // Energy after the step.
749 double en1 = en;
750 // Determine the timestep.
751 double dt = 0.;
752 // Parameters for B-field stepping.
753 double cphi = 1., sphi = 0.;
754 double a3 = 0., a4 = 0.;
755 while (1) {
756 // Sample the flight time.
757 const double r = RndmUniformPos();
758 dt += -log(r) * fInv;
759 // Calculate the energy after the proposed step.
760 if (m_useBfield) {
761 en1 = en + (a1 + a2 * dt) * dt;
762 if (omega > Small) {
763 cphi = cos(omega * dt);
764 sphi = sin(omega * dt);
765 a3 = sphi / omega;
766 a4 = (1. - cphi) / omega;
767 en1 += ez * (vz * a3 - vy * a4);
768 }
769 } else if (sc) {
770 const double cdt = dt * SpeedOfLight;
771 const double kx1 = kx + ex * cdt;
772 const double ky1 = ky + ey * cdt;
773 const double kz1 = kz + ez * cdt;
774 double vx1 = 0., vy1 = 0., vz1 = 0.;
775 en1 = medium->GetElectronEnergy(kx1, ky1, kz1, vx1, vy1, vz1, band);
776 } else {
777 en1 = en + (a1 + a2 * dt) * dt;
778 }
779 en1 = std::max(en1, Small);
780 // Get the real collision rate at the updated energy.
781 double fReal = medium->GetElectronCollisionRate(en1, band);
782 if (fReal <= 0.) {
783 std::cerr << hdr << "Got collision rate <= 0 at " << en1
784 << " eV (band " << band << ").\n";
785 return false;
786 }
787 if (fReal > fLim) {
788 // Real collision rate is higher than null-collision rate.
789 dt += log(r) * fInv;
790 // Increase the null collision rate and try again.
791 std::cerr << hdr << "Increasing null-collision rate by 5%.\n";
792 if (sc) std::cerr << " Band " << band << "\n";
793 fLim *= 1.05;
794 fInv = 1. / fLim;
795 continue;
796 }
797 // Check for real or null collision.
798 if (RndmUniform() <= fReal * fInv) break;
799 if (m_useNullCollisionSteps) {
800 isNullCollision = true;
801 break;
802 }
803 }
804 if (!ok) break;
805
806 // Increase the collision counter.
807 ++nCollTemp;
808
809 // Calculate the direction at the instant before the collision
810 // and the proposed new position.
811 double kx1 = 0., ky1 = 0., kz1 = 0.;
812 double dx = 0., dy = 0., dz = 0.;
813 if (m_useBfield) {
814 // Calculate the new velocity.
815 double vx1 = vx + 2. * c2 * ex * dt;
816 double vy1 = vy * cphi + vz * sphi + ezovb;
817 double vz1 = vz * cphi - vy * sphi;
818 if (omega < Small) vz1 += 2. * c2 * ez * dt;
819 // Rotate back to the global frame and normalise.
820 ToGlobal(rot, vx1, vy1, vz1, kx1, ky1, kz1);
821 const double scale = 1. / Mag(kx1, ky1, kz1);
822 kx1 *= scale;
823 ky1 *= scale;
824 kz1 *= scale;
825 // Calculate the step in coordinate space.
826 dx = vx * dt + c2 * ex * dt * dt;
827 if (omega > Small) {
828 dy = vy * a3 + vz * a4 + ezovb * dt;
829 dz = vz * a3 - vy * a4;
830 } else {
831 dy = vy * dt;
832 dz = vz * dt + c2 * ez * dt * dt;
833 }
834 // Rotate back to the global frame.
835 ToGlobal(rot, dx, dy, dz, dx, dy, dz);
836 } else if (sc) {
837 // Update the wave-vector.
838 const double cdt = dt * SpeedOfLight;
839 kx1 = kx + ex * cdt;
840 ky1 = ky + ey * cdt;
841 kz1 = kz + ez * cdt;
842 double vx1 = 0., vy1 = 0, vz1 = 0.;
843 en1 = medium->GetElectronEnergy(kx1, ky1, kz1, vx1, vy1, vz1, band);
844 dx = 0.5 * (vx + vx1) * dt;
845 dy = 0.5 * (vy + vy1) * dt;
846 dz = 0.5 * (vz + vz1) * dt;
847 } else {
848 // Update the direction.
849 const double b1 = sqrt(en / en1);
850 const double b2 = 0.5 * c1 * dt / sqrt(en1);
851 kx1 = kx * b1 + ex * b2;
852 ky1 = ky * b1 + ey * b2;
853 kz1 = kz * b1 + ez * b2;
854
855 // Calculate the step in coordinate space.
856 const double b3 = dt * dt * c2;
857 dx = vx * dt + ex * b3;
858 dy = vy * dt + ey * b3;
859 dz = vz * dt + ez * b3;
860 }
861 double x1 = x + dx;
862 double y1 = y + dy;
863 double z1 = z + dz;
864 double t1 = t + dt;
865 // Get the electric field and medium at the proposed new position.
866 m_sensor->ElectricField(x1, y1, z1, ex, ey, ez, medium, status);
867 if (!hole) {
868 ex = -ex;
869 ey = -ey;
870 ez = -ez;
871 }
872 // Check if the electron is still inside a drift medium/the drift area.
873 if (status != 0 || !m_sensor->IsInArea(x1, y1, z1)) {
874 // Try to terminate the drift line close to the boundary (endpoint
875 // outside the drift medium/drift area) using iterative bisection.
876 Terminate(x, y, z, t, x1, y1, z1, t1);
877 if (m_doSignal) {
878 const int q = hole ? 1 : -1;
879 m_sensor->AddSignal(q, t, t1, x, y, z, x1, y1, z1,
880 m_integrateWeightingField,
881 m_useWeightingPotential);
882 }
883 Update(it, x1, y1, z1, t1, en, kx1, ky1, kz1, band);
884 if (status != 0) {
885 (*it).status = StatusLeftDriftMedium;
886 if (m_debug)
887 PrintStatus(hdr, "left the drift medium", x1, y1, z1, hole);
888 } else {
889 (*it).status = StatusLeftDriftArea;
890 if (m_debug)
891 PrintStatus(hdr, "left the drift area", x1, y1, z1, hole);
892 }
893 AddToEndPoints(*it, hole);
894 ok = false;
895 break;
896 }
897
898 // Check if the electron/hole has crossed a wire.
899 double xc = x, yc = y, zc = z;
900 double rc = 0.;
901 if (m_sensor->IsWireCrossed(x, y, z, x1, y1, z1, xc, yc, zc, false,
902 rc)) {
903 const double dc = Mag(xc - x, yc - y, zc - z);
904 const double tc = t + dt * dc / Mag(dx, dy, dz);
905 // If switched on, calculated the induced signal over this step.
906 if (m_doSignal) {
907 const int q = hole ? 1 : -1;
908 m_sensor->AddSignal(q, t, tc, x, y, z, xc, yc, zc,
909 m_integrateWeightingField,
910 m_useWeightingPotential);
911 }
912 Update(it, xc, yc, zc, tc, en, kx1, ky1, kz1, band);
913 (*it).status = StatusLeftDriftMedium;
914 AddToEndPoints(*it, hole);
915 ok = false;
916 if (m_debug) PrintStatus(hdr, "hit a wire", x, y, z, hole);
917 break;
918 }
919
920 // If switched on, calculate the induced signal.
921 if (m_doSignal) {
922 const int q = hole ? 1 : -1;
923 m_sensor->AddSignal(q, t, t + dt, x, y, z, x1, y1, z1,
924 m_integrateWeightingField,
925 m_useWeightingPotential);
926 }
927
928 // Update the coordinates.
929 x = x1;
930 y = y1;
931 z = z1;
932 t = t1;
933
934 // If switched on, get the magnetic field at the new location.
935 if (m_useBfield) {
936 m_sensor->MagneticField(x, y, z, bx, by, bz, status);
937 const double scale = hole ? Tesla2Internal : -Tesla2Internal;
938 bx *= scale;
939 by *= scale;
940 bz *= scale;
941 const double bmag = Mag(bx, by, bz);
942 // Update the rotation matrix.
943 RotationMatrix(bx, by, bz, bmag, ex, ey, ez, rot);
944 omega = OmegaCyclotronOverB * bmag;
945 // Calculate the electric field in the local frame.
946 ToLocal(rot, ex, ey, ez, ex, ey, ez);
947 ezovb = bmag > Small ? ez / bmag : 0.;
948 }
949
950 if (isNullCollision) {
951 en = en1;
952 kx = kx1;
953 ky = ky1;
954 kz = kz1;
955 continue;
956 }
957
958 // Get the collision type and parameters.
959 int cstype = 0;
960 int level = 0;
961 int ndxc = 0;
962 medium->GetElectronCollision(en1, cstype, level, en, kx1, ky1, kz1,
963 secondaries, ndxc, band);
964 // If activated, histogram the distance with respect to the
965 // last collision.
966 if (m_histDistance && !m_distanceHistogramType.empty()) {
967 for (const auto& htype : m_distanceHistogramType) {
968 if (htype != cstype) continue;
969 if (m_debug) {
970 std::cout << m_className << "::TransportElectron: Collision type "
971 << cstype << ". Fill distance histogram.\n";
972 getchar();
973 }
974 switch (m_distanceOption) {
975 case 'x':
976 m_histDistance->Fill((*it).xLast - x);
977 break;
978 case 'y':
979 m_histDistance->Fill((*it).yLast - y);
980 break;
981 case 'z':
982 m_histDistance->Fill((*it).zLast - z);
983 break;
984 case 'r':
985 m_histDistance->Fill(
986 Mag((*it).xLast - x, (*it).yLast - y, (*it).zLast - z));
987 break;
988 }
989 (*it).xLast = x;
990 (*it).yLast = y;
991 (*it).zLast = z;
992 break;
993 }
994 }
995
996 if (m_userHandleCollision) {
997 m_userHandleCollision(x, y, z, t, cstype, level, medium, en1, en, kx,
998 ky, kz, kx1, ky1, kz1);
999 }
1000 switch (cstype) {
1001 // Elastic collision
1002 case ElectronCollisionTypeElastic:
1003 break;
1004 // Ionising collision
1005 case ElectronCollisionTypeIonisation:
1006 if (m_viewer && m_plotIonisations) {
1007 m_viewer->AddIonisation(x, y, z);
1008 }
1009 if (m_userHandleIonisation) {
1010 m_userHandleIonisation(x, y, z, t, cstype, level, medium);
1011 }
1012 for (const auto& secondary : secondaries) {
1013 if (secondary.first == IonProdTypeElectron) {
1014 const double esec = std::max(secondary.second, Small);
1015 if (m_histSecondary) m_histSecondary->Fill(esec);
1016 // Increment the electron counter.
1017 ++m_nElectrons;
1018 if (!aval) continue;
1019 // Add the secondary electron to the stack.
1020 if (sc) {
1021 double kxs = 0., kys = 0., kzs = 0.;
1022 int bs = -1;
1023 medium->GetElectronMomentum(esec, kxs, kys, kzs, bs);
1024 AddToStack(x, y, z, t, esec, kxs, kys, kzs, bs, false,
1025 stackNew);
1026 } else {
1027 AddToStack(x, y, z, t, esec, false, stackNew);
1028 }
1029 } else if (secondary.first == IonProdTypeHole) {
1030 const double esec = std::max(secondary.second, Small);
1031 // Increment the hole counter.
1032 ++m_nHoles;
1033 if (!aval) continue;
1034 // Add the secondary hole to the stack.
1035 if (sc) {
1036 double kxs = 0., kys = 0., kzs = 0.;
1037 int bs = -1;
1038 medium->GetElectronMomentum(esec, kxs, kys, kzs, bs);
1039 AddToStack(x, y, z, t, esec, kxs, kys, kzs, bs, true,
1040 stackNew);
1041 } else {
1042 AddToStack(x, y, z, t, esec, true, stackNew);
1043 }
1044 } else if (secondary.first == IonProdTypeIon) {
1045 ++m_nIons;
1046 }
1047 }
1048 secondaries.clear();
1049 if (m_debug) PrintStatus(hdr, "ionised", x, y, z, hole);
1050 break;
1051 // Attachment
1052 case ElectronCollisionTypeAttachment:
1053 if (m_viewer && m_plotAttachments) {
1054 m_viewer->AddAttachment(x, y, z);
1055 }
1056 if (m_userHandleAttachment) {
1057 m_userHandleAttachment(x, y, z, t, cstype, level, medium);
1058 }
1059 Update(it, x, y, z, t, en, kx1, ky1, kz1, band);
1060 (*it).status = StatusAttached;
1061 if (hole) {
1062 m_endpointsHoles.push_back(*it);
1063 --m_nHoles;
1064 } else {
1065 m_endpointsElectrons.push_back(*it);
1066 --m_nElectrons;
1067 }
1068 ok = false;
1069 break;
1070 // Inelastic collision
1071 case ElectronCollisionTypeInelastic:
1072 if (m_userHandleInelastic) {
1073 m_userHandleInelastic(x, y, z, t, cstype, level, medium);
1074 }
1075 break;
1076 // Excitation
1077 case ElectronCollisionTypeExcitation:
1078 if (m_viewer && m_plotExcitations) {
1079 m_viewer->AddExcitation(x, y, z);
1080 }
1081 if (m_userHandleInelastic) {
1082 m_userHandleInelastic(x, y, z, t, cstype, level, medium);
1083 }
1084 if (ndxc <= 0) break;
1085 // Get the electrons/photons produced in the deexcitation cascade.
1086 stackPhotons.clear();
1087 for (int j = ndxc; j--;) {
1088 double tdx = 0., sdx = 0., edx = 0.;
1089 int typedx = 0;
1090 if (!medium->GetDeexcitationProduct(j, tdx, sdx, typedx, edx)) {
1091 std::cerr << hdr << "Cannot retrieve deexcitation product " << j
1092 << "/" << ndxc << ".\n";
1093 break;
1094 }
1095
1096 if (typedx == DxcProdTypeElectron) {
1097 // Penning ionisation
1098 double xp = x, yp = y, zp = z;
1099 if (sdx > Small) {
1100 // Randomise the point of creation.
1101 double dxp = 0., dyp = 0., dzp = 0.;
1102 RndmDirection(dxp, dyp, dzp);
1103 xp += sdx * dxp;
1104 yp += sdx * dyp;
1105 zp += sdx * dzp;
1106 }
1107 // Get the electric field and medium at this location.
1108 Medium* med = nullptr;
1109 double fx = 0., fy = 0., fz = 0.;
1110 m_sensor->ElectricField(xp, yp, zp, fx, fy, fz, med, status);
1111 // Check if this location is inside a drift medium/area.
1112 if (status != 0 || !m_sensor->IsInArea(xp, yp, zp)) continue;
1113 // Make sure we haven't jumped across a wire.
1114 if (m_sensor->IsWireCrossed(x, y, z, xp, yp, zp, xc, yc, zc,
1115 false, rc)) {
1116 continue;
1117 }
1118 // Increment the electron and ion counters.
1119 ++m_nElectrons;
1120 ++m_nIons;
1121 if (m_userHandleIonisation) {
1122 m_userHandleIonisation(xp, yp, zp, t, cstype, level, medium);
1123 }
1124 if (!aval) continue;
1125 // Add the Penning electron to the list.
1126 AddToStack(xp, yp, zp, t + tdx, std::max(edx, Small), false,
1127 stackNew);
1128 } else if (typedx == DxcProdTypePhoton && m_usePhotons &&
1129 edx > m_gammaCut) {
1130 // Radiative de-excitation
1131 stackPhotons.emplace_back(std::make_pair(t + tdx, edx));
1132 }
1133 }
1134
1135 // Transport the photons (if any)
1136 if (aval) {
1137 for (const auto& ph : stackPhotons) {
1138 TransportPhoton(x, y, z, ph.first, ph.second, stackNew);
1139 }
1140 }
1141 break;
1142 // Super-elastic collision
1143 case ElectronCollisionTypeSuperelastic:
1144 break;
1145 // Virtual/null collision
1146 case ElectronCollisionTypeVirtual:
1147 break;
1148 // Acoustic phonon scattering (intravalley)
1149 case ElectronCollisionTypeAcousticPhonon:
1150 break;
1151 // Optical phonon scattering (intravalley)
1152 case ElectronCollisionTypeOpticalPhonon:
1153 break;
1154 // Intervalley scattering (phonon assisted)
1155 case ElectronCollisionTypeIntervalleyG:
1156 case ElectronCollisionTypeIntervalleyF:
1157 case ElectronCollisionTypeInterbandXL:
1158 case ElectronCollisionTypeInterbandXG:
1159 case ElectronCollisionTypeInterbandLG:
1160 break;
1161 // Coulomb scattering
1162 case ElectronCollisionTypeImpurity:
1163 break;
1164 default:
1165 std::cerr << hdr << "Unknown collision type.\n";
1166 ok = false;
1167 break;
1168 }
1169
1170 // Continue with the next electron/hole?
1171 if (!ok || nCollTemp > m_nCollSkip ||
1172 cstype == ElectronCollisionTypeIonisation ||
1173 (m_plotExcitations && cstype == ElectronCollisionTypeExcitation) ||
1174 (m_plotAttachments && cstype == ElectronCollisionTypeAttachment)) {
1175 break;
1176 }
1177 kx = kx1;
1178 ky = ky1;
1179 kz = kz1;
1180 }
1181
1182 if (!ok) continue;
1183
1184 if (!sc) {
1185 // Normalise the direction vector.
1186 const double scale = 1. / Mag(kx, ky, kz);
1187 kx *= scale;
1188 ky *= scale;
1189 kz *= scale;
1190 }
1191 // Update the stack.
1192 Update(it, x, y, z, t, en, kx, ky, kz, band);
1193 // Add a new point to the drift line (if enabled).
1194 if (m_storeDriftLines) {
1195 point newPoint;
1196 newPoint.x = x;
1197 newPoint.y = y;
1198 newPoint.z = z;
1199 newPoint.t = t;
1200 (*it).driftLine.push_back(std::move(newPoint));
1201 }
1202 }
1203 }
1204
1205 // Calculate the induced charge.
1206 if (m_doInducedCharge) {
1207 for (const auto& ep : m_endpointsElectrons) {
1208 m_sensor->AddInducedCharge(-1, ep.x0, ep.y0, ep.z0, ep.x, ep.y, ep.z);
1209 }
1210 for (const auto& ep : m_endpointsHoles) {
1211 m_sensor->AddInducedCharge(+1, ep.x0, ep.y0, ep.z0, ep.x, ep.y, ep.z);
1212 }
1213 }
1214
1215 // Plot the drift paths and photon tracks.
1216 if (m_viewer) {
1217 // Electrons
1218 const size_t nElectronEndpoints = m_endpointsElectrons.size();
1219 for (size_t i = 0; i < nElectronEndpoints; ++i) {
1220 const size_t np = GetNumberOfElectronDriftLinePoints(i);
1221 if (np <= 0) continue;
1222 size_t k;
1223 const Electron& p = m_endpointsElectrons[i];
1224 m_viewer->NewElectronDriftLine(np, k, p.x0, p.y0, p.z0);
1225 for (size_t j = 0; j < np; ++j) {
1226 double x = 0., y = 0., z = 0., t = 0.;
1227 GetElectronDriftLinePoint(x, y, z, t, j, i);
1228 m_viewer->SetDriftLinePoint(k, j, x, y, z);
1229 }
1230 }
1231 // Holes
1232 const size_t nHoleEndpoints = m_endpointsHoles.size();
1233 for (size_t i = 0; i < nHoleEndpoints; ++i) {
1234 const size_t np = GetNumberOfHoleDriftLinePoints(i);
1235 if (np <= 0) continue;
1236 size_t k;
1237 const Electron& p = m_endpointsHoles[i];
1238 m_viewer->NewHoleDriftLine(np, k, p.x0, p.y0, p.z0);
1239 for (size_t j = 0; j < np; ++j) {
1240 double x = 0., y = 0., z = 0., t = 0.;
1241 GetHoleDriftLinePoint(x, y, z, t, j, i);
1242 m_viewer->SetDriftLinePoint(k, j, x, y, z);
1243 }
1244 }
1245 // Photons
1246 for (const auto& ph : m_photons) {
1247 m_viewer->AddPhoton(ph.x0, ph.y0, ph.z0, ph.x1, ph.y1, ph.z1);
1248 }
1249 }
1250 return true;
1251}
1252
1253void AvalancheMicroscopic::TransportPhoton(const double x0, const double y0,
1254 const double z0, const double t0,
1255 const double e0,
1256 std::vector<Electron>& stack) {
1257 // Make sure that the sensor is defined.
1258 if (!m_sensor) {
1259 std::cerr << m_className << "::TransportPhoton: Sensor is not defined.\n";
1260 return;
1261 }
1262
1263 // Make sure that the starting point is inside a medium.
1264 Medium* medium;
1265 if (!m_sensor->GetMedium(x0, y0, z0, medium)) {
1266 std::cerr << m_className << "::TransportPhoton:\n"
1267 << " No medium at initial position.\n";
1268 return;
1269 }
1270
1271 // Make sure that the medium is "driftable" and microscopic.
1272 if (!medium->IsDriftable() || !medium->IsMicroscopic()) {
1273 std::cerr << m_className << "::TransportPhoton:\n"
1274 << " Medium at initial position does not provide "
1275 << " microscopic tracking data.\n";
1276 return;
1277 }
1278
1279 // Get the id number of the drift medium.
1280 int id = medium->GetId();
1281
1282 // Position
1283 double x = x0, y = y0, z = z0;
1284 double t = t0;
1285 // Initial direction (randomised).
1286 double dx = 0., dy = 0., dz = 0.;
1287 RndmDirection(dx, dy, dz);
1288 // Energy
1289 double e = e0;
1290
1291 // Photon collision rate
1292 double f = medium->GetPhotonCollisionRate(e);
1293 if (f <= 0.) return;
1294 // Timestep
1295 double dt = -log(RndmUniformPos()) / f;
1296 t += dt;
1297 dt *= SpeedOfLight;
1298 x += dt * dx;
1299 y += dt * dy;
1300 z += dt * dz;
1301
1302 // Check if the photon is still inside a medium.
1303 if (!m_sensor->GetMedium(x, y, z, medium) || medium->GetId() != id) {
1304 // Try to terminate the photon track close to the boundary
1305 // by means of iterative bisection.
1306 dx *= dt;
1307 dy *= dt;
1308 dz *= dt;
1309 x -= dx;
1310 y -= dy;
1311 z -= dz;
1312 double delta = Mag(dx, dy, dz);
1313 if (delta > 0) {
1314 dx /= delta;
1315 dy /= delta;
1316 dz /= delta;
1317 }
1318 // Mid-point
1319 double xM = x, yM = y, zM = z;
1320 while (delta > BoundaryDistance) {
1321 delta *= 0.5;
1322 dt *= 0.5;
1323 xM = x + delta * dx;
1324 yM = y + delta * dy;
1325 zM = z + delta * dz;
1326 // Check if the mid-point is inside the drift medium.
1327 if (m_sensor->GetMedium(xM, yM, zM, medium) && medium->GetId() == id) {
1328 x = xM;
1329 y = yM;
1330 z = zM;
1331 t += dt;
1332 }
1333 }
1334 photon newPhoton;
1335 newPhoton.x0 = x0;
1336 newPhoton.y0 = y0;
1337 newPhoton.z0 = z0;
1338 newPhoton.x1 = x;
1339 newPhoton.y1 = y;
1340 newPhoton.z1 = z;
1341 newPhoton.energy = e0;
1342 newPhoton.status = StatusLeftDriftMedium;
1343 m_photons.push_back(std::move(newPhoton));
1344 return;
1345 }
1346
1347 int type, level;
1348 double e1;
1349 double ctheta = 0.;
1350 int nsec = 0;
1351 double esec = 0.;
1352 if (!medium->GetPhotonCollision(e, type, level, e1, ctheta, nsec, esec))
1353 return;
1354
1355 if (type == PhotonCollisionTypeIonisation) {
1356 // Add the secondary electron (random direction) to the stack.
1357 if (m_sizeCut == 0 || stack.size() < m_sizeCut) {
1358 AddToStack(x, y, z, t, std::max(esec, Small), false, stack);
1359 }
1360 // Increment the electron and ion counters.
1361 ++m_nElectrons;
1362 ++m_nIons;
1363 } else if (type == PhotonCollisionTypeExcitation) {
1364 double tdx = 0.;
1365 double sdx = 0.;
1366 int typedx = 0;
1367 std::vector<double> tPhotons;
1368 std::vector<double> ePhotons;
1369 for (int j = nsec; j--;) {
1370 if (!medium->GetDeexcitationProduct(j, tdx, sdx, typedx, esec)) continue;
1371 if (typedx == DxcProdTypeElectron) {
1372 // Ionisation.
1373 AddToStack(x, y, z, t + tdx, std::max(esec, Small), false, stack);
1374 // Increment the electron and ion counters.
1375 ++m_nElectrons;
1376 ++m_nIons;
1377 } else if (typedx == DxcProdTypePhoton && m_usePhotons &&
1378 esec > m_gammaCut) {
1379 // Radiative de-excitation
1380 tPhotons.push_back(t + tdx);
1381 ePhotons.push_back(esec);
1382 }
1383 }
1384 // Transport the photons (if any).
1385 const int nSizePhotons = tPhotons.size();
1386 for (int k = nSizePhotons; k--;) {
1387 TransportPhoton(x, y, z, tPhotons[k], ePhotons[k], stack);
1388 }
1389 }
1390
1391 photon newPhoton;
1392 newPhoton.x0 = x0;
1393 newPhoton.y0 = y0;
1394 newPhoton.z0 = z0;
1395 newPhoton.x1 = x;
1396 newPhoton.y1 = y;
1397 newPhoton.z1 = z;
1398 newPhoton.energy = e0;
1399 newPhoton.status = -2;
1400 m_photons.push_back(std::move(newPhoton));
1401}
1402
1403void AvalancheMicroscopic::Update(std::vector<Electron>::iterator it,
1404 const double x, const double y,
1405 const double z, const double t,
1406 const double energy, const double kx,
1407 const double ky, const double kz,
1408 const int band) {
1409 (*it).x = x;
1410 (*it).y = y;
1411 (*it).z = z;
1412 (*it).t = t;
1413 (*it).energy = energy;
1414 (*it).kx = kx;
1415 (*it).ky = ky;
1416 (*it).kz = kz;
1417 (*it).band = band;
1418}
1419
1420void AvalancheMicroscopic::AddToStack(const double x, const double y,
1421 const double z, const double t,
1422 const double energy, const bool hole,
1423 std::vector<Electron>& container) const {
1424 // Randomise the direction.
1425 double dx = 0., dy = 0., dz = 1.;
1426 RndmDirection(dx, dy, dz);
1427 AddToStack(x, y, z, t, energy, dx, dy, dz, 0, hole, container);
1428}
1429
1430void AvalancheMicroscopic::AddToStack(const double x, const double y,
1431 const double z, const double t,
1432 const double energy, const double dx,
1433 const double dy, const double dz,
1434 const int band, const bool hole,
1435 std::vector<Electron>& container) const {
1436 Electron electron;
1437 electron.status = 0;
1438 electron.hole = hole;
1439 electron.x0 = x;
1440 electron.y0 = y;
1441 electron.z0 = z;
1442 electron.t0 = t;
1443 electron.e0 = energy;
1444 electron.x = x;
1445 electron.y = y;
1446 electron.z = z;
1447 electron.t = t;
1448 electron.energy = energy;
1449 electron.kx = dx;
1450 electron.ky = dy;
1451 electron.kz = dz;
1452 electron.band = band;
1453 // Previous coordinates for distance histogramming.
1454 electron.xLast = x;
1455 electron.yLast = y;
1456 electron.zLast = z;
1457 electron.driftLine.reserve(1000);
1458 container.push_back(std::move(electron));
1459}
1460
1461void AvalancheMicroscopic::Terminate(double x0, double y0, double z0, double t0,
1462 double& x1, double& y1, double& z1,
1463 double& t1) {
1464 const double dx = x1 - x0;
1465 const double dy = y1 - y0;
1466 const double dz = z1 - z0;
1467 double d = Mag(dx, dy, dz);
1468 while (d > BoundaryDistance) {
1469 d *= 0.5;
1470 const double xm = 0.5 * (x0 + x1);
1471 const double ym = 0.5 * (y0 + y1);
1472 const double zm = 0.5 * (z0 + z1);
1473 const double tm = 0.5 * (t0 + t1);
1474 // Check if the mid-point is inside the drift medium.
1475 double ex = 0., ey = 0., ez = 0.;
1476 Medium* medium = nullptr;
1477 int status = 0;
1478 m_sensor->ElectricField(xm, ym, zm, ex, ey, ez, medium, status);
1479 if (status == 0 && m_sensor->IsInArea(xm, ym, zm)) {
1480 x0 = xm;
1481 y0 = ym;
1482 z0 = zm;
1483 t0 = tm;
1484 } else {
1485 x1 = xm;
1486 y1 = ym;
1487 z1 = zm;
1488 t1 = tm;
1489 }
1490 }
1491}
1492
1493} // namespace Garfield
void EnableDistanceHistogramming(const int type)
Fill distance distribution histograms for a given collision type.
void SetDistanceHistogram(TH1 *histo, const char opt='r')
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 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))
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.)
unsigned int GetNumberOfElectronDriftLinePoints(const unsigned int i=0) const
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 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 GetElectronDriftLinePoint(double &x, double &y, double &z, double &t, const int ip, const unsigned int iel=0) const
void EnableElectronEnergyHistogramming(TH1 *histo)
Fill a histogram with the electron energy distribution.
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 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.
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 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
void MagneticField(const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
Get the magnetic field at (x, y, z).
Definition: Sensor.cc:116
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&medium, int &status)
Get the drift field and potential at (x, y, z).
Definition: Sensor.cc:65
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:254
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:159
void AddSignal(const double q, const double t0, const double t1, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const bool integrateWeightingField, const bool useWeightingPotential=false)
Add the signal from a charge-carrier step.
Definition: Sensor.cc:452
bool IsWireCrossed(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc, const bool centre, double &rc)
Definition: Sensor.cc:276
void AddInducedCharge(const double q, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1)
Add the induced charge from a charge carrier drift between two points.
Definition: Sensor.cc:798
Visualize drift lines and tracks.
Definition: ViewDrift.hh:18
void AddIonisation(const float x, const float y, const float z)
Definition: ViewDrift.cc:173
void NewHoleDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition: ViewDrift.cc:81
void SetDriftLinePoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
Definition: ViewDrift.cc:124
void AddAttachment(const float x, const float y, const float z)
Definition: ViewDrift.cc:179
void AddPhoton(const float x0, const float y0, const float z0, const float x1, const float y1, const float z1)
Definition: ViewDrift.cc:103
void NewElectronDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
Definition: ViewDrift.cc:68
void AddExcitation(const float x, const float y, const float z)
Definition: ViewDrift.cc:167
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:107
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition: Random.hh:17
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314