Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
DriftLineRKF.cc
Go to the documentation of this file.
1#include <iostream>
2#include <cstdio>
3#include <cmath>
4#include <numeric>
5
9#include "Garfield/Random.hh"
10
11namespace {
12
13std::string PrintVec(const std::array<double, 3>& x) {
14
15 return "(" + std::to_string(x[0]) + ", " + std::to_string(x[1]) + ", " +
16 std::to_string(x[2]) + ")";
17}
18
19double Mag(const std::array<double, 3>& x) {
20
21 return sqrt(x[0] * x[0] + x[1] * x[1] + x[2] * x[2]);
22}
23
24double Mag(const double x, const double y, const double z) {
25
26 return sqrt(x * x + y * y + z * z);
27}
28
29double Mag(const double x, const double y) {
30
31 return sqrt(x * x + y * y);
32}
33
34double Mag2(const double x, const double y) {
35
36 return x * x + y * y;
37}
38
39}
40
41namespace Garfield {
42
43using Vec = std::array<double, 3>;
44
46 m_t.reserve(1000);
47 m_x.reserve(1000);
48}
49
51 if (!s) {
52 std::cerr << m_className << "::SetSensor: Null pointer.\n";
53 return;
54 }
55 m_sensor = s;
56}
57
59 if (eps > 0.) {
60 m_accuracy = eps;
61 } else {
62 std::cerr << m_className << "::SetIntegrationAccuracy:\n"
63 << " Accuracy must be greater than zero.\n";
64 }
65}
66
67void DriftLineRKF::SetMaximumStepSize(const double ms) {
68 if (ms > 0.) {
69 m_maxStepSize = ms;
70 m_useStepSizeLimit = true;
71 } else {
72 std::cerr << m_className << "::SetMaximumStepSize:\n"
73 << " Step size must be greater than zero.\n";
74 }
75}
76
78 if (!view) {
79 std::cerr << m_className << "::EnablePlotting: Null pointer.\n";
80 return;
81 }
82 m_view = view;
83}
84
85void DriftLineRKF::DisablePlotting() { m_view = nullptr; }
86
88
89 if (gain > 1.) {
90 std::cout << m_className << "::SetGainFluctuationsFixed: "
91 << "Avalanche size set to " << gain << ".\n";
92 } else {
93 std::cout << m_className << "::SetGainFluctuationsFixed:\n "
94 << "Avalanche size will be given by "
95 << "the integrated Townsend coefficient.\n";
96 }
97 m_gain = gain;
98 m_gainFluctuations = GainFluctuations::None;
99}
100
102 const double mean) {
103
104 if (theta < 0.) {
105 std::cerr << m_className << "::SetGainFluctuationsPolya: "
106 << "Shape parameter must be >= 0.\n";
107 return;
108 }
109 if (mean > 1.) {
110 std::cout << m_className << "::SetGainFluctuationsPolya: "
111 << "Mean avalanche size set to " << mean << ".\n";
112 } else {
113 std::cout << m_className << "::SetGainFluctuationsPolya:\n "
114 << "Mean avalanche size will be given by "
115 << "the integrated Townsend coefficient.\n";
116 }
117 m_gain = mean;
118 m_theta = theta;
119 m_gainFluctuations = GainFluctuations::Polya;
120}
121
122bool DriftLineRKF::DriftElectron(const double x0, const double y0,
123 const double z0, const double t0) {
124 m_particle = Particle::Electron;
125 if (!DriftLine(x0, y0, z0, t0, Particle::Electron, m_t, m_x, m_status)) {
126 return false;
127 }
128 std::vector<double> ne(m_t.size(), 1.);
129 std::vector<double> ni(m_t.size(), 0.);
130 double scale = 1.;
131 if (m_doAvalanche) Avalanche(Particle::Electron, m_x, ne, ni, scale);
132 if (m_doSignal) {
133 ComputeSignal(Particle::Electron, scale * m_scaleE, m_t, m_x, ne);
134 }
135 if (m_doAvalanche && m_doIonTail) AddIonTail(m_t, m_x, ni, scale);
136 return true;
137}
138
139bool DriftLineRKF::AddIonTail(const std::vector<double>& te,
140 const std::vector<std::array<double, 3> >& xe,
141 const std::vector<double>& ni,
142 const double scale) {
143 // SIGETR, SIGIOR
144 const unsigned int nPoints = te.size();
145 if (nPoints < 2) return false;
146 if (ni.size() != nPoints) return false;
147 // Loop over the electron track.
148 for (unsigned int i = 1; i < nPoints; ++i) {
149 // Skip points where there are no ions yet.
150 if (scale * ni[i] < 1.) continue;
151 // Skip also points with a negligible contribution.
152 // if (scale * ni[i] < threshold * m_nI) continue;
153 // Compute the ion drift line.
154 const auto& x0 = xe[i];
155 std::vector<double> ti;
156 std::vector<std::array<double, 3> > xi;
157 int stat = 0;
158 if (!DriftLine(x0[0], x0[1], x0[2], te[i], Particle::Ion, ti, xi, stat)) {
159 std::cerr << m_className << "::AddIonTail:\n"
160 << " Unable to obtain an ion tail; tail not added.\n";
161 return false;
162 }
163 if (m_debug) {
164 std::cout << m_className << "::AddIonTail: Origin = " << PrintVec(x0)
165 << ", n = " << ti.size() << ", status = " << stat << "\n";
166 }
167 // Compute the contribution of the drift line to the signal.
168 ComputeSignal(Particle::Ion, scale * m_scaleI * ni[i], ti, xi, {});
169 }
170 return true;
171}
172
173bool DriftLineRKF::DriftPositron(const double x0, const double y0,
174 const double z0, const double t0) {
175
176 m_particle = Particle::Positron;
177 if (!DriftLine(x0, y0, z0, t0, Particle::Positron, m_t, m_x, m_status)) {
178 return false;
179 }
180 if (m_doSignal) {
181 ComputeSignal(Particle::Positron, m_scaleE, m_t, m_x, {});
182 }
183 return true;
184}
185
186bool DriftLineRKF::DriftHole(const double x0, const double y0, const double z0,
187 const double t0) {
188 m_particle = Particle::Hole;
189 if (!DriftLine(x0, y0, z0, t0, Particle::Hole, m_t, m_x, m_status)) {
190 return false;
191 }
192 if (m_doSignal) ComputeSignal(Particle::Hole, m_scaleH, m_t, m_x, {});
193 return true;
194}
195
196bool DriftLineRKF::DriftIon(const double x0, const double y0, const double z0,
197 const double t0) {
198 m_particle = Particle::Ion;
199 if (!DriftLine(x0, y0, z0, t0, Particle::Ion, m_t, m_x, m_status)) {
200 return false;
201 }
202 if (m_doSignal) ComputeSignal(Particle::Ion, m_scaleI, m_t, m_x, {});
203 return true;
204}
205
206bool DriftLineRKF::DriftNegativeIon(const double x0, const double y0,
207 const double z0, const double t0) {
208
209 m_particle = Particle::NegativeIon;
210 if (!DriftLine(x0, y0, z0, t0, Particle::NegativeIon, m_t, m_x, m_status)) {
211 return false;
212 }
213 if (m_doSignal) {
214 ComputeSignal(Particle::NegativeIon, m_scaleI, m_t, m_x, {});
215 }
216 return true;
217}
218
219bool DriftLineRKF::DriftLine(const double xi, const double yi, const double zi,
220 const double ti, const Particle particle,
221 std::vector<double>& ts,
222 std::vector<Vec>& xs, int& flag) {
223
224 // -----------------------------------------------------------------------
225 // DLCALC - Subroutine doing the actual drift line calculations.
226 // The calculations are based on a Runge-Kutta-Fehlberg method
227 // which has the advantage of controlling the stepsize and the
228 // error while needing only relatively few calls to EFIELD.
229 // Full details are given in the reference quoted below.
230 // VARIABLES : H : Current stepsize (it is in fact a delta t).
231 // HPREV : Stores the previous value of H (comparison)
232 // INITCH : Used for checking initial stepsize (1 = ok)
233 // Other variables such as F0, F1, F2, F3, PHII, PHIII,
234 // CI. ,CII. , BETA.. etc are explained in the reference.
235 // REFERENCE : Stoer + Bulirsch, Einfuhrung in die Numerische
236 // Mathematic II, chapter 7, page 122, 1978, HTB, Springer.
237 // -----------------------------------------------------------------------
238
239
240 // Reset the drift line.
241 ts.clear();
242 xs.clear();
243 // Reset the status flag.
244 flag = StatusAlive;
245
246 // Check if the sensor is defined.
247 if (!m_sensor) {
248 std::cerr << m_className << "::DriftLine: Sensor is not defined.\n";
249 flag = StatusCalculationAbandoned;
250 return false;
251 }
252
253 // Get the sensor's bounding box.
254 double xmin = 0., xmax = 0.;
255 double ymin = 0., ymax = 0.;
256 double zmin = 0., zmax = 0.;
257 bool bbox = m_sensor->GetArea(xmin, ymin, zmin, xmax, ymax, zmax);
258
259 // Initialise the current position.
260 Vec x0 = {xi, yi, zi};
261 // Make sure the initial position is at a valid location.
262 double ex = 0., ey = 0., ez = 0.;
263 double bx = 0., by = 0., bz = 0.;
264 Medium* medium = nullptr;
265 if (GetField(x0, ex, ey, ez, bx, by, bz, medium) != 0) {
266 std::cerr << m_className << "::DriftLine:\n"
267 << " No valid field at initial position.\n";
268 flag = StatusLeftDriftMedium;
269 return false;
270 }
271
272 // Start plotting a new line if requested.
273 int iLine = 0;
274 if (m_view) {
275 if (particle == Particle::Ion || particle == Particle::NegativeIon) {
276 m_view->NewIonDriftLine(1, iLine, xi, yi, zi);
277 } else if (particle == Particle::Electron ||
278 particle == Particle::Positron) {
279 m_view->NewElectronDriftLine(1, iLine, xi, yi, zi);
280 } else if (particle == Particle::Hole) {
281 m_view->NewHoleDriftLine(1, iLine, xi, yi, zi);
282 }
283 }
284
285 // Set the numerical constants for the RKF integration.
286 constexpr double c10 = 214. / 891.;
287 constexpr double c11 = 1. / 33.;
288 constexpr double c12 = 650. / 891.;
289 constexpr double c20 = 533. / 2106.;
290 constexpr double c22 = 800. / 1053.;
291 constexpr double c23 = -1. / 78.;
292
293 constexpr double b10 = 1. / 4.;
294 constexpr double b20 = -189. / 800.;
295 constexpr double b21 = 729. / 800.;
296 constexpr double b30 = 214. / 891.;
297 constexpr double b31 = 1. / 33.;
298 constexpr double b32 = 650. / 891.;
299
300 // Set the charge of the drifting particle.
301 const double charge = particle == Particle::Electron ? -1 : 1;
302
303 // Initialise the particle velocity.
304 Vec v0 = {0., 0., 0.};
305 if (!GetVelocity(ex, ey, ez, bx, by, bz, medium, particle, v0)) {
306 std::cerr << m_className << "::DriftLine:\n"
307 << " Cannot retrieve drift velocity.\n";
308 return false;
309 }
310
311 const double speed0 = Mag(v0);
312 if (speed0 < Small) {
313 std::cerr << m_className << "::DriftLine:\n"
314 << " Zero velocity at initial position.\n";
315 return false;
316 }
317
318 // Initialise time step and previous time step.
319 double h = m_accuracy / speed0;
320 double hprev = h;
321 double t0 = ti;
322
323 // Set the initial point.
324 ts.push_back(t0);
325 xs.push_back(x0);
326
327 int initCycle = 3;
328 bool ok = true;
329 while (ok) {
330 // Get the velocity at the first probe point.
331 Vec x1 = x0;
332 for (unsigned int i = 0; i < 3; ++i) {
333 x1[i] += h * b10 * v0[i];
334 }
335 Vec v1;
336 int stat = 0;
337 if (!GetVelocity(x1, particle, v1, stat)) {
338 flag = StatusCalculationAbandoned;
339 break;
340 } else if (stat != 0) {
341 if (m_debug) {
342 std::cout << m_className << "::DriftLine: Point 1 outside.\n";
343 }
344 if (!Terminate(x0, x1, particle, ts, xs)) {
345 flag = StatusCalculationAbandoned;
346 } else {
347 flag = stat;
348 }
349 break;
350 }
351 // Get the velocity at the second probe point.
352 Vec x2 = x0;
353 for (unsigned int i = 0; i < 3; ++i) {
354 x2[i] += h * (b20 * v0[i] + b21 * v1[i]);
355 }
356 Vec v2;
357 if (!GetVelocity(x2, particle, v2, stat)) {
358 flag = StatusCalculationAbandoned;
359 break;
360 } else if (stat != 0) {
361 if (m_debug) {
362 std::cout << m_className << "::DriftLine: Point 2 outside.\n";
363 }
364 if (!Terminate(x0, x2, particle, ts, xs)) {
365 flag = StatusCalculationAbandoned;
366 } else {
367 flag = stat;
368 }
369 break;
370 }
371 // Get the velocity at the third probe point.
372 Vec x3 = x0;
373 for (unsigned int i = 0; i < 3; ++i) {
374 x3[i] += h * (b30 * v0[i] + b31 * v1[i] + b32 * v2[i]);
375 }
376 Vec v3;
377 if (!GetVelocity(x3, particle, v3, stat)) {
378 flag = StatusCalculationAbandoned;
379 break;
380 } else if (stat != 0) {
381 if (m_debug) {
382 std::cout << m_className << "::DriftLine: Point 3 outside.\n";
383 }
384 if (!Terminate(x0, x3, particle, ts, xs)) {
385 flag = StatusCalculationAbandoned;
386 } else {
387 flag = stat;
388 }
389 break;
390 }
391 // Check if we crossed a wire.
392 double xw = 0., yw = 0., zw = 0., rw = 0.;
393 if (m_sensor->IsWireCrossed(x0[0], x0[1], x0[2],
394 x1[0], x1[1], x1[2], xw, yw, zw, true, rw) ||
395 m_sensor->IsWireCrossed(x0[0], x0[1], x0[2],
396 x2[0], x2[1], x2[2], xw, yw, zw, true, rw) ||
397 m_sensor->IsWireCrossed(x0[0], x0[1], x0[2],
398 x3[0], x3[1], x3[2], xw, yw, zw, true, rw)) {
399 if (m_debug) {
400 std::cout << m_className << "::DriftLine: Crossed wire.\n";
401 }
402 if (DriftToWire(xw, yw, rw, particle, ts, xs, stat)) {
403 flag = stat;
404 } else if (h > Small) {
405 h *= 0.5;
406 continue;
407 } else {
408 std::cerr << m_className << "::DriftLine: Step size too small. Stop.\n";
409 flag = StatusCalculationAbandoned;
410 }
411 break;
412 }
413 // Check if we are inside the trap radius of a wire.
414 if (particle != Particle::Ion) {
415 if (m_sensor->IsInTrapRadius(charge, x1[0], x1[1], x1[2], xw, yw, rw)) {
416 if (!DriftToWire(xw, yw, rw, particle, ts, xs, flag)) {
417 flag = StatusCalculationAbandoned;
418 }
419 break;
420 }
421 if (m_sensor->IsInTrapRadius(charge, x2[0], x2[1], x2[2], xw, yw, rw)) {
422 if (!DriftToWire(xw, yw, rw, particle, ts, xs, flag)) {
423 flag = StatusCalculationAbandoned;
424 }
425 break;
426 }
427 if (m_sensor->IsInTrapRadius(charge, x3[0], x3[1], x3[2], xw, yw, rw)) {
428 if (!DriftToWire(xw, yw, rw, particle, ts, xs, flag)) {
429 flag = StatusCalculationAbandoned;
430 }
431 break;
432 }
433 }
434 // Calculate the correction terms.
435 Vec phi1 = {0., 0., 0.};
436 Vec phi2 = {0., 0., 0.};
437 for (unsigned int i = 0; i < 3; ++i) {
438 phi1[i] = c10 * v0[i] + c11 * v1[i] + c12 * v2[i];
439 phi2[i] = c20 * v0[i] + c22 * v2[i] + c23 * v3[i];
440 }
441 // Check if the step length is valid.
442 const double phi1mag = Mag(phi1);
443 if (phi1mag < Small) {
444 std::cerr << m_className << "::DriftLine: Step has zero length. Stop.\n";
445 flag = StatusCalculationAbandoned;
446 break;
447 } else if (m_useStepSizeLimit && h * phi1mag > m_maxStepSize) {
448 if (m_debug) {
449 std::cout << m_className << "::DriftLine: Step is considered too long. "
450 << "H is reduced.\n";
451 }
452 h = 0.5 * m_maxStepSize / phi1mag;
453 continue;
454 } else if (bbox) {
455 // Don't allow h to become too large in view of the time resolution.
456 if (h * fabs(phi1[0]) > 0.1 * fabs(xmax - xmin) ||
457 h * fabs(phi1[1]) > 0.1 * fabs(ymax - ymin)) {
458 h *= 0.5;
459 if (m_debug) {
460 std::cout << m_className << "::DriftLine: Step is considered too long. "
461 << "H is divided by two.\n";
462 }
463 continue;
464 }
465 } else if (m_rejectKinks && xs.size() > 1) {
466 const unsigned int np = xs.size();
467 const auto& x = xs[np - 1];
468 const auto& xp = xs[np - 2];
469 if (phi1[0] * (x[0] - xp[0]) + phi1[1] * (x[1] - xp[1]) +
470 phi1[2] * (x[2] - xp[2]) < 0.) {
471 std::cerr << m_className << "::DriftLine: Bending angle > 90 degree.\n";
472 flag = StatusSharpKink;
473 break;
474 }
475 }
476 if (m_debug) std::cout << m_className << "::DriftLine: Step size ok.\n";
477 // Update the position and time.
478 for (unsigned int i = 0; i < 3; ++i) x0[i] += h * phi1[i];
479 t0 += h;
480 // Check the new position.
481 m_sensor->ElectricField(x0[0], x0[1], x0[2], ex, ey, ez, medium, stat);
482 if (stat != 0) {
483 // The new position is not inside a valid drift medium.
484 // Terminate the drift line.
485 if (m_debug) {
486 std::cout << m_className << "::DriftLine: Point outside. Terminate.\n";
487 }
488 if (!Terminate(xs.back(), x0, particle, ts, xs)) {
489 flag = StatusCalculationAbandoned;
490 }
491 break;
492 }
493 // Add the new point to the drift line.
494 ts.push_back(t0);
495 xs.push_back(x0);
496 // Adjust the step size according to the accuracy of the two estimates.
497 hprev = h;
498 const double dphi = fabs(phi1[0] - phi2[0]) + fabs(phi1[1] - phi2[1]) +
499 fabs(phi1[2] - phi2[2]);
500 if (dphi > 0) {
501 h = sqrt(h * m_accuracy / dphi);
502 if (m_debug) {
503 std::cout << m_className << "::DriftLine: Adapting H to " << h << ".\n";
504 }
505 } else {
506 h *= 2;
507 if (m_debug) {
508 std::cout << m_className << "::DriftLine: H increased by factor two.\n";
509 }
510 }
511 // Make sure that H is different from zero; this should always be ok.
512 if (h < Small) {
513 std::cerr << m_className << "::DriftLine: Step size is zero. Stop.\n";
514 flag = StatusCalculationAbandoned;
515 break;
516 }
517 // Check the initial step size.
518 if (initCycle > 0 && h < 0.2 * hprev) {
519 if (m_debug) {
520 std::cout << m_className << "::DriftLine: Reinitialise step size.\n";
521 }
522 --initCycle;
523 t0 = ti;
524 x0 = {xi, yi, zi};
525 ts = {t0};
526 xs = {x0};
527 continue;
528 }
529 initCycle = 0;
530 // Don't allow H to grow too quickly
531 if (h > 10 * hprev) {
532 h = 10 * hprev;
533 if (m_debug) {
534 std::cout << m_className << "::DriftLine: H restricted to 10 times "
535 << "the previous value.\n";
536 }
537 }
538 // Stop in case H tends to become too small.
539 if (h * (fabs(phi1[0]) + fabs(phi1[1]) + fabs(phi1[2])) < m_accuracy) {
540 std::cerr << m_className << "::DriftLine: Step size has become smaller "
541 << "than int. accuracy. Stop.\n";
542 flag = StatusCalculationAbandoned;
543 break;
544 }
545 // Update the velocity.
546 v0 = v3;
547 }
548 if (m_view) {
549 for (const auto& x : m_x) {
550 m_view->AddDriftLinePoint(iLine, x[0], x[1], x[2]);
551 }
552 }
553 if (flag == StatusCalculationAbandoned) return false;
554 return true;
555}
556
557bool DriftLineRKF::Avalanche(const Particle particle,
558 const std::vector<Vec>& xs,
559 std::vector<double>& ne,
560 std::vector<double>& ni, double& scale) {
561
562 // SIGETR
563 const unsigned int nPoints = xs.size();
564 if (nPoints < 2) return true;
565 // Locations and weights for 6-point Gaussian integration
566 constexpr double tg[6] = {-0.932469514203152028, -0.661209386466264514,
567 -0.238619186083196909, 0.238619186083196909,
568 0.661209386466264514, 0.932469514203152028};
569 constexpr double wg[6] = {0.171324492379170345, 0.360761573048138608,
570 0.467913934572691047, 0.467913934572691047,
571 0.360761573048138608, 0.171324492379170345};
572
573 ne.assign(nPoints, 1.);
574 ni.assign(nPoints, 0.);
575 bool start = false;
576 bool overflow = false;
577 // Loop over the drift line.
578 for (unsigned int i = 1; i < nPoints; ++i) {
579 const auto& xp = xs[i - 1];
580 const auto& x = xs[i];
581 const Vec dx = {x[0] - xp[0], x[1] - xp[1], x[2] - xp[2]};
582 // Calculate the integrated Townsend and attachment coefficients.
583 double alpsum = 0.;
584 double etasum = 0.;
585 for (unsigned int j = 0; j < 6; ++j) {
586 const double f = 0.5 * (1. + tg[j]);
587 Vec xj = xp;
588 for (unsigned int k = 0; k < 3; ++k) xj[k] += f * dx[k];
589 double ex = 0., ey = 0., ez = 0.;
590 double bx = 0., by = 0., bz = 0.;
591 Medium* medium = nullptr;
592 if (GetField(xj, ex, ey, ez, bx, by, bz, medium) != 0) {
593 std::cerr << m_className << "::Avalanche:\n Invalid field at "
594 << "drift line point " << i << ", segment " << j << ".\n";
595 continue;
596 }
597 double alp = 0.;
598 if (!GetAlpha(ex, ey, ez, bx, by, bz, medium, particle, alp)) {
599 std::cerr << m_className << "::Avalanche:\n Cannot retrieve alpha at "
600 << "drift line point " << i << ", segment " << j << ".\n";
601 continue;
602 }
603 double eta = 0.;
604 if (!GetEta(ex, ey, ez, bx, by, bz, medium, particle, eta)) {
605 std::cerr << m_className << "::Avalanche:\n Cannot retrieve eta at "
606 << "drift line point " << i << ", segment " << j << ".\n";
607 continue;
608 }
609 alpsum += wg[j] * alp;
610 etasum += wg[j] * eta;
611 }
612 alpsum *= 0.5;
613 etasum *= 0.5;
614 if (alpsum > 1.e-6 && !start) {
615 if (m_debug) {
616 std::cout << m_className << "::Avalanche: Avalanche starts at step "
617 << i << ".\n";
618 }
619 start = true;
620 }
621 const double d = Mag(dx);
622 // Update the numbers of electrons and ions.
623 constexpr double expmax = 30.;
624 const double logp = log(std::max(1., ne[i - 1]));
625 if (logp + d * (alpsum - etasum) > expmax) {
626 overflow = true;
627 ne[i] = exp(expmax);
628 } else {
629 ne[i] = ne[i - 1] * exp(d * (alpsum - etasum));
630 }
631 if (logp + d * alpsum > expmax) {
632 overflow = true;
633 ni[i] = exp(expmax);
634 } else {
635 ni[i] = ne[i - 1] * (exp(d * alpsum) - 1);
636 }
637 }
638 if (overflow) {
639 std::cerr << m_className << "::Avalanche:\n "
640 << "Warning: Integrating the Townsend coefficients "
641 << "would lead to exponential overflow.\n "
642 << "Avalanche truncated.\n";
643 }
644 const double qe = ne.back();
645 const double qi = std::accumulate(ni.begin(), ni.end(), 0.);
646 scale = 1.;
647 if (qi > 1. &&
648 !(m_gainFluctuations == GainFluctuations::None && m_gain < 1.)) {
649 const double gain = m_gain > 1. ? m_gain : GetGain();
650 double q1 = gain;
651 if (m_gainFluctuations == GainFluctuations::Polya) {
652 for (unsigned int i = 0; i < 100; ++i) {
653 q1 = gain * RndmPolya(m_theta);
654 if (q1 >= 1.) break;
655 }
656 q1 = std::max(q1, 1.);
657 }
658 q1 *= GetLoss();
659 scale = (q1 + 1.) / (qi + 1.);
660 }
661 if (m_debug) {
662 std::cout << m_className << "::Avalanche:\n "
663 << "Final number of electrons: " << qe << "\n "
664 << "Number of ions: " << qi << "\n "
665 << "Charge scaling factor: " << scale << "\n "
666 << "Avalanche development:\n Step Electrons Ions\n";
667 for (unsigned int i = 0; i < nPoints; ++i) {
668 std::printf("%6d %15.7f %15.7f\n", i, scale * ne[i], scale * ni[i]);
669 }
670 }
671 m_nE = scale * qe;
672 m_nI = scale * qi;
673 return true;
674}
675
676double DriftLineRKF::GetArrivalTimeSpread(const double eps) {
677
678 // -----------------------------------------------------------------------
679 // DLCDF1 - Routine returning the integrated diffusion coefficient of
680 // the current drift line. The routine uses an adaptive
681 // Simpson integration.
682 // -----------------------------------------------------------------------
683
684 const unsigned int nPoints = m_x.size();
685 // Return straight away if there is only one point.
686 if (nPoints < 2) return 0.;
687 const Particle particle = m_particle;
688
689 // First get a rough estimate of the result.
690 double crude = 0.;
691 double varPrev = 0.;
692 for (unsigned int i = 0; i < nPoints; ++i) {
693 // Get the drift velocity and diffusion coefficients at this step.
694 double ex = 0., ey = 0., ez = 0.;
695 double bx = 0., by = 0., bz = 0.;
696 Medium* medium = nullptr;
697 if (GetField(m_x[i], ex, ey, ez, bx, by, bz, medium) != 0) {
698 std::cerr << m_className << "::GetArrivalTimeSpread:\n"
699 << " Invalid drift line point " << i << ".\n";
700 continue;
701 }
702 Vec v;
703 if (!GetVelocity(ex, ey, ez, bx, by, bz, medium, particle, v)) {
704 std::cerr << m_className << "::GetArrivalTimeSpread:\n"
705 << " Cannot retrieve drift velocity at point " << i << ".\n";
706 continue;
707 }
708 const double speed = Mag(v);
709 if (speed < Small) {
710 std::cerr << m_className << "::GetArrivalTimeSpread:\n"
711 << " Zero drift velocity at point " << i << ".\n";
712 continue;
713 }
714 double dl = 0., dt = 0.;
715 if (!GetDiffusion(ex, ey, ez, bx, by, bz, medium, particle, dl, dt)) {
716 std::cerr << m_className << "::GetArrivalTimeSpread:\n"
717 << " Cannot retrieve diffusion at point " << i << ".\n";
718 continue;
719 }
720 const double sigma = dl / speed;
721 const double var = sigma * sigma;
722 if (i > 0) {
723 const auto& x = m_x[i];
724 const auto& xPrev = m_x[i - 1];
725 const double d = Mag(x[0] - xPrev[0], x[1] - xPrev[1], x[2] - xPrev[2]);
726 crude += 0.5 * d * (var + varPrev);
727 }
728 varPrev = var;
729 }
730 crude = sqrt(crude);
731
732 const double tol = eps * crude;
733 double sum = 0.;
734 for (unsigned int i = 0; i < nPoints - 1; ++i) {
735 sum += IntegrateDiffusion(m_x[i], m_x[i + 1], particle, tol);
736 }
737 return sqrt(sum);
738}
739
740double DriftLineRKF::GetGain(const double eps) {
741
742 // -----------------------------------------------------------------------
743 // DLCTW1 - Routine returning the multiplication factor for the current
744 // drift line. The routine uses an adaptive Simpson style
745 // integration.
746 // -----------------------------------------------------------------------
747
748 const unsigned int nPoints = m_x.size();
749 // Return straight away if there is only one point.
750 if (nPoints < 2) return 1.;
751 const Particle particle = m_particle;
752 if (particle == Particle::Ion) return 1.;
753 if (m_status == StatusCalculationAbandoned) return 1.;
754
755 // First get a rough estimate of the result.
756 double crude = 0.;
757 double alphaPrev = 0.;
758 for (unsigned int i = 0; i < nPoints; ++i) {
759 // Get the Townsend coefficient at this step.
760 double ex = 0., ey = 0., ez = 0.;
761 double bx = 0., by = 0., bz = 0.;
762 Medium* medium = nullptr;
763 if (GetField(m_x[i], ex, ey, ez, bx, by, bz, medium) != 0) {
764 std::cerr << m_className << "::GetGain:\n"
765 << " Invalid drift line point " << i << ".\n";
766 continue;
767 }
768 double alpha = 0.;
769 if (!GetAlpha(ex, ey, ez, bx, by, bz, medium, particle, alpha)) {
770 std::cerr << m_className << "::GetGain:\n"
771 << " Cannot retrieve alpha at point " << i << ".\n";
772 continue;
773 }
774 if (i > 0) {
775 const auto& x = m_x[i];
776 const auto& xPrev = m_x[i - 1];
777 const double d = Mag(x[0] - xPrev[0], x[1] - xPrev[1], x[2] - xPrev[2]);
778 crude += 0.5 * d * (alpha + alphaPrev);
779 }
780 alphaPrev = alpha;
781 }
782 // Stop if the rough estimate is negligibly small.
783 if (crude < Small) return 1.;
784
785 // Calculate the integration tolerance based on the rough estimate.
786 const double tol = eps * crude;
787 double sum = 0.;
788 for (unsigned int i = 0; i < nPoints - 1; ++i) {
789 sum += IntegrateAlpha(m_x[i], m_x[i + 1], particle, tol);
790 }
791 return exp(sum);
792}
793
794double DriftLineRKF::GetLoss(const double eps) {
795
796 // -----------------------------------------------------------------------
797 // DLCAT1 - Routine returning the attachment losses for the current
798 // drift line. The routine uses an adaptive Simpson style
799 // integration.
800 // -----------------------------------------------------------------------
801
802 const unsigned int nPoints = m_x.size();
803 // Return straight away if there is only one point.
804 if (nPoints < 2) return 1.;
805 const Particle particle = m_particle;
806 if (particle == Particle::Ion) return 1.;
807 if (m_status == StatusCalculationAbandoned) return 1.;
808
809 // First get a rough estimate of the result.
810 double crude = 0.;
811 double etaPrev = 0.;
812 for (unsigned int i = 0; i < nPoints; ++i) {
813 // Get the Townsend coefficient at this step.
814 double ex = 0., ey = 0., ez = 0.;
815 double bx = 0., by = 0., bz = 0.;
816 Medium* medium = nullptr;
817 if (GetField(m_x[i], ex, ey, ez, bx, by, bz, medium) != 0) {
818 std::cerr << m_className << "::GetLoss:\n"
819 << " Invalid drift line point " << i << ".\n";
820 continue;
821 }
822 double eta = 0.;
823 if (!GetEta(ex, ey, ez, bx, by, bz, medium, particle, eta)) {
824 std::cerr << m_className << "::GetLoss:\n"
825 << " Cannot retrieve eta at point " << i << ".\n";
826 continue;
827 }
828 if (i > 0) {
829 const auto& x = m_x[i];
830 const auto& xPrev = m_x[i - 1];
831 const double d = Mag(x[0] - xPrev[0], x[1] - xPrev[1], x[2] - xPrev[2]);
832 crude += 0.5 * d * (eta + etaPrev);
833 }
834 etaPrev = eta;
835 }
836
837 // Calculate the integration tolerance based on the rough estimate.
838 const double tol = eps * crude;
839 double sum = 0.;
840 for (unsigned int i = 0; i < nPoints - 1; ++i) {
841 sum += IntegrateEta(m_x[i], m_x[i + 1], particle, tol);
842 }
843 return exp(-sum);
844}
845
846int DriftLineRKF::GetField(const std::array<double, 3>& x,
847 double& ex, double& ey, double& ez,
848 double& bx, double& by, double& bz,
849 Medium*& medium) const {
850 int status = 0;
851 m_sensor->MagneticField(x[0], x[1], x[2], bx, by, bz, status);
852 m_sensor->ElectricField(x[0], x[1], x[2], ex, ey, ez, medium, status);
853 return status;
854}
855
856bool DriftLineRKF::GetVelocity(const std::array<double, 3>& x,
857 const Particle particle,
858 std::array<double, 3>& v, int& status) const {
859 v.fill(0.);
860 double ex = 0., ey = 0., ez = 0.;
861 double bx = 0., by = 0., bz = 0.;
862 Medium* medium = nullptr;
863 // Stop if we are outside a valid drift medium.
864 status = GetField(x, ex, ey, ez, bx, by, bz, medium);
865 if (status != 0) return true;
866 if (particle == Particle::Electron) {
867 return medium->ElectronVelocity(ex, ey, ez, bx, by, bz, v[0], v[1], v[2]);
868 } else if (particle == Particle::Ion) {
869 return medium->IonVelocity(ex, ey, ez, bx, by, bz, v[0], v[1], v[2]);
870 } else if (particle == Particle::Hole) {
871 return medium->HoleVelocity(ex, ey, ez, bx, by, bz, v[0], v[1], v[2]);
872 } else if (particle == Particle::Positron) {
873 const bool ok = medium->ElectronVelocity(ex, ey, ez, bx, by, bz,
874 v[0], v[1], v[2]);
875 for (unsigned int i = 0; i < 3; ++i) v[i] *= -1;
876 return ok;
877 } else if (particle == Particle::NegativeIon) {
878 const bool ok = medium->IonVelocity(ex, ey, ez, bx, by, bz,
879 v[0], v[1], v[2]);
880 for (unsigned int i = 0; i < 3; ++i) v[i] *= -1;
881 return ok;
882 }
883 std::cerr << m_className << "::GetVelocity:\n"
884 << " Cannot retrieve drift velocity at " << PrintVec(x) << ".\n";
885 return false;
886}
887
888bool DriftLineRKF::GetVelocity(const double ex, const double ey,
889 const double ez, const double bx,
890 const double by, const double bz,
891 Medium* medium, const Particle particle,
892 std::array<double, 3>& v) const {
893 v.fill(0.);
894 if (particle == Particle::Electron) {
895 return medium->ElectronVelocity(ex, ey, ez, bx, by, bz, v[0], v[1], v[2]);
896 } else if (particle == Particle::Ion) {
897 return medium->IonVelocity(ex, ey, ez, bx, by, bz, v[0], v[1], v[2]);
898 } else if (particle == Particle::Hole) {
899 return medium->HoleVelocity(ex, ey, ez, bx, by, bz, v[0], v[1], v[2]);
900 } else if (particle == Particle::Positron) {
901 const bool ok = medium->ElectronVelocity(ex, ey, ez, bx, by, bz,
902 v[0], v[1], v[2]);
903 for (unsigned int i = 0; i < 3; ++i) v[i] *= -1;
904 return ok;
905 } else if (particle == Particle::NegativeIon) {
906 const bool ok = medium->IonVelocity(ex, ey, ez, bx, by, bz,
907 v[0], v[1], v[2]);
908 for (unsigned int i = 0; i < 3; ++i) v[i] *= -1;
909 return ok;
910 }
911 return false;
912}
913
914bool DriftLineRKF::GetDiffusion(const double ex, const double ey,
915 const double ez, const double bx,
916 const double by, const double bz,
917 Medium* medium, const Particle particle,
918 double& dl, double& dt) const {
919 if (particle == Particle::Electron || particle == Particle::Positron) {
920 return medium->ElectronDiffusion(ex, ey, ez, bx, by, bz, dl, dt);
921 } else if (particle == Particle::Ion || particle == Particle::NegativeIon) {
922 return medium->IonDiffusion(ex, ey, ez, bx, by, bz, dl, dt);
923 } else if (particle == Particle::Hole) {
924 return medium->HoleDiffusion(ex, ey, ez, bx, by, bz, dl, dt);
925 }
926 return false;
927}
928
929bool DriftLineRKF::GetAlpha(const double ex, const double ey, const double ez,
930 const double bx, const double by, const double bz,
931 Medium* medium, const Particle particle,
932 double& alpha) const {
933 if (particle == Particle::Electron || particle == Particle::Positron) {
934 return medium->ElectronTownsend(ex, ey, ez, bx, by, bz, alpha);
935 } else if (particle == Particle::Hole) {
936 return medium->HoleTownsend(ex, ey, ez, bx, by, bz, alpha);
937 }
938 return false;
939}
940
941bool DriftLineRKF::GetEta(const double ex, const double ey, const double ez,
942 const double bx, const double by, const double bz,
943 Medium* medium, const Particle particle,
944 double& eta) const {
945 if (particle == Particle::Electron) {
946 return medium->ElectronAttachment(ex, ey, ez, bx, by, bz, eta);
947 } else if (particle == Particle::Hole) {
948 return medium->HoleAttachment(ex, ey, ez, bx, by, bz, eta);
949 }
950 return false;
951}
952
953bool DriftLineRKF::Terminate(const std::array<double, 3>& xx0,
954 const std::array<double, 3>& xx1,
955 const Particle particle,
956 std::vector<double>& ts,
957 std::vector<Vec>& xs) {
958
959 // -----------------------------------------------------------------------
960 // DLCFMP - Terminates drift line calculation by making a last step
961 // to the boundary of the mesh or the drift medium.
962 // VARIABLES : XX0: Last point in drift medium.
963 // XX1: Estimated step, outside drift medium.
964 // -----------------------------------------------------------------------
965
966 // Check the validity of the initial point.
967 Vec vv0 = {0., 0., 0.};
968 int status = 0;
969 if (!GetVelocity(xx0, particle, vv0, status)) {
970 std::cerr << m_className << "::Terminate:\n"
971 << " Cannot retrieve initial drift velocity.\n";
972 return false;
973 } else if (status != 0) {
974 std::cerr << m_className << "::Terminate:\n"
975 << " No valid field at initial point. Program bug!\n";
976 return false;
977 }
978 double speed = Mag(vv0);
979 if (speed < Small) {
980 std::cerr << m_className << "::Terminate:\n"
981 << " Zero velocity at initial position.\n";
982 return false;
983 }
984
985 // Final point just inside the medium.
986 Vec x0 = xx0;
987 // Final point just outside the medium.
988 Vec x1 = xx1;
989 // Perform some bisections.
990 const unsigned int nBisections = 20;
991 for (unsigned int i = 0; i < nBisections; ++i) {
992 // Quit bisection when interval becomes too small.
993 bool small = true;
994 for (unsigned int j = 0; j < 3; ++j) {
995 if (fabs(x1[j] - x0[j]) > 1.e-6 * (fabs(x0[j]) + fabs(x1[j]))) {
996 small = false;
997 break;
998 }
999 }
1000 if (small) {
1001 if (m_debug) {
1002 std::cout << m_className << "::Terminate:\n"
1003 << " Bisection ended at cycle " << i << ".\n";
1004 }
1005 break;
1006 }
1007 // Calculate the mid point and evaluate the field.
1008 Vec xm;
1009 for (unsigned int j = 0; j < 3; ++j) xm[j] = 0.5 * (x0[j] + x1[j]);
1010 double ex = 0., ey = 0., ez = 0.;
1011 Medium* medium = nullptr;
1012 m_sensor->ElectricField(xm[0], xm[1], xm[2], ex, ey, ez, medium, status);
1013 if (status == 0) {
1014 x0 = xm;
1015 } else {
1016 x1 = xm;
1017 }
1018 }
1019
1020 // Compute drift velocity at the end of the step.
1021 Vec v0;
1022 if (!GetVelocity(x0, particle, v0, status) || status != 0) {
1023 std::cerr << m_className << "::Terminate:\n"
1024 << " Warning: Unable to compute mean velocity at last step.\n";
1025 } else {
1026 speed = 0.5 * (speed + Mag(v0));
1027 }
1028 // Calculate the time step.
1029 const double dt = Mag(x0[0] - xx0[0], x0[1] - xx0[1], x0[2] - xx0[2]) / speed;
1030 // Add the last point, just inside the drift area.
1031 ts.push_back(ts.back() + dt);
1032 xs.push_back(x0);
1033 m_status = StatusLeftDriftMedium;
1034 return true;
1035}
1036
1037bool DriftLineRKF::DriftToWire(const double xw, const double yw,
1038 const double rw, const Particle particle,
1039 std::vector<double>& ts,
1040 std::vector<Vec>& xs, int& stat) {
1041
1042 // -----------------------------------------------------------------------
1043 // DLCWIR - Terminates drift line calculation by making some last steps
1044 // towards the surface of the wire on which it is supposed to
1045 // end. The precision is controlled in order to obtain a
1046 // good estimate of the total remaining drift-time.
1047 // -----------------------------------------------------------------------
1048
1049 // Get the starting point.
1050 Vec x0 = xs.back();
1051 double t0 = ts.back() - ts.front();
1052 if (m_debug) {
1053 std::cout << m_className << "::DriftToWire:\n Drifting from ("
1054 << x0[0] << ", " << x0[1] << ") to wire at ("
1055 << xw << ", " << yw << ") with radius " << rw << " cm.\n";
1056 }
1057
1058 // Get the initial drift velocity.
1059 Vec v0;
1060 int status = 0;
1061 if (!GetVelocity(x0, particle, v0, status)) {
1062 std::cerr << m_className << "::DriftToWire:\n"
1063 << " Cannot retrieve initial drift velocity.\n";
1064 return false;
1065 } else if (status != 0) {
1066 std::cerr << m_className << "::DriftToWire:\n"
1067 << " No valid field at initial point. Program bug!\n";
1068 return false;
1069 }
1070
1071 // Estimate the time needed to reach the wire
1072 // assuming a straight-line trajectory and constant velocity.
1073 double dt = (Mag(xw - x0[0], yw - x0[1]) - rw) / Mag(v0[0], v0[1]);
1074 if (m_debug) {
1075 std::cout << m_className << "::DriftToWire: "
1076 << "Estimated time needed to reach the wire: " << dt << " ns.\n";
1077 }
1078
1079 constexpr unsigned int nMaxSplit = 10;
1080 unsigned int nSplit = 0;
1081 // Move towards the wire.
1082 bool onwire = false;
1083 const double r2 = rw * rw;
1084 while (!onwire && dt > 1.e-6 * t0) {
1085 // Calculate the estimated end point.
1086 Vec x1 = x0;
1087 for (unsigned int j = 0; j < 3; ++j) x1[j] += dt * v0[j];
1088 // Make sure we are not moving away from the wire.
1089 const double xinp0 = (x1[0] - x0[0]) * (xw - x0[0]) +
1090 (x1[1] - x0[1]) * (yw - x0[1]);
1091 if (xinp0 < 0.) {
1092 if (m_debug) {
1093 std::cerr << m_className << "::DriftToWire:\n"
1094 << " Particle moves away from the wire. Quit.\n";
1095 }
1096 return false;
1097 }
1098 // Check if the end point is inside the wire or the wire was crossed.
1099 const double xinp1 = (x0[0] - x1[0]) * (xw - x1[0]) +
1100 (x0[1] - x1[1]) * (yw - x1[1]);
1101 if (xinp1 < 0.) {
1102 if (Mag2(xw - x1[0], yw - x1[1]) <= r2) {
1103 onwire = true;
1104 if (m_debug) std::cout << m_className << "::DriftToWire: Inside.\n";
1105 }
1106 } else {
1107 if (m_debug) std::cout << m_className << "::DriftToWire: Wire crossed.\n";
1108 onwire = true;
1109 }
1110 if (onwire) {
1111 const double dw0 = Mag(xw - x0[0], yw - x0[1]);
1112 x1[0] = xw - (rw + BoundaryDistance) * (xw - x0[0]) / dw0;
1113 x1[1] = yw - (rw + BoundaryDistance) * (yw - x0[1]) / dw0;
1114 dt = Mag(x1[0] - x0[0], x1[1] - x0[1]) / Mag(v0[0], v0[1]);
1115 x1[2] = x0[2] + dt * v0[2];
1116 }
1117 // Calculate the drift velocity at the end point.
1118 Vec v1;
1119 if (!GetVelocity(x1, particle, v1, status)) {
1120 std::cerr << m_className << "::DriftToWire:\n"
1121 << " Cannot retrieve drift velocity at end point. Quit.\n";
1122 return false;
1123 } else if (status != 0) {
1124 std::cerr << m_className << "::DriftToWire:\n"
1125 << " End point is not in a valid drift medium. Quit.\n";
1126 return false;
1127 }
1128 // Get a point halfway between for an accuracy check.
1129 Vec xm;
1130 for (unsigned int j = 0; j < 3; ++j) xm[j] = 0.5 * (x0[j] + x1[j]);
1131 // Calculate the drift velocity at the mid point.
1132 Vec vm;
1133 if (!GetVelocity(xm, particle, vm, status)) {
1134 std::cerr << m_className << "::DriftToWire:\n"
1135 << " Cannot retrieve drift velocity at mid point. Quit.\n";
1136 return false;
1137 } else if (status != 0) {
1138 std::cerr << m_className << "::DriftToWire:\n"
1139 << " Mid point is not in a valid drift medium. Quit.\n";
1140 return false;
1141 }
1142 // Make sure the velocities are non-zero.
1143 const double speed0 = Mag(v0[0], v0[1]);
1144 const double speed1 = Mag(v1[0], v1[1]);
1145 const double speedm = Mag(vm[0], vm[1]);
1146 if (speed0 < Small || speed1 < Small || speedm < Small) {
1147 std::cerr << m_className << "::DriftToWire: Zero velocity. Stop.\n";
1148 return false;
1149 }
1150 const double p0 = 1. / speed0;
1151 const double p1 = 1. / speed1;
1152 const double pm = 1. / speedm;
1153 // Compare first and second order estimates.
1154 const double d = Mag(x0[0] - x1[0], x0[1] - x1[1]);
1155 const double tol = 1.e-4 * (1. + fabs(t0));
1156 if (d * fabs(p0 - 2 * pm + p1) / 3. > tol && nSplit < nMaxSplit) {
1157 // Accuracy was not good enough so halve the step time.
1158 if (m_debug) {
1159 std::cout << m_className << "::DriftToWire: Reducing step size.\n";
1160 }
1161 dt *= 0.5;
1162 onwire = false;
1163 ++nSplit;
1164 continue;
1165 }
1166 const double t1 = t0 + d * (p0 + 4 * pm + p1) / 6.;
1167 // Add a new point to the drift line.
1168 ts.push_back(ts[0] + t1);
1169 xs.push_back(x1);
1170 // Proceed to the next step.
1171 x0 = x1;
1172 t0 = t1;
1173 v0 = v1;
1174 }
1175 // Get the wire index (status code inside the wire).
1176 double ex = 0., ey = 0., ez = 0.;
1177 Medium* medium = nullptr;
1178 m_sensor->ElectricField(xw, yw, 0., ex, ey, ez, medium, stat);
1179 return true;
1180}
1181
1183
1184 std::cout << m_className << "::PrintDriftLine:\n";
1185 if (m_x.empty()) {
1186 std::cout << " No drift line present.\n";
1187 return;
1188 }
1189 if (m_particle == Particle::Electron) {
1190 std::cout << " Particle: electron\n";
1191 } else if (m_particle == Particle::Ion) {
1192 std::cout << " Particle: ion\n";
1193 } else if (m_particle == Particle::Hole) {
1194 std::cout << " Particle: hole\n";
1195 } else {
1196 std::cout << " Particle: unknown\n";
1197 }
1198 std::cout << " Status: " << m_status << "\n"
1199 << " Step time [ns] "
1200 << "x [cm] y [cm] z [cm]\n";
1201 const unsigned int nPoints = m_x.size();
1202 for (unsigned int i = 0; i < nPoints; ++i) {
1203 std::printf("%6d %15.7f %15.7f %15.7f %15.7f\n",
1204 i, m_t[i], m_x[i][0], m_x[i][1], m_x[i][2]);
1205 }
1206
1207}
1208
1209void DriftLineRKF::GetEndPoint(double& x, double& y, double& z, double& t,
1210 int& stat) const {
1211 if (m_x.empty()) {
1212 x = y = z = t = 0.;
1213 stat = m_status;
1214 return;
1215 }
1216 const auto& p = m_x.back();
1217 x = p[0];
1218 y = p[1];
1219 z = p[2];
1220 t = m_t.back();
1221 stat = m_status;
1222}
1223
1224void DriftLineRKF::GetDriftLinePoint(const unsigned int i, double& x, double& y,
1225 double& z, double& t) const {
1226 if (i >= m_x.size()) {
1227 std::cerr << m_className << "::GetDriftLinePoint: Index out of range.\n";
1228 return;
1229 }
1230
1231 const auto& p = m_x[i];
1232 x = p[0];
1233 y = p[1];
1234 z = p[2];
1235 t = m_t[i];
1236}
1237
1238double DriftLineRKF::IntegrateDiffusion(const std::array<double, 3>& xi,
1239 const std::array<double, 3>& xe,
1240 const Particle particle,
1241 const double tol) {
1242
1243 double ex = 0., ey = 0., ez = 0.;
1244 double bx = 0., by = 0., bz = 0.;
1245 Medium* medium = nullptr;
1246
1247 // Make sure the starting point is valid.
1248 Vec x0 = xi;
1249 if (GetField(x0, ex, ey, ez, bx, by, bz, medium) != 0) {
1250 std::cerr << m_className << "::IntegrateDiffusion: Invalid starting point "
1251 << PrintVec(xi) << ".\n";
1252 return 0.;
1253 }
1254
1255 // Determine the drift velocity at the starting point.
1256 Vec v0;
1257 if (!GetVelocity(ex, ey, ez, bx, by, bz, medium, particle, v0)) {
1258 std::cerr << m_className << "::IntegrateDiffusion:\n"
1259 << " Cannot retrieve drift velocity at initial point.\n";
1260 return 0.;
1261 }
1262 double speed0 = Mag(v0);
1263 if (speed0 < Small) {
1264 std::cerr << m_className << "::IntegrateDiffusion:\n"
1265 << " Zero velocity at starting point.\n";
1266 return 0.;
1267 }
1268 // Determine the diffusion coefficients at the initial point.
1269 double dl0 = 0., dt0 = 0.;
1270 if (!GetDiffusion(ex, ey, ez, bx, by, bz, medium, particle, dl0, dt0)) {
1271 std::cerr << m_className << "::IntegrateDiffusion:\n"
1272 << " Cannot retrieve diffusion at initial point.\n";
1273 return 0.;
1274 }
1275 const double sigma0 = dl0 / speed0;
1276 double var0 = sigma0 * sigma0;
1277
1278 // Make sure the end point is valid.
1279 Vec x1 = xe;
1280 if (GetField(x1, ex, ey, ez, bx, by, bz, medium) != 0) {
1281 std::cerr << m_className << "::IntegrateDiffusion: Invalid end point "
1282 << PrintVec(xe) << ".\n";
1283 return 0.;
1284 }
1285 // Determine the drift velocity at the end point.
1286 Vec v1;
1287 if (!GetVelocity(ex, ey, ez, bx, by, bz, medium, particle, v1)) {
1288 std::cerr << m_className << "::IntegrateDiffusion:\n"
1289 << " Cannot retrieve drift velocity at end point.\n";
1290 return 0.;
1291 }
1292 double speed1 = Mag(v1);
1293 if (speed1 < Small) {
1294 std::cerr << m_className << "::IntegrateDiffusion:\n"
1295 << " Zero velocity at end point.\n";
1296 return 0.;
1297 }
1298 // Determine the diffusion coefficients at the end point.
1299 double dl1 = 0., dt1 = 0.;
1300 if (!GetDiffusion(ex, ey, ez, bx, by, bz, medium, particle, dl1, dt1)) {
1301 std::cerr << m_className << "::IntegrateDiffusion:\n"
1302 << " Cannot retrieve diffusion at initial point.\n";
1303 return 0.;
1304 }
1305 double sigma1 = dl1 / speed1;
1306 double var1 = sigma1 * sigma1;
1307
1308 double integral = 0.;
1309 bool done = false;
1310 while (!done) {
1311 // Check if we are close to the end point.
1312 if (Mag(xe[0] - x0[0], xe[1] - x0[1], xe[2] - x0[2]) < 1.e-6) {
1313 done = true;
1314 break;
1315 }
1316 const double d = Mag(x1[0] - x0[0], x1[1] - x0[1], x1[2] - x0[2]);
1317 if (d < 1.e-6) {
1318 // Step length has become very small.
1319 if (m_debug) {
1320 std::cout << m_className << "::IntegrateDiffusion: Small step.\n";
1321 }
1322 integral += var0 * d;
1323 // Proceed with the next step.
1324 x0 = x1;
1325 x1 = xe;
1326 continue;
1327 }
1328 // Calculate drift velocity and diffusion at the end point of the step.
1329 if (GetField(x1, ex, ey, ez, bx, by, bz, medium) != 0) {
1330 std::cerr << m_className << "::IntegrateDiffusion: Invalid end point.\n";
1331 break;
1332 }
1333 if (!GetVelocity(ex, ey, ez, bx, by, bz, medium, particle, v1)) {
1334 std::cerr << m_className << "::IntegrateDiffusion:\n"
1335 << " Cannot retrieve drift velocity at end point.\n";
1336 break;
1337 }
1338 speed1 = Mag(v1);
1339 if (speed1 < Small) {
1340 std::cerr << m_className << "::IntegrateDiffusion:\n"
1341 << " Zero drift velocity at end point.\n";
1342 break;
1343 }
1344 if (!GetDiffusion(ex, ey, ez, bx, by, bz, medium, particle, dl1, dt1)) {
1345 std::cerr << m_className << "::IntegrateDiffusion:\n"
1346 << " Cannot retrieve diffusion at end point.\n";
1347 break;
1348 }
1349 // Determine drift velocity and diffusion at the mid point of the step.
1350 Vec xm;
1351 for (unsigned int i = 0; i < 3; ++i) xm[i] = 0.5 * (x0[i] + x1[i]);
1352 if (GetField(xm, ex, ey, ez, bx, by, bz, medium) != 0) {
1353 std::cerr << m_className << "::IntegrateDiffusion: Invalid mid point.\n";
1354 break;
1355 }
1356 Vec vm;
1357 if (!GetVelocity(ex, ey, ez, bx, by, bz, medium, particle, vm)) {
1358 std::cerr << m_className << "::IntegrateDiffusion:\n"
1359 << " Cannot retrieve drift velocity at mid point.\n";
1360 break;
1361 }
1362 const double speedm = Mag(vm);
1363 if (speedm < Small) {
1364 std::cerr << m_className << "::IntegrateDiffusion:\n"
1365 << " Zero drift velocity at mid point.\n";
1366 break;
1367 }
1368 double dlm = 0., dtm = 0.;
1369 if (!GetDiffusion(ex, ey, ez, bx, by, bz, medium, particle, dlm, dtm)) {
1370 std::cerr << m_className << "::IntegrateDiffusion:\n"
1371 << " Cannot retrieve diffusion at mid point.\n";
1372 break;
1373 }
1374 sigma1 = dl1 / speed1;
1375 var1 = sigma1 * sigma1;
1376 const double sigmam = dlm / speedm;
1377 const double varm = sigmam * sigmam;
1378 // Compare first and second order estimates
1379 // (integrals calculated using trapezoidal and Simpson's rule).
1380 if (fabs(var0 - 2 * varm + var1) * sqrt(d * 2 / (var0 + var1)) / 6. < tol) {
1381 // Accuracy is good enough.
1382 integral += d * (var0 + 4 * varm + var1) / 6.;
1383 // Proceed to the next step.
1384 x0 = x1;
1385 x1 = xe;
1386 var0 = var1;
1387 } else {
1388 // Accuracy is not good enough, so halve the step.
1389 x1 = xm;
1390 var1 = varm;
1391 }
1392 }
1393 return integral;
1394}
1395
1396double DriftLineRKF::IntegrateAlpha(const std::array<double, 3>& xi,
1397 const std::array<double, 3>& xe,
1398 const Particle particle, const double tol) {
1399
1400 double ex = 0., ey = 0., ez = 0.;
1401 double bx = 0., by = 0., bz = 0.;
1402 Medium* medium = nullptr;
1403
1404 // Make sure the starting point is valid.
1405 Vec x0 = xi;
1406 if (GetField(x0, ex, ey, ez, bx, by, bz, medium) != 0) {
1407 std::cerr << m_className << "::IntegrateAlpha: Invalid starting point "
1408 << PrintVec(xi) << ".\n";
1409 return 0.;
1410 }
1411 // Determine the Townsend coefficient at the initial point.
1412 double alpha0 = 0.;
1413 if (!GetAlpha(ex, ey, ez, bx, by, bz, medium, particle, alpha0)) {
1414 std::cerr << m_className << "::IntegrateAlpha:\n"
1415 << " Cannot retrieve Townsend coefficient at initial point.\n";
1416 return 0.;
1417 }
1418 // Make sure the end point is valid.
1419 Vec x1 = xe;
1420 if (GetField(x1, ex, ey, ez, bx, by, bz, medium) != 0) {
1421 std::cerr << m_className << "::IntegrateAlpha: Invalid end point "
1422 << PrintVec(xe) << ".\n";
1423 return 0.;
1424 }
1425 // Determine the Townsend coefficient at the end point.
1426 double alpha1 = 0.;
1427 if (!GetAlpha(ex, ey, ez, bx, by, bz, medium, particle, alpha1)) {
1428 std::cerr << m_className << "::IntegrateAlpha:\n"
1429 << " Cannot retrieve Townsend coefficient at end point.\n";
1430 return 0.;
1431 }
1432 double integral = 0.;
1433 bool done = false;
1434 while (!done) {
1435 // Check if we are close to the end point.
1436 if (Mag(xe[0] - x0[0], xe[1] - x0[1], xe[2] - x0[2]) < 1.e-6) {
1437 done = true;
1438 break;
1439 }
1440 const double d = Mag(x1[0] - x0[0], x1[1] - x0[1], x1[2] - x0[2]);
1441 if (d < 1.e-6) {
1442 // Step length has become very small.
1443 if (m_debug) {
1444 std::cout << m_className << "::IntegrateAlpha: Small step.\n";
1445 }
1446 integral += alpha0 * d;
1447 // Proceed with the next step.
1448 x0 = x1;
1449 x1 = xe;
1450 continue;
1451 }
1452 // Calculate the Townsend coefficient at the end point of the step.
1453 if (GetField(x1, ex, ey, ez, bx, by, bz, medium) != 0) {
1454 std::cerr << m_className << "::IntegrateAlpha: Invalid end point.\n";
1455 break;
1456 }
1457 if (!GetAlpha(ex, ey, ez, bx, by, bz, medium, particle, alpha1)) {
1458 std::cerr << m_className << "::IntegrateAlpha:\n"
1459 << " Cannot retrieve Townsend coefficient at end point.\n";
1460 break;
1461 }
1462 // Calculate the Townsend coefficient at the mid point of the step.
1463 Vec xm;
1464 for (unsigned int i = 0; i < 3; ++i) xm[i] = 0.5 * (x0[i] + x1[i]);
1465 if (GetField(xm, ex, ey, ez, bx, by, bz, medium) != 0) {
1466 std::cerr << m_className << "::IntegrateAlpha: Invalid mid point.\n";
1467 break;
1468 }
1469 double alpham = 0.;
1470 if (!GetAlpha(ex, ey, ez, bx, by, bz, medium, particle, alpham)) {
1471 std::cerr << m_className << "::IntegrateAlpha:\n"
1472 << " Cannot retrieve Townsend coefficient at mid point.\n";
1473 break;
1474 }
1475 // Compare first and second order estimates.
1476 if (d * fabs(alpha0 - 2 * alpham + alpha1) / 3. < tol) {
1477 // Accuracy is good enough.
1478 integral += d * (alpha0 + 4 * alpham + alpha1) / 6.;
1479 // Proceed to the next step.
1480 x0 = x1;
1481 x1 = xe;
1482 alpha0 = alpha1;
1483 } else {
1484 // Accuracy is not good enough, so halve the step.
1485 x1 = xm;
1486 alpha1 = alpham;
1487 }
1488 }
1489 return integral;
1490}
1491
1492double DriftLineRKF::IntegrateEta(const std::array<double, 3>& xi,
1493 const std::array<double, 3>& xe,
1494 const Particle particle, const double tol) {
1495
1496 double ex = 0., ey = 0., ez = 0.;
1497 double bx = 0., by = 0., bz = 0.;
1498 Medium* medium = nullptr;
1499
1500 // Make sure the starting point is valid.
1501 Vec x0 = xi;
1502 if (GetField(x0, ex, ey, ez, bx, by, bz, medium) != 0) {
1503 std::cerr << m_className << "::IntegrateEta: Invalid starting point "
1504 << PrintVec(xi) << ".\n";
1505 return 0.;
1506 }
1507 // Determine the attachment coefficient at the initial point.
1508 double eta0 = 0.;
1509 if (!GetEta(ex, ey, ez, bx, by, bz, medium, particle, eta0)) {
1510 std::cerr << m_className << "::IntegrateEta:\n"
1511 << " Cannot retrieve att. coefficient at initial point.\n";
1512 return 0.;
1513 }
1514 // Make sure the end point is valid.
1515 Vec x1 = xe;
1516 if (GetField(x1, ex, ey, ez, bx, by, bz, medium) != 0) {
1517 std::cerr << m_className << "::IntegrateEta: Invalid end point "
1518 << PrintVec(xe) << ".\n";
1519 return 0.;
1520 }
1521 // Determine the attachment coefficient at the end point.
1522 double eta1 = 0.;
1523 if (!GetEta(ex, ey, ez, bx, by, bz, medium, particle, eta1)) {
1524 std::cerr << m_className << "::IntegrateEta:\n"
1525 << " Cannot retrieve att. coefficient at end point.\n";
1526 return 0.;
1527 }
1528 double integral = 0.;
1529 bool done = false;
1530 while (!done) {
1531 // Check if we are close to the end point.
1532 if (Mag(xe[0] - x0[0], xe[1] - x0[1], xe[2] - x0[2]) < 1.e-6) {
1533 done = true;
1534 break;
1535 }
1536 const double d = Mag(x1[0] - x0[0], x1[1] - x0[1], x1[2] - x0[2]);
1537 if (d < 1.e-6) {
1538 // Step length has become very small.
1539 if (m_debug) {
1540 std::cout << m_className << "::IntegrateEta: Small step.\n";
1541 }
1542 integral += eta0 * d;
1543 // Proceed with the next step.
1544 x0 = x1;
1545 x1 = xe;
1546 continue;
1547 }
1548 // Calculate the attachment coefficient at the end point of the step.
1549 if (GetField(x1, ex, ey, ez, bx, by, bz, medium) != 0) {
1550 std::cerr << m_className << "::IntegrateEta: Invalid end point.\n";
1551 break;
1552 }
1553 if (!GetEta(ex, ey, ez, bx, by, bz, medium, particle, eta1)) {
1554 std::cerr << m_className << "::IntegrateEta:\n"
1555 << " Cannot retrieve att. coefficient at end point.\n";
1556 break;
1557 }
1558 // Calculate the attachment coefficient at the mid point of the step.
1559 Vec xm;
1560 for (unsigned int i = 0; i < 3; ++i) xm[i] = 0.5 * (x0[i] + x1[i]);
1561 if (GetField(xm, ex, ey, ez, bx, by, bz, medium) != 0) {
1562 std::cerr << m_className << "::IntegrateEta: Invalid mid point.\n";
1563 break;
1564 }
1565 double etam = 0.;
1566 if (!GetEta(ex, ey, ez, bx, by, bz, medium, particle, etam)) {
1567 std::cerr << m_className << "::IntegrateEta:\n"
1568 << " Cannot retrieve att. coefficient at mid point.\n";
1569 break;
1570 }
1571 // Compare first and second order estimates.
1572 if (d * fabs(eta0 - 2 * etam + eta1) / 3. < tol) {
1573 // Accuracy is good enough.
1574 integral += d * (eta0 + 4 * etam + eta1) / 6.;
1575 // Proceed to the next step.
1576 x0 = x1;
1577 x1 = xe;
1578 eta0 = eta1;
1579 } else {
1580 // Accuracy is not good enough, so halve the step.
1581 x1 = xm;
1582 eta1 = etam;
1583 }
1584 }
1585 return integral;
1586}
1587
1588void DriftLineRKF::ComputeSignal(const Particle particle, const double scale,
1589 const std::vector<double>& ts,
1590 const std::vector<Vec>& xs,
1591 const std::vector<double>& ne) const {
1592
1593 const unsigned int nPoints = ts.size();
1594 if (nPoints < 2) return;
1595 const double q = particle == Particle::Electron ? -1 * scale : scale;
1596
1597 // Get the drift velocity at each point.
1598 std::vector<std::array<double, 3> > vs;
1599 for (const auto& x : xs) {
1600 int stat = 0;
1601 Vec v;
1602 if (!GetVelocity(x, particle, v, stat)) {
1603 std::cerr << m_className << "::ComputeSignal:\n"
1604 << " Cannot retrieve velocity at " << PrintVec(x) << "\n";
1605 }
1606 vs.push_back(std::move(v));
1607 }
1608 m_sensor->AddSignal(q, ts, xs, vs, ne, m_navg);
1609}
1610
1611}
void DisablePlotting()
Switch off drift line plotting.
Definition: DriftLineRKF.cc:85
bool DriftNegativeIon(const double x0, const double y0, const double z0, const double t0)
double GetLoss(const double eps=1.e-4)
Compute the attachment loss factor for the current drift line.
bool DriftHole(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of a hole with a given starting point.
bool DriftIon(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an ion with a given starting point.
void SetSensor(Sensor *s)
Set the sensor.
Definition: DriftLineRKF.cc:50
void SetIntegrationAccuracy(const double eps)
Set the accuracy of the Runge Kutta Fehlberg drift line integration.
Definition: DriftLineRKF.cc:58
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an electron with a given starting point.
void SetGainFluctuationsPolya(const double theta, const double mean=-1.)
double GetGain(const double eps=1.e-4)
Compute the multiplication factor for the current drift line.
void SetGainFluctuationsFixed(const double gain=-1.)
Do not randomize the avalanche size.
Definition: DriftLineRKF.cc:87
bool DriftPositron(const double x0, const double y0, const double z0, const double t0)
void GetDriftLinePoint(const unsigned int i, double &x, double &y, double &z, double &t) const
Get the coordinates and time of a point along the most recent drift line.
void PrintDriftLine() const
Print the trajectory of the most recent drift line.
void EnablePlotting(ViewDrift *view)
Switch on drift line plotting.
Definition: DriftLineRKF.cc:77
DriftLineRKF()
Constructor.
Definition: DriftLineRKF.cc:45
double GetArrivalTimeSpread(const double eps=1.e-4)
void GetEndPoint(double &x, double &y, double &z, double &t, int &st) const
Get the end point and status flag of the most recent drift line.
void SetMaximumStepSize(const double ms)
Set the maximum step size that is allowed. By default, there is no limit.
Definition: DriftLineRKF.cc:67
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:124
bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
Definition: Sensor.cc:293
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:75
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:409
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:280
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
Definition: Sensor.cc:232
Visualize drift lines and tracks.
Definition: ViewDrift.hh:19
void NewElectronDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:74
void AddDriftLinePoint(const unsigned int iL, const double x, const double y, const double z)
Definition: ViewDrift.cc:174
void NewHoleDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:99
void NewIonDriftLine(const unsigned int np, int &id, const double x0, const double y0, const double z0)
Definition: ViewDrift.cc:119
std::array< double, 3 > Vec
Definition: DriftLineRKF.cc:43
double RndmPolya(const double theta)
Draw a Polya distributed random number.
Definition: Random.hh:70
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314