13double Mag(
const double x,
const double y,
const double z) {
14 return sqrt(x * x + y * y + z * z);
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,
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;
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,
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;
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) {
41 std::array<std::array<double, 3>, 3> rB = {{{1, 0, 0}, {0, 1, 0}, {0, 0, 1}}};
42 if (bmag > Garfield::Small) {
46 const double bt = by * by + bz * bz;
47 if (bt > Garfield::Small) {
48 const double btInv = 1. / bt;
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;
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) {
66 rX[1][1] = rX[2][2] = fz / ft;
70 for (
unsigned int i = 0; i < 3; ++i) {
71 for (
unsigned int j = 0; j < 3; ++j) {
73 for (
unsigned int k = 0; k < 3; ++k) {
74 rot[i][j] += rX[i][k] * rB[k][j];
80void PrintStatus(
const std::string& hdr,
const std::string& status,
81 const double x,
const double y,
const double z,
83 const std::string eh = hole ?
"Hole " :
"Electron ";
84 std::cout << hdr << eh << status <<
" at " <<
x <<
", " <<
y <<
", " <<
z
92 m_endpointsElectrons.reserve(10000);
93 m_endpointsHoles.reserve(10000);
94 m_photons.reserve(1000);
99 std::cerr << m_className <<
"::SetSensor: Null pointer.\n";
107 std::cerr << m_className <<
"::EnablePlotting: Null pointer.\n";
112 if (!m_storeDriftLines) {
113 std::cout << m_className <<
"::EnablePlotting:\n"
114 <<
" Enabling storage of drift line.\n";
121 std::cerr << m_className <<
"::EnableElectronEnergyHistogramming:\n"
122 <<
" Null pointer.\n";
126 m_histElectronEnergy = histo;
131 std::cerr << m_className <<
"::EnableHoleEnergyHistogramming:\n"
132 <<
" Null pointer.\n";
136 m_histHoleEnergy = histo;
141 std::cerr << m_className <<
"::SetDistanceHistogram: Null pointer.\n";
145 m_histDistance = histo;
147 if (opt ==
'x' || opt ==
'y' || opt ==
'z' || opt ==
'r') {
148 m_distanceOption = opt;
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';
157 if (m_distanceHistogramType.empty()) {
158 std::cout << m_className <<
"::SetDistanceHistogram:\n";
159 std::cout <<
" Don't forget to call EnableDistanceHistogramming.\n";
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";
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";
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";
193 m_distanceHistogramType.erase(
194 std::remove(m_distanceHistogramType.begin(),
195 m_distanceHistogramType.end(), type),
196 m_distanceHistogramType.end());
200 m_histDistance =
nullptr;
201 m_distanceHistogramType.clear();
206 std::cerr << m_className <<
"::EnableSecondaryEnergyHistogramming:\n"
207 <<
" Null pointer.\n";
211 m_histSecondary = histo;
215 if (fabs(t1 - t0) < Small) {
216 std::cerr << m_className <<
"::SetTimeWindow:\n";
217 std::cerr <<
" Time interval must be greater than zero.\n";
221 m_tMin = std::min(t0, t1);
222 m_tMax = std::max(t0, t1);
223 m_hasTimeWindow =
true;
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.;
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;
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.;
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;
283 double& y0,
double& z0,
double& t0,
284 double& e0,
double& x1,
double& y1,
285 double& z1,
double& t1,
double& e1,
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.;
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;
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";
316 if (!m_storeDriftLines)
return 2;
318 return m_endpointsElectrons[i].driftLine.size() + 2;
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";
329 if (!m_storeDriftLines)
return 2;
331 return m_endpointsHoles[i].driftLine.size() + 2;
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";
344 x = m_endpointsElectrons[iel].x0;
345 y = m_endpointsElectrons[iel].y0;
346 z = m_endpointsElectrons[iel].z0;
347 t = m_endpointsElectrons[iel].t0;
351 const int np = m_endpointsElectrons[iel].driftLine.size();
353 x = m_endpointsElectrons[iel].x;
354 y = m_endpointsElectrons[iel].y;
355 z = m_endpointsElectrons[iel].z;
356 t = m_endpointsElectrons[iel].t;
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;
367 double& z,
double& t,
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";
377 x = m_endpointsHoles[ih].x0;
378 y = m_endpointsHoles[ih].y0;
379 z = m_endpointsHoles[ih].z0;
380 t = m_endpointsHoles[ih].t0;
384 const int np = m_endpointsHoles[ih].driftLine.size();
386 x = m_endpointsHoles[ih].x;
387 y = m_endpointsHoles[ih].y;
388 z = m_endpointsHoles[ih].z;
389 t = m_endpointsHoles[ih].t;
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;
400 double& x0,
double& y0,
double& z0,
401 double& t0,
double& x1,
double& y1,
402 double& z1,
double& t1,
404 if (i >= m_photons.size()) {
405 std::cerr << m_className <<
"::GetPhoton: Index out of range.\n";
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;
422 void (*f)(
double x,
double y,
double z,
double t,
double e,
double dx,
423 double dy,
double dz,
bool hole)) {
425 std::cerr << m_className <<
"::SetUserHandleStep: Null pointer.\n";
428 m_userHandleStep = f;
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;
439 double x,
double y,
double z,
double t,
int type,
int level,
Medium* m)) {
440 m_userHandleAttachment = f;
444 double x,
double y,
double z,
double t,
int type,
int level,
Medium* m)) {
445 m_userHandleInelastic = f;
449 double x,
double y,
double z,
double t,
int type,
int level,
Medium* m)) {
450 m_userHandleIonisation = f;
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);
463 const double z0,
const double t0,
464 const double e0,
const double dx0,
467 std::vector<Electron> stack;
468 AddToStack(x0, y0, z0, t0, e0, dx0, dy0, dz0, 0,
false, stack);
469 return TransportElectrons(stack,
true);
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);
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);
486 return TransportElectrons(stack,
true);
489bool AvalancheMicroscopic::TransportElectrons(std::vector<Electron>& stack,
493 m_endpointsElectrons.clear();
494 m_endpointsHoles.clear();
498 m_nElectrons = m_nHoles = m_nIons = 0;
500 const std::string hdr = m_className +
"::TransportElectron: ";
503 std::cerr << hdr <<
"Sensor is not defined.\n";
508 for (
auto& p : stack) {
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";
516 if (!medium->IsDriftable() || !medium->IsMicroscopic()) {
517 std::cerr << hdr <<
"Medium does not have cross-section data.\n";
521 p.e0 = std::max(p.e0, Small);
522 if (medium->IsSemiconductor() && m_useBandStructure) {
525 medium->GetElectronMomentum(p.e0, p.kx, p.ky, p.kz, p.band);
529 const double kmag = Mag(p.kx, p.ky, p.kz);
530 if (
fabs(kmag) < Small) {
548 const double c1 = SpeedOfLight *
sqrt(2. / ElectronMass);
549 const double c2 = c1 * c1 / 4.;
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;
558 stack.erase(std::remove_if(stack.begin(), stack.end(), IsInactive),
561 if (aval && m_sizeCut > 0) {
563 if (stack.size() > m_sizeCut) {
565 }
else if (stack.size() + stackNew.size() > m_sizeCut) {
566 stackNew.resize(m_sizeCut - stack.size());
569 stack.insert(stack.end(), std::make_move_iterator(stackNew.begin()),
570 std::make_move_iterator(stackNew.end()));
573 if (stack.empty())
break;
575 for (
auto it = stack.begin(), end = stack.end(); it != end; ++it) {
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;
591 unsigned int nCollTemp = 0;
594 double ex = 0., ey = 0., ez = 0.;
595 Medium* medium =
nullptr;
597 m_sensor->
ElectricField(x, y, z, ex, ey, ez, medium, status);
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";
613 if (status != 0 || !medium || !medium->IsDriftable() ||
614 !medium->IsMicroscopic()) {
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);
623 auto id = medium->GetId();
624 bool sc = (medium->IsSemiconductor() && m_useBandStructure);
626 double fLim = medium->GetElectronNullCollisionRate(band);
628 std::cerr << hdr <<
"Got null-collision rate <= 0.\n";
631 double fInv = 1. / fLim;
634 double bx = 0., by = 0., bz = 0.;
639 std::array<std::array<double, 3>, 3> rot;
642 const double scale = hole ? Tesla2Internal : -Tesla2Internal;
646 const double bmag = Mag(bx, by, bz);
649 RotationMatrix(bx, by, bz, bmag, ex, ey, ez, rot);
651 omega = OmegaCyclotronOverB *
bmag;
653 ToLocal(rot, ex, ey, ez, ex, ey, ez);
654 ezovb =
bmag > Small ? ez /
bmag : 0.;
659 bool isNullCollision =
false;
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);
667 std::cout << hdr <<
"Kinetic energy (" << en
668 <<
") below transport cut.\n";
675 if (hole && m_histHoleEnergy) {
676 m_histHoleEnergy->Fill(en);
677 }
else if (!hole && m_histElectronEnergy) {
678 m_histElectronEnergy->Fill(en);
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);
691 if (medium->GetId() !=
id) {
693 if (!medium->IsMicroscopic()) {
695 Update(it, x, y, z, t, en, kx, ky, kz, band);
696 (*it).status = StatusLeftDriftMedium;
697 AddToEndPoints(*it, hole);
700 std::cout << hdr <<
"\n Medium at " <<
x <<
", " <<
y <<
", "
701 <<
z <<
" does not have cross-section data.\n";
705 id = medium->GetId();
706 sc = (medium->IsSemiconductor() && m_useBandStructure);
708 fLim = medium->GetElectronNullCollisionRate(band);
710 std::cerr << hdr <<
"Got null-collision rate <= 0.\n";
716 double a1 = 0., a2 = 0.;
718 double vx = 0., vy = 0., vz = 0.;
721 const double vmag = c1 *
sqrt(en);
722 ToLocal(rot, vmag * kx, vmag * ky, vmag * kz, vx, vy, vz);
732 en = medium->GetElectronEnergy(kx, ky, kz, vx, vy, vz, band);
736 const double vmag = c1 *
sqrt(en);
740 a1 = vx * ex + vy * ey + vz * ez;
741 a2 = c2 * (ex * ex + ey * ey + ez * ez);
744 if (m_userHandleStep) {
745 m_userHandleStep(x, y, z, t, en, kx, ky, kz, hole);
753 double cphi = 1., sphi = 0.;
754 double a3 = 0., a4 = 0.;
758 dt += -log(r) * fInv;
761 en1 = en + (a1 + a2 * dt) * dt;
763 cphi =
cos(omega * dt);
764 sphi =
sin(omega * dt);
766 a4 = (1. - cphi) / omega;
767 en1 += ez * (vz * a3 - vy * a4);
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);
777 en1 = en + (a1 + a2 * dt) * dt;
779 en1 = std::max(en1, Small);
781 double fReal = medium->GetElectronCollisionRate(en1, band);
783 std::cerr << hdr <<
"Got collision rate <= 0 at " << en1
784 <<
" eV (band " << band <<
").\n";
791 std::cerr << hdr <<
"Increasing null-collision rate by 5%.\n";
792 if (sc) std::cerr <<
" Band " << band <<
"\n";
799 if (m_useNullCollisionSteps) {
800 isNullCollision =
true;
811 double kx1 = 0., ky1 = 0., kz1 = 0.;
812 double dx = 0., dy = 0., dz = 0.;
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;
820 ToGlobal(rot, vx1, vy1, vz1, kx1, ky1, kz1);
821 const double scale = 1. / Mag(kx1, ky1, kz1);
826 dx = vx * dt + c2 * ex * dt * dt;
828 dy = vy * a3 + vz * a4 + ezovb * dt;
829 dz = vz * a3 - vy * a4;
832 dz = vz * dt + c2 * ez * dt * dt;
835 ToGlobal(rot, dx, dy, dz, dx, dy, dz);
838 const double cdt = dt * SpeedOfLight;
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;
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;
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;
866 m_sensor->
ElectricField(x1, y1, z1, ex, ey, ez, medium, status);
873 if (status != 0 || !m_sensor->
IsInArea(x1, y1, z1)) {
876 Terminate(x, y, z, t, x1, y1, z1, t1);
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);
883 Update(it, x1, y1, z1, t1, en, kx1, ky1, kz1, band);
885 (*it).status = StatusLeftDriftMedium;
887 PrintStatus(hdr,
"left the drift medium", x1, y1, z1, hole);
889 (*it).status = StatusLeftDriftArea;
891 PrintStatus(hdr,
"left the drift area", x1, y1, z1, hole);
893 AddToEndPoints(*it, hole);
899 double xc =
x, yc =
y, zc =
z;
901 if (m_sensor->
IsWireCrossed(x, y, z, x1, y1, z1, xc, yc, zc,
false,
903 const double dc = Mag(xc - x, yc - y, zc - z);
904 const double tc = t + dt * dc / Mag(dx, dy, dz);
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);
912 Update(it, xc, yc, zc, tc, en, kx1, ky1, kz1, band);
913 (*it).status = StatusLeftDriftMedium;
914 AddToEndPoints(*it, hole);
916 if (m_debug) PrintStatus(hdr,
"hit a wire", x, y, z, hole);
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);
937 const double scale = hole ? Tesla2Internal : -Tesla2Internal;
941 const double bmag = Mag(bx, by, bz);
943 RotationMatrix(bx, by, bz, bmag, ex, ey, ez, rot);
944 omega = OmegaCyclotronOverB *
bmag;
946 ToLocal(rot, ex, ey, ez, ex, ey, ez);
947 ezovb =
bmag > Small ? ez /
bmag : 0.;
950 if (isNullCollision) {
962 medium->GetElectronCollision(en1, cstype, level, en, kx1, ky1, kz1,
963 secondaries, ndxc, band);
966 if (m_histDistance && !m_distanceHistogramType.empty()) {
967 for (
const auto& htype : m_distanceHistogramType) {
968 if (htype != cstype)
continue;
970 std::cout << m_className <<
"::TransportElectron: Collision type "
971 << cstype <<
". Fill distance histogram.\n";
974 switch (m_distanceOption) {
976 m_histDistance->Fill((*it).xLast - x);
979 m_histDistance->Fill((*it).yLast - y);
982 m_histDistance->Fill((*it).zLast - z);
985 m_histDistance->Fill(
986 Mag((*it).xLast - x, (*it).yLast - y, (*it).zLast - z));
996 if (m_userHandleCollision) {
997 m_userHandleCollision(x, y, z, t, cstype, level, medium, en1, en, kx,
998 ky, kz, kx1, ky1, kz1);
1002 case ElectronCollisionTypeElastic:
1005 case ElectronCollisionTypeIonisation:
1006 if (m_viewer && m_plotIonisations) {
1009 if (m_userHandleIonisation) {
1010 m_userHandleIonisation(x, y, z, t, cstype, level, medium);
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);
1018 if (!aval)
continue;
1021 double kxs = 0., kys = 0., kzs = 0.;
1023 medium->GetElectronMomentum(esec, kxs, kys, kzs, bs);
1024 AddToStack(x, y, z, t, esec, kxs, kys, kzs, bs,
false,
1027 AddToStack(x, y, z, t, esec,
false, stackNew);
1029 }
else if (secondary.first == IonProdTypeHole) {
1030 const double esec = std::max(secondary.second, Small);
1033 if (!aval)
continue;
1036 double kxs = 0., kys = 0., kzs = 0.;
1038 medium->GetElectronMomentum(esec, kxs, kys, kzs, bs);
1039 AddToStack(x, y, z, t, esec, kxs, kys, kzs, bs,
true,
1042 AddToStack(x, y, z, t, esec,
true, stackNew);
1044 }
else if (secondary.first == IonProdTypeIon) {
1048 secondaries.clear();
1049 if (m_debug) PrintStatus(hdr,
"ionised", x, y, z, hole);
1052 case ElectronCollisionTypeAttachment:
1053 if (m_viewer && m_plotAttachments) {
1056 if (m_userHandleAttachment) {
1057 m_userHandleAttachment(x, y, z, t, cstype, level, medium);
1059 Update(it, x, y, z, t, en, kx1, ky1, kz1, band);
1060 (*it).status = StatusAttached;
1062 m_endpointsHoles.push_back(*it);
1065 m_endpointsElectrons.push_back(*it);
1071 case ElectronCollisionTypeInelastic:
1072 if (m_userHandleInelastic) {
1073 m_userHandleInelastic(x, y, z, t, cstype, level, medium);
1077 case ElectronCollisionTypeExcitation:
1078 if (m_viewer && m_plotExcitations) {
1081 if (m_userHandleInelastic) {
1082 m_userHandleInelastic(x, y, z, t, cstype, level, medium);
1084 if (ndxc <= 0)
break;
1086 stackPhotons.clear();
1087 for (
int j = ndxc; j--;) {
1088 double tdx = 0., sdx = 0., edx = 0.;
1090 if (!medium->GetDeexcitationProduct(j, tdx, sdx, typedx, edx)) {
1091 std::cerr << hdr <<
"Cannot retrieve deexcitation product " << j
1092 <<
"/" << ndxc <<
".\n";
1096 if (typedx == DxcProdTypeElectron) {
1098 double xp =
x, yp =
y, zp =
z;
1101 double dxp = 0., dyp = 0., dzp = 0.;
1108 Medium* med =
nullptr;
1109 double fx = 0., fy = 0., fz = 0.;
1110 m_sensor->
ElectricField(xp, yp, zp, fx, fy, fz, med, status);
1112 if (status != 0 || !m_sensor->
IsInArea(xp, yp, zp))
continue;
1114 if (m_sensor->
IsWireCrossed(x, y, z, xp, yp, zp, xc, yc, zc,
1121 if (m_userHandleIonisation) {
1122 m_userHandleIonisation(xp, yp, zp, t, cstype, level, medium);
1124 if (!aval)
continue;
1126 AddToStack(xp, yp, zp, t + tdx, std::max(edx, Small),
false,
1128 }
else if (typedx == DxcProdTypePhoton && m_usePhotons &&
1131 stackPhotons.emplace_back(std::make_pair(t + tdx, edx));
1137 for (
const auto& ph : stackPhotons) {
1138 TransportPhoton(x, y, z, ph.first, ph.second, stackNew);
1143 case ElectronCollisionTypeSuperelastic:
1146 case ElectronCollisionTypeVirtual:
1149 case ElectronCollisionTypeAcousticPhonon:
1152 case ElectronCollisionTypeOpticalPhonon:
1155 case ElectronCollisionTypeIntervalleyG:
1156 case ElectronCollisionTypeIntervalleyF:
1157 case ElectronCollisionTypeInterbandXL:
1158 case ElectronCollisionTypeInterbandXG:
1159 case ElectronCollisionTypeInterbandLG:
1162 case ElectronCollisionTypeImpurity:
1165 std::cerr << hdr <<
"Unknown collision type.\n";
1171 if (!ok || nCollTemp > m_nCollSkip ||
1172 cstype == ElectronCollisionTypeIonisation ||
1173 (m_plotExcitations && cstype == ElectronCollisionTypeExcitation) ||
1174 (m_plotAttachments && cstype == ElectronCollisionTypeAttachment)) {
1186 const double scale = 1. / Mag(kx, ky, kz);
1192 Update(it, x, y, z, t, en, kx, ky, kz, band);
1194 if (m_storeDriftLines) {
1200 (*it).driftLine.push_back(std::move(newPoint));
1206 if (m_doInducedCharge) {
1207 for (
const auto& ep : m_endpointsElectrons) {
1210 for (
const auto& ep : m_endpointsHoles) {
1218 const size_t nElectronEndpoints = m_endpointsElectrons.size();
1219 for (
size_t i = 0; i < nElectronEndpoints; ++i) {
1221 if (np <= 0)
continue;
1223 const Electron& p = m_endpointsElectrons[i];
1225 for (
size_t j = 0; j < np; ++j) {
1226 double x = 0.,
y = 0.,
z = 0., t = 0.;
1232 const size_t nHoleEndpoints = m_endpointsHoles.size();
1233 for (
size_t i = 0; i < nHoleEndpoints; ++i) {
1235 if (np <= 0)
continue;
1237 const Electron& p = m_endpointsHoles[i];
1239 for (
size_t j = 0; j < np; ++j) {
1240 double x = 0.,
y = 0.,
z = 0., t = 0.;
1246 for (
const auto& ph : m_photons) {
1247 m_viewer->
AddPhoton(ph.x0, ph.y0, ph.z0, ph.x1, ph.y1, ph.z1);
1253void AvalancheMicroscopic::TransportPhoton(
const double x0,
const double y0,
1254 const double z0,
const double t0,
1256 std::vector<Electron>& stack) {
1259 std::cerr << m_className <<
"::TransportPhoton: Sensor is not defined.\n";
1265 if (!m_sensor->
GetMedium(x0, y0, z0, medium)) {
1266 std::cerr << m_className <<
"::TransportPhoton:\n"
1267 <<
" No medium at initial position.\n";
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";
1280 int id = medium->GetId();
1283 double x = x0,
y = y0,
z = z0;
1286 double dx = 0., dy = 0., dz = 0.;
1292 double f = medium->GetPhotonCollisionRate(e);
1293 if (f <= 0.)
return;
1303 if (!m_sensor->
GetMedium(x, y, z, medium) || medium->GetId() !=
id) {
1312 double delta = Mag(dx, dy, dz);
1319 double xM =
x, yM =
y, zM =
z;
1320 while (delta > BoundaryDistance) {
1323 xM =
x + delta * dx;
1324 yM =
y + delta * dy;
1325 zM =
z + delta * dz;
1327 if (m_sensor->
GetMedium(xM, yM, zM, medium) && medium->GetId() ==
id) {
1341 newPhoton.energy = e0;
1342 newPhoton.status = StatusLeftDriftMedium;
1343 m_photons.push_back(std::move(newPhoton));
1352 if (!medium->GetPhotonCollision(e, type, level, e1, ctheta, nsec, esec))
1355 if (type == PhotonCollisionTypeIonisation) {
1357 if (m_sizeCut == 0 || stack.size() < m_sizeCut) {
1358 AddToStack(x, y, z, t, std::max(esec, Small),
false, stack);
1363 }
else if (type == PhotonCollisionTypeExcitation) {
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) {
1373 AddToStack(x, y, z, t + tdx, std::max(esec, Small),
false, stack);
1377 }
else if (typedx == DxcProdTypePhoton && m_usePhotons &&
1378 esec > m_gammaCut) {
1380 tPhotons.push_back(t + tdx);
1381 ePhotons.push_back(esec);
1385 const int nSizePhotons = tPhotons.size();
1386 for (
int k = nSizePhotons; k--;) {
1387 TransportPhoton(x, y, z, tPhotons[k], ePhotons[k], stack);
1398 newPhoton.energy = e0;
1399 newPhoton.status = -2;
1400 m_photons.push_back(std::move(newPhoton));
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,
1413 (*it).energy = energy;
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 {
1425 double dx = 0., dy = 0., dz = 1.;
1427 AddToStack(x, y, z, t, energy, dx, dy, dz, 0, hole, container);
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 {
1437 electron.status = 0;
1438 electron.hole = hole;
1443 electron.e0 = energy;
1448 electron.energy = energy;
1452 electron.band = band;
1457 electron.driftLine.reserve(1000);
1458 container.push_back(std::move(electron));
1461void AvalancheMicroscopic::Terminate(
double x0,
double y0,
double z0,
double t0,
1462 double& x1,
double& y1,
double& z1,
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) {
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);
1475 double ex = 0., ey = 0., ez = 0.;
1476 Medium* medium =
nullptr;
1478 m_sensor->
ElectricField(xm, ym, zm, ex, ey, ez, medium, status);
1479 if (status == 0 && m_sensor->
IsInArea(xm, ym, zm)) {
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
AvalancheMicroscopic()
Constructor.
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.
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).
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).
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
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.
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)
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.
Visualize drift lines and tracks.
void AddIonisation(const float x, const float y, const float z)
void NewHoleDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
void SetDriftLinePoint(const size_t iL, const size_t iP, const float x, const float y, const float z)
void AddAttachment(const float x, const float y, const float z)
void AddPhoton(const float x0, const float y0, const float z0, const float x1, const float y1, const float z1)
void NewElectronDriftLine(const size_t np, size_t &id, const float x0, const float y0, const float z0)
void AddExcitation(const float x, const float y, const float z)
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
DoubleAc cos(const DoubleAc &f)
DoubleAc fabs(const DoubleAc &f)
DoubleAc sin(const DoubleAc &f)
DoubleAc sqrt(const DoubleAc &f)