Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
TrackPAI.cc
Go to the documentation of this file.
1#include <iostream>
2#include <fstream>
3#include <cmath>
4
5#include "Sensor.hh"
6#include "TrackPAI.hh"
9#include "Random.hh"
10
11namespace Garfield {
12
14 : m_ready(false),
15 m_x(0.),
16 m_y(0.),
17 m_z(0.),
18 m_t(0.),
19 m_dx(0.),
20 m_dy(0),
21 m_dz(1.),
22 m_e(0.),
23 m_speed(0.),
24 m_emax(0.),
25 m_imfp(0.),
26 m_dedx(0.),
27 m_nSteps(1000),
28 m_mediumName(""),
29 m_mediumDensity(0.),
30 m_electronDensity(0.) {
31
32 m_className = "TrackPAI";
33
34}
35
36bool TrackPAI::NewTrack(const double x0, const double y0, const double z0,
37 const double t0, const double dx0, const double dy0,
38 const double dz0) {
39
40 m_ready = false;
41
42 // Make sure the sensor has been set.
43 if (!m_sensor) {
44 std::cerr << m_className << "::NewTrack: Sensor is not defined.\n";
45 return false;
46 }
47
48 // Get the medium at this location and check if it is "ionisable".
49 Medium* medium = NULL;
50 if (!m_sensor->GetMedium(x0, y0, z0, medium)) {
51 std::cerr << m_className << "::NewTrack: No medium at initial position.\n";
52 return false;
53 }
54 if (!medium->IsIonisable()) {
55 std::cerr << m_className << "::NewTrack:\n"
56 << " Medium at initial position is not ionisable.\n";
57 return false;
58 }
59
60 if (medium->GetName() != m_mediumName ||
61 medium->GetNumberDensity() != m_mediumDensity) {
62 m_isChanged = true;
63 if (!SetupMedium(medium)) {
64 std::cerr << m_className << "::NewTrack:\n Properties of medium "
65 << medium->GetName() << " are not available.\n";
66 return false;
67 }
68 m_mediumName = medium->GetName();
69 m_mediumDensity = medium->GetNumberDensity();
70 }
71
72 m_ready = true;
73
74 if (m_isChanged) {
75 if (!SetupCrossSectionTable()) {
76 std::cerr << m_className << "::NewTrack:\n"
77 << " Calculation of ionisation cross-section failed.\n";
78 m_ready = false;
79 return false;
80 }
81 m_isChanged = false;
82 }
83
84 m_x = x0;
85 m_y = y0;
86 m_z = z0;
87 m_t = t0;
88 const double d = sqrt(dx0 * dx0 + dy0 * dy0 + dz0 * dz0);
89 if (d < Small) {
90 if (m_debug) {
91 std::cout << m_className << "::NewTrack:\n"
92 << " Direction vector has zero norm.\n"
93 << " Initial direction is randomized.\n";
94 }
95 RndmDirection(m_dx, m_dy, m_dz);
96 } else {
97 // Normalize the direction vector.
98 m_dx = dx0 / d;
99 m_dy = dy0 / d;
100 m_dz = dz0 / d;
101 }
102 return true;
103}
104
105bool TrackPAI::GetCluster(double& xcls, double& ycls, double& zcls,
106 double& tcls, int& ncls, double& edep,
107 double& extra) {
108
109 ncls = 0;
110 edep = extra = 0.;
111
112 // Clear the stack.
113 m_electrons.clear();
114 m_holes.clear();
115
116 if (!m_ready) {
117 std::cerr << m_className << "::GetCluster:\n";
118 std::cerr << " Track not initialized. Call NewTrack first.\n";
119 return false;
120 }
121
122 if (m_isChanged) {
123 if (SetupCrossSectionTable()) {
124 m_isChanged = false;
125 } else {
126 std::cerr << m_className << "::GetCluster:\n";
127 std::cerr << " Calculation of ionisation cross-section failed.\n";
128 return false;
129 }
130 }
131
132 // Draw a step length and propagate the particle.
133 const double d = -m_imfp * log(RndmUniformPos());
134 m_x += d * m_dx;
135 m_y += d * m_dy;
136 m_z += d * m_dz;
137 m_t += d / m_speed;
138
139 // Check the medium at this location.
140 Medium* medium = NULL;
141 if (!m_sensor->GetMedium(m_x, m_y, m_z, medium)) {
142 m_ready = false;
143 return false;
144 }
145 if (medium->GetName() != m_mediumName ||
146 medium->GetNumberDensity() != m_mediumDensity || !medium->IsIonisable()) {
147 m_ready = false;
148 return false;
149 }
150
151 // Check if the particle is still inside the drift area.
152 if (!m_sensor->IsInArea(m_x, m_y, m_z)) {
153 m_ready = false;
154 return false;
155 }
156
157 xcls = m_x;
158 ycls = m_y;
159 zcls = m_z;
160 tcls = m_t;
161
162 // Sample the energy deposition.
163 double f = 0.;
164 edep = SampleEnergyDeposit(RndmUniform(), f);
165 // Update the particle energy.
166 m_e -= edep;
167
168 // Number of electron/hole (or electron/ion pairs) produced.
169 ncls = 1;
170
171 if (m_debug) {
172 std::cout << m_className << "::GetCluster:\n";
173 std::cout << " Fraction of Rutherford scattering: " << f << "\n";
174 }
175 return true;
176}
177
178double TrackPAI::SampleEnergyDeposit(const double u, double& f) const {
179
180 if (u > m_cdf.back()) {
181 // Use the free-electron differential cross-section.
182 f = 1.;
183 return SampleAsymptoticCs(u);
184 }
185
186 if (u <= m_cdf[0]) return m_energies[0];
187 if (u >= 1.) return m_energies.back();
188
189 // Find the energy loss by interpolation
190 // from the cumulative distribution table
191 int iLow = 0;
192 int iUp = m_cdf.size();
193 while (iUp - iLow > 1) {
194 const int iMid = (iUp + iLow) >> 1;
195 if (u >= m_cdf[iMid]) {
196 iLow = iMid;
197 } else {
198 iUp = iMid;
199 }
200 }
201 const double e0 = m_energies[iLow];
202 const double e1 = m_energies[iUp];
203 const double c0 = m_cdf[iLow];
204 const double c1 = m_cdf[iUp];
205 const double r0 = m_rutherford[iLow];
206 const double r1 = m_rutherford[iUp];
207 if (e0 < 100.) {
208 const double edep = e0 + (u - c0) * (e1 - e0) / (c1 - c0);
209 f = r0 + (edep - e0) * (r1 - r0) / (e1 - e0);
210 return edep;
211 }
212 const double loge0 = log(e0);
213 const double loge1 = log(e1);
214 const double logc0 = log(c0);
215 const double logc1 = log(c1);
216 double edep = loge0 + (log(u) - logc0) * (loge1 - loge0) / (logc1 - logc0);
217 f = r0 + (log(edep) - loge0) * (r1 - r0) / (loge1 - loge0);
218 edep = exp(edep);
219 return edep;
220}
221
222bool TrackPAI::SetupMedium(Medium* medium) {
223
224 // Make sure that the medium is defined.
225 if (!medium) {
226 std::cerr << m_className << "::SetupMedium: Null pointer.\n";
227 return false;
228 }
229
230 // Get the density and effective Z.
231 m_electronDensity = medium->GetNumberDensity() * medium->GetAtomicNumber();
232 if (m_electronDensity <= 0.) {
233 std::cerr << m_className << "::SetupMedium:\n"
234 << " Unphysical electron density ("
235 << m_electronDensity << ")\n";
236 return false;
237 }
238
239 // Get the dielectric function.
240 double emin, emax;
241 if (!medium->GetOpticalDataRange(emin, emax)) {
242 std::cerr << m_className << "::SetupMedium:\n";
243 std::cerr << " Could not load optical data for medium " << m_mediumName
244 << ".\n";
245 return false;
246 }
247
248 // Make sure the minimum energy is positive.
249 if (emin < Small) emin = Small;
250
251 // Reset the arrays.
252 m_energies.clear();
253 m_opticalDataTable.clear();
254 opticalData newEpsilon;
255
256 // Use logarithmically spaced energy steps.
257 const double r = pow(emax / emin, 1. / double(m_nSteps));
258 double eps1, eps2;
259
260 double eC = 0.5 * emin * (1. + r);
261 for (int i = 0; i < m_nSteps; ++i) {
262 medium->GetDielectricFunction(eC, eps1, eps2);
263 newEpsilon.eps1 = eps1;
264 newEpsilon.eps2 = eps2;
265 m_opticalDataTable.push_back(newEpsilon);
266 m_energies.push_back(eC);
267 eC *= r;
268 }
269
270 // Compute the integral of loss function times energy.
271 m_opticalDataTable[0].integral = 0.;
272 double integral = 0.;
273 double f1 = m_energies[0] * LossFunction(m_opticalDataTable[0].eps1,
274 m_opticalDataTable[0].eps2);
275 double f2 = f1;
276 for (int i = 1; i < m_nSteps; ++i) {
277 f2 = m_energies[i] *
278 LossFunction(m_opticalDataTable[i].eps1, m_opticalDataTable[i].eps2);
279 const double eM = 0.5 * (m_energies[i - 1] + m_energies[i]);
280 medium->GetDielectricFunction(eM, eps1, eps2);
281 const double fM = eM * LossFunction(eps1, eps2);
282 // Simpson's rule
283 integral += (f1 + 4 * fM + f2) * (m_energies[i] - m_energies[i - 1]) / 6.;
284 m_opticalDataTable[i].integral = integral;
285 f1 = f2;
286 }
287
288 // Check the consistency of the optical data by means of the TRK sum rule
289 const double trk = 2 * Pi2 * FineStructureConstant * pow(HbarC, 3) *
290 m_electronDensity / ElectronMass;
291 if (fabs(integral - trk) > 0.2 * trk) {
292 std::cerr << m_className << "::SetupMedium:\n";
293 std::cerr << " Deviation from Thomas-Reiche-Kuhn sum rule by > 20%.\n";
294 std::cerr << " Optical data are probably incomplete or erroneous!\n";
295 }
296
297 return true;
298}
299
301
302 if (!m_ready) {
303 std::cerr << m_className << "::GetClusterDensity:\n";
304 std::cerr << " Track has not been initialized.\n";
305 return 0.;
306 }
307
308 if (m_isChanged) {
309 if (SetupCrossSectionTable()) {
310 m_isChanged = false;
311 } else {
312 std::cerr << m_className << "::GetClusterDensity:\n";
313 std::cerr << " Ionisation cross-section could not be calculated.\n";
314 return 0.;
315 }
316 }
317
318 return 1. / m_imfp;
319}
320
322
323 if (!m_ready) {
324 std::cerr << m_className << "::GetStoppingPower:\n";
325 std::cerr << " Track has not been initialised.\n";
326 return 0.;
327 }
328
329 if (m_isChanged) {
330 if (SetupCrossSectionTable()) {
331 m_isChanged = false;
332 } else {
333 std::cerr << m_className << "::GetStoppingPower:\n";
334 std::cerr << " Ionisation cross-section could not be calculated.\n";
335 return 0.;
336 }
337 }
338
339 return m_dedx;
340}
341
342bool TrackPAI::SetupCrossSectionTable() {
343
344 if (!m_ready) {
345 std::cerr << m_className << "::SetupCrossSectionTable:\n"
346 << " Medium not set up.\n";
347 return false;
348 }
349
350 const double c1 = 2. * Pi2 * FineStructureConstant * pow(HbarC, 3) *
351 m_electronDensity / ElectronMass;
352 const double c2 = m_q * m_q * FineStructureConstant / (m_beta2 * Pi * HbarC);
353
354 // Get the max. allowed energy transfer.
355 m_emax = ComputeMaxTransfer();
356
357 std::ofstream outfile;
358 if (m_debug) outfile.open("dcs.txt", std::ios::out);
359
360 // Compute the differential cross-section.
361 std::vector<double> dcs;
362 m_rutherford.clear();
363
364 for (int i = 0; i < m_nSteps; ++i) {
365 // Define shorthand variables for photon energy and dielectric function.
366 const double egamma = m_energies[i];
367 const double eps1 = m_opticalDataTable[i].eps1;
368 const double eps2 = m_opticalDataTable[i].eps2;
369 const double integral = m_opticalDataTable[i].integral;
370
371 // First, calculate the distant-collision terms.
372 double dcsLog = 0., dcsDensity = 0., dcsCher = 0.;
373 if (eps2 > 0.) {
374 // Normal case (loss function > 0).
375 const double lf = LossFunction(eps1, eps2);
376 // Non-relativistic logarithmic term.
377 dcsLog = lf * log(2 * ElectronMass * m_beta2 / egamma);
378 // Relativistic logarithmic term (density effect)
379 const double u = 1. - m_beta2 * eps1;
380 const double v = m_beta2 * eps2;
381 dcsDensity = -0.5 * lf * log(u * u + v * v);
382 // "Cerenkov" term
383 dcsCher =
384 (m_beta2 - eps1 / (eps1 * eps1 + eps2 * eps2)) * (HalfPi - atan(u / v));
385 } else if (eps1 > 1. / m_beta2) {
386 // Imaginary part is zero, only the Cerenkov term contributes.
387 dcsCher = Pi * (m_beta2 - 1. / eps1);
388 }
389
390 // Calculate the close-collision term (quasi-free scattering)
391 double dcsRuth = 0.;
392 double f = 0.;
393 if (egamma > 0. && integral > 0.) {
394 dcsRuth = integral / (egamma * egamma);
395 f = dcsRuth / (dcsLog + dcsDensity + dcsCher);
396 }
397 m_rutherford.push_back(f);
398 dcs.push_back(dcsLog + dcsDensity + dcsCher + dcsRuth);
399 // If requested, write the cross-section terms to file.
400 if (m_debug) {
401 outfile << egamma << " " << eps1 << " " << eps2 << " " << dcsLog* c2
402 << " " << dcsDensity* c2 << " " << dcsCher* c2 << " "
403 << dcsRuth* c2 << "\n";
404 }
405 }
406 if (m_debug) outfile.close();
407
408 // Compute the cumulative distribution,
409 // total cross-section and stopping power.
410 m_cdf.clear();
411 m_cdf.push_back(0.);
412 m_dedx = 0.;
413 double cs = 0.;
414 for (int i = 1; i < m_nSteps; ++i) {
415 const double e0 = m_energies[i - 1];
416 const double e1 = m_energies[i];
417 const double de = e1 - e0;
418 const double dcs0 = dcs[i - 1];
419 const double dcs1 = dcs[i];
420 cs += 0.5 * (dcs0 + dcs1) * de;
421 m_cdf.push_back(cs);
422 m_dedx += 0.5 * (dcs0 * e0 + dcs1 * e1) * de;
423 }
424
425 // Add the contribution of high energy transfers to the stopping power
426 // and the total cross-section
427 const double elim = m_energies.back();
428 if (elim < m_emax) {
429 cs += c1 * ComputeCsTail(elim, m_emax);
430 m_dedx += c1 * ComputeDeDxTail(elim, m_emax);
431 } else {
432 std::cerr << m_className << "::SetupCrossSectionTable:\n";
433 std::cerr << " Max. energy transfer lower than optical data range.\n";
434 }
435
436 if (cs <= 0.) {
437 std::cerr << m_className << "::SetupCrossSectionTable:\n";
438 std::cerr << " Total cross-section <= 0.\n";
439 return false;
440 }
441
442 // Normalise the cumulative distribution.
443 for (int i = m_nSteps; i--;) m_cdf[i] /= cs;
444
445 cs *= c2;
446 m_dedx *= c2;
447
448 // Compute the inelastic mean free path
449 m_imfp = 1. / cs;
450
451 return true;
452}
453
454double TrackPAI::ComputeMaxTransfer() const {
455
456 if (m_isElectron) {
457 // Max. transfer for electrons is half the kinetic energy.
458 return 0.5 * (m_energy - m_mass);
459 }
460
461 // Other particles.
462 const double bg2 = m_beta2 / (1. - m_beta2);
463 const double mass2 = m_mass * m_mass;
464
465 return 2. * mass2 * ElectronMass * bg2 /
466 (mass2 + ElectronMass * ElectronMass + 2. * m_energy * ElectronMass);
467}
468
469double TrackPAI::ComputeCsTail(const double emin, const double emax) {
470
471 if (m_isElectron) {
472 // Electrons
473 const double ek = m_energy - m_mass;
474 return 1. / emin - 1. / emax - 2 * emin / (ek * ek) -
475 emin * emin / ((ek - emin) * ek * ek);
476 } else if (m_mass == ElectronMass) {
477 // Positrons
478 const double ek = m_energy - m_mass;
479 return 1. / emin - 1. / emax + 3 * (emax - emin) / (ek * ek) -
480 (emax - emin) * (ek * (emax + emin) +
481 (emin * emin + emin * emax + emax * emax) / 3.) /
482 pow(ek, 4) -
483 (2. / ek) * log(emax / emin);
484 }
485
486 switch (m_spin) {
487 case 0:
488 // Spin 0
489 return 1. / emin - 1. / emax - m_beta2 * log(emax / emin) / emax;
490 break;
491 case 1:
492 // Spin 1/2
493 return 1. / emin - 1. / emax - m_beta2 * log(emax / emin) / emax +
494 (emax - emin) / (2 * m_energy * m_energy);
495 break;
496 case 2: {
497 // Spin 1
498 const double e2 = 2 * m_energy * m_energy;
499 const double ec = m_mass * m_mass / ElectronMass;
500 const double a = 1. / (3 * ec);
501 const double b = (emax - emin);
502 return 1. / emin - 1. / emax + a * b * (emin + e2 + emax) / e2 -
503 m_beta2 * a * b / emax + (a - m_beta2 / emax) * log(emax / emin);
504 break;
505 }
506 default:
507 break;
508 }
509 // Rutherford type cross-section
510 return 1. / emin - 1. / emax;
511}
512
513double TrackPAI::ComputeDeDxTail(const double emin, const double emax) {
514
515 if (m_isElectron) {
516 const double ek = m_energy - m_mass;
517 return -log(emin * (ek - emin) / (ek * ek)) +
518 (1. / (8 * (emin - ek) * ek * ek)) *
519 (-4 * pow(emin, 3) + 4 * emin * emin * ek +
520 emin * ek * ek * (17. - 16. * CLog2) +
521 pow(ek, 3) * (-9. + 16. * CLog2));
522 } else if (m_mass == ElectronMass) {
523 // Positron
524 const double ek = m_energy - m_mass;
525 return log(ek / emin) -
526 (ek - emin) * (ek - emin) *
527 (3. * emin * emin - 2. * emin * ek + 11. * ek * ek) /
528 (12. * pow(ek, 4));
529 }
530
531 switch (m_spin) {
532 case 0:
533 return log(emax / emin) - m_beta2 * (emax - emin) / emax;
534 break;
535 case 1:
536 // Spin 1/2
537 return log(emax / emin) - m_beta2 * (emax - emin) / emax +
538 (emax * emax - emin * emin) / (2 * m_energy * m_energy);
539 break;
540 case 2: {
541 // Spin 1
542 double e2 = m_energy * m_energy;
543 double ec = m_mass * m_mass / ElectronMass;
544 return log(emax / emin) + (pow(emax, 3) - pow(emin, 3)) / (9. * e2 * ec) +
545 (emax * emax - emin * emin) / (6. * e2) + (emax - emin) *
546 (2. - (1. + emin / emax + 6 * ec / emax) * m_beta2) / (6. * ec);
547 break;
548 }
549 default:
550 break;
551 }
552
553 // Rutherford cross-section
554 return log(emax / emin);
555}
556
557double TrackPAI::SampleAsymptoticCs(double u) const {
558
559 const double emin = m_energies.back();
560 // Rescale the random number
561 u = (u - m_cdf.back()) / (1. - m_cdf.back());
562
563 if (m_isElectron) {
564 return SampleAsymptoticCsElectron(emin, u);
565 } else if (m_mass == ElectronMass) {
566 return SampleAsymptoticCsPositron(emin, u);
567 }
568
569 switch (m_spin) {
570 case 0:
571 // Spin 0
572 return SampleAsymptoticCsSpinZero(emin, u);
573 break;
574 case 1:
575 // Spin 1/2
576 return SampleAsymptoticCsSpinHalf(emin, u);
577 break;
578 case 2:
579 // Spin 1
580 return SampleAsymptoticCsSpinOne(emin, u);
581 break;
582 default:
583 break;
584 }
585 // Rutherford cross-section (analytic inversion)
586 return emin * m_emax / (u * emin + (1. - u) * m_emax);
587}
588
589double TrackPAI::SampleAsymptoticCsSpinZero(const double emin, double u) const {
590
591 const double a = emin / m_emax;
592 const double b = m_beta2 * a;
593 u *= (1. - a + b * log(a));
594 double eLow = emin, eUp = m_emax;
595 while (eUp - eLow > 1.) {
596 const double eM = 0.5 * (eUp + eLow);
597 if (u >= 1. - emin / eM - b * log(eM / emin)) {
598 eLow = eM;
599 } else {
600 eUp = eM;
601 }
602 }
603
604 return 0.5 * (eLow + eUp);
605}
606
607double TrackPAI::SampleAsymptoticCsSpinHalf(const double emin, double u) const {
608
609 const double a = emin / m_emax;
610 const double b = m_beta2 * a;
611 const double c = emin / (2. * m_energy * m_energy);
612 u *= 1. - a + b * log(a) + (m_emax - emin) * c;
613 double eLow = emin, eUp = m_emax;
614 while (eUp - eLow > 1.) {
615 const double eM = 0.5 * (eUp + eLow);
616 if (u >= 1. - emin / eM - b * log(eM / emin) + (eM - emin) * c) {
617 eLow = eM;
618 } else {
619 eUp = eM;
620 }
621 }
622
623 return 0.5 * (eLow + eUp);
624}
625
626double TrackPAI::SampleAsymptoticCsSpinOne(const double emin, double u) const {
627
628 const double e2 = 2 * m_energy * m_energy;
629 const double ec = m_mass * m_mass / ElectronMass;
630 const double a = 2 * ec / e2 - m_beta2 / m_emax;
631 const double b = 1.5 * ec / emin;
632 const double c = 1. - 1.5 * ec * m_beta2 / m_emax;
633 u *= (m_emax - emin) * (0.5 * (emin + m_emax) / e2 + a + b / m_emax) +
634 c * log(m_emax / emin);
635 double eLow = emin, eUp = m_emax;
636 while (eUp - eLow > 1.) {
637 const double eM = 0.5 * (eUp + eLow);
638 if (u >=
639 (eM - emin) * ((emin + eM) / e2 + a + b / eM) + c * log(eM / emin)) {
640 eLow = eM;
641 } else {
642 eUp = eM;
643 }
644 }
645
646 return 0.5 * (eLow + eUp);
647}
648
649double TrackPAI::SampleAsymptoticCsElectron(const double emin, double u) const {
650
651 const double ek = m_energy - m_mass;
652 const double ek2 = ek * ek;
653 const double a = ek / (emin * (ek - emin));
654 const double norm = 1. / emin - 0.5 / ek - emin * emin / ((ek - emin) * ek2) -
655 2. * emin / ek2;
656 u *= norm;
657 double eLow = emin, eUp = m_emax;
658 while (eUp - eLow > 1.) {
659 const double eM = 0.5 * (eUp + eLow);
660 if (u >= a - 1. / eM + (eM - emin) / ek2 + 1. / (ek - eM)) {
661 eLow = eM;
662 } else {
663 eUp = eM;
664 }
665 }
666 return 0.5 * (eLow + eUp);
667}
668
669double TrackPAI::SampleAsymptoticCsPositron(const double emin, double u) const {
670
671 const double ek = m_energy - m_mass;
672 const double ek2 = ek * ek;
673 const double ek3 = ek2 * ek;
674 const double ek4 = 3 * ek3 * ek;
675 const double emin2 = emin * emin;
676 const double a = 1. / emin;
677 const double b = 3. / ek2;
678 const double c = 2. / ek;
679 u *= 1. / emin - 1. / m_emax + 3 * (m_emax - emin) / ek2 -
680 (m_emax - emin) * (m_emax + emin) / ek3 +
681 (m_emax - emin) * (emin * emin + emin * m_emax + m_emax * m_emax) / ek4 -
682 (2. / ek) * log(m_emax / emin);
683 double eLow = emin, eUp = m_emax;
684 while (eUp - eLow > 1.) {
685 const double eM = 0.5 * (eUp + eLow);
686 const double eM2 = eM * eM;
687 if (u >= a - 1. / eM + b * (eM - emin) - (eM2 - emin2) / ek3 +
688 (eM - emin) * (emin2 + emin * eM + eM2) / ek4 -
689 c * log(eM / emin)) {
690 eLow = eM;
691 } else {
692 eUp = eM;
693 }
694 }
695
696 return 0.5 * (eLow + eUp);
697}
698
699}
Abstract base class for media.
Definition: Medium.hh:11
virtual double GetNumberDensity() const
Definition: Medium.hh:47
const std::string & GetName() const
Definition: Medium.hh:22
bool IsIonisable() const
Definition: Medium.hh:59
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:264
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:150
virtual bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &ncls, double &ecls, double &extra)
Definition: TrackPAI.cc:105
virtual double GetStoppingPower()
Get the stopping power (mean energy loss [eV] per cm).
Definition: TrackPAI.cc:321
virtual bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0)
Definition: TrackPAI.cc:36
virtual double GetClusterDensity()
Definition: TrackPAI.cc:300
Sensor * m_sensor
Definition: Track.hh:90
bool m_debug
Definition: Track.hh:97
bool m_isElectron
Definition: Track.hh:87
bool m_isChanged
Definition: Track.hh:92
double m_q
Definition: Track.hh:82
double m_beta2
Definition: Track.hh:86
std::string m_className
Definition: Track.hh:80
double m_mass
Definition: Track.hh:84
double m_energy
Definition: Track.hh:85
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition: Random.hh:14
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:106
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition: Random.hh:17
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:337
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615