Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
TrackSrim.cc
Go to the documentation of this file.
1#include <fstream>
2#include <iostream>
3#include <algorithm>
4
5#include <TCanvas.h>
6#include <TGraph.h>
7#include <TH1.h>
8#include <TLegend.h>
9#include <TLatex.h>
10
13#include "Garfield/Numerics.hh"
14#include "Garfield/Random.hh"
15#include "Garfield/Sensor.hh"
16#include "Garfield/ViewBase.hh"
17#include "Garfield/TrackSrim.hh"
18
19namespace {
20
21double StepVavilov(const double rkappa) {
22 double xlmin = -3.7;
23 if (rkappa < 0.1) {
24 xlmin = -2.7;
25 } else if (rkappa < 1) {
26 xlmin = -2.9;
27 } else if (rkappa < 2) {
28 xlmin = -3.0;
29 } else if (rkappa < 3) {
30 xlmin = -3.1;
31 } else if (rkappa < 4) {
32 xlmin = -3.2;
33 } else if (rkappa < 5) {
34 xlmin = -3.3;
35 } else if (rkappa < 6) {
36 xlmin = -3.4;
37 } else if (rkappa < 7) {
38 xlmin = -3.5;
39 } else if (rkappa < 8) {
40 xlmin = -3.6;
41 }
42 return xlmin;
43}
44
45double Interpolate(const double x, const std::vector<double>& xtab,
46 const std::vector<double>& ytab) {
47 if (x < xtab[0]) {
48 return ytab[0];
49 } else if (x > xtab.back()) {
50 return ytab.back();
51 }
52 return Garfield::Numerics::Divdif(ytab, xtab, xtab.size(), x, 2);
53}
54
55void PrintSettings(const std::string& hdr, const double de, const double step,
56 const double ekin, const double beta2, const double gamma,
57 const double agas, const double zgas, const double density,
58 const double qp, const double mp, const double emax,
59 const double xi, const double kappa) {
60 std::cout << hdr << "Settings:\n"
61 << " dE = " << de << " MeV,\n"
62 << " step = " << step << " cm.\n"
63 << " Ekin = " << ekin << " MeV,\n"
64 << " beta2 = " << beta2 << ",\n"
65 << " gamma = " << gamma << ".\n"
66 << " Agas = " << agas << ", Zgas = " << zgas << ",\n"
67 << " density = " << density << " g/cm3.\n"
68 << " Qpart = " << qp << ", mpart = " << 1.e-6 * mp << " MeV.\n"
69 << " Emax = " << emax << " MeV,\n"
70 << " xi = " << xi << " MeV,\n"
71 << " kappa = " << kappa << ".\n";
72}
73}
74
75namespace Garfield {
76
77TrackSrim::TrackSrim() : Track() { m_className = "TrackSrim"; }
78
79bool TrackSrim::ReadFile(const std::string& file) {
80 // SRMREA
81
82 const std::string hdr = m_className + "::ReadFile:\n ";
83 // Open the material list.
84 std::ifstream fsrim;
85 fsrim.open(file.c_str(), std::ios::in);
86 if (fsrim.fail()) {
87 std::cerr << hdr << "Could not open SRIM file " << file
88 << " for reading.\n The file perhaps does not exist.\n";
89 return false;
90 }
91 unsigned int nread = 0;
92
93 // Read the header
94 if (m_debug) {
95 std::cout << hdr << "SRIM header records from file " << file << "\n";
96 }
97 constexpr size_t size = 100;
98 char line[size];
99 while (fsrim.getline(line, size, '\n')) {
100 nread++;
101 if (strstr(line, "SRIM version") != NULL) {
102 if (m_debug) std::cout << "\t" << line << "\n";
103 } else if (strstr(line, "Calc. date") != NULL) {
104 if (m_debug) std::cout << "\t" << line << "\n";
105 } else if (strstr(line, "Ion =") != NULL) {
106 break;
107 }
108 }
109
110 // Identify the ion
111 char* token = NULL;
112 token = strtok(line, " []=");
113 token = strtok(NULL, " []=");
114 token = strtok(NULL, " []=");
115 // Set the ion charge.
116 m_qion = std::atof(token);
117 m_chargeset = true;
118 token = strtok(NULL, " []=");
119 token = strtok(NULL, " []=");
120 token = strtok(NULL, " []=");
121 // Set the ion mass (convert amu to eV).
122 m_mion = std::atof(token) * AtomicMassUnitElectronVolt;
123
124 // Find the target density
125 if (!fsrim.getline(line, size, '\n')) {
126 std::cerr << hdr << "Premature EOF looking for target density (line "
127 << nread << ").\n";
128 return false;
129 }
130 nread++;
131 if (!fsrim.getline(line, size, '\n')) {
132 std::cerr << hdr << "Premature EOF looking for target density (line "
133 << nread << ").\n";
134 return false;
135 }
136 nread++;
137 const bool pre2013 = (strstr(line, "Target Density") != NULL);
138 token = strtok(line, " ");
139 token = strtok(NULL, " ");
140 token = strtok(NULL, " ");
141 if (pre2013) token = strtok(NULL, " ");
142 SetDensity(std::atof(token));
143
144 // Check the stopping units
145 while (fsrim.getline(line, size, '\n')) {
146 nread++;
147 if (strstr(line, "Stopping Units") == NULL) continue;
148 if (strstr(line, "Stopping Units = MeV / (mg/cm2)") != NULL ||
149 strstr(line, "Stopping Units = MeV/(mg/cm2)") != NULL) {
150 if (m_debug) {
151 std::cout << hdr << "Stopping units: MeV / (mg/cm2) as expected.\n";
152 }
153 break;
154 }
155 std::cerr << hdr << "Unknown stopping units. Aborting (line " << nread
156 << ").\n";
157 return false;
158 }
159
160 // Skip to the table
161 while (fsrim.getline(line, size, '\n')) {
162 nread++;
163 if (strstr(line, "-----------") != NULL) break;
164 }
165
166 // Read the table line by line
167 m_ekin.clear();
168 m_emloss.clear();
169 m_hdloss.clear();
170 m_range.clear();
171 m_transstraggle.clear();
172 m_longstraggle.clear();
173 unsigned int ntable = 0;
174 while (fsrim.getline(line, size, '\n')) {
175 nread++;
176 if (strstr(line, "-----------") != NULL) break;
177 // Energy
178 token = strtok(line, " ");
179 m_ekin.push_back(atof(token));
180 token = strtok(NULL, " ");
181 if (strcmp(token, "eV") == 0) {
182 m_ekin[ntable] *= 1.0e-6;
183 } else if (strcmp(token, "keV") == 0) {
184 m_ekin[ntable] *= 1.0e-3;
185 } else if (strcmp(token, "GeV") == 0) {
186 m_ekin[ntable] *= 1.0e3;
187 } else if (strcmp(token, "MeV") != 0) {
188 std::cerr << hdr << "Unknown energy unit " << token << "; aborting\n";
189 return false;
190 }
191 // EM loss
192 token = strtok(NULL, " ");
193 m_emloss.push_back(atof(token));
194 // HD loss
195 token = strtok(NULL, " ");
196 m_hdloss.push_back(atof(token));
197 // Projected range
198 token = strtok(NULL, " ");
199 m_range.push_back(atof(token));
200 token = strtok(NULL, " ");
201 if (strcmp(token, "A") == 0) {
202 m_range[ntable] *= 1.0e-8;
203 } else if (strcmp(token, "um") == 0) {
204 m_range[ntable] *= 1.0e-4;
205 } else if (strcmp(token, "mm") == 0) {
206 m_range[ntable] *= 1.0e-1;
207 } else if (strcmp(token, "m") == 0) {
208 m_range[ntable] *= 1.0e2;
209 } else if (strcmp(token, "km") == 0) {
210 m_range[ntable] *= 1.0e5;
211 } else if (strcmp(token, "cm") != 0) {
212 std::cerr << hdr << "Unknown distance unit " << token << "; aborting\n";
213 return false;
214 }
215 // Longitudinal straggling
216 token = strtok(NULL, " ");
217 m_longstraggle.push_back(atof(token));
218 token = strtok(NULL, " ");
219 if (strcmp(token, "A") == 0) {
220 m_longstraggle[ntable] *= 1.0e-8;
221 } else if (strcmp(token, "um") == 0) {
222 m_longstraggle[ntable] *= 1.0e-4;
223 } else if (strcmp(token, "mm") == 0) {
224 m_longstraggle[ntable] *= 1.0e-1;
225 } else if (strcmp(token, "m") == 0) {
226 m_longstraggle[ntable] *= 1.0e2;
227 } else if (strcmp(token, "km") == 0) {
228 m_longstraggle[ntable] *= 1.0e5;
229 } else if (strcmp(token, "cm") != 0) {
230 std::cerr << hdr << "Unknown distance unit " << token << "; aborting\n";
231 return false;
232 }
233 // Transverse straggling
234 token = strtok(NULL, " ");
235 m_transstraggle.push_back(atof(token));
236 token = strtok(NULL, " ");
237 if (strcmp(token, "A") == 0) {
238 m_transstraggle[ntable] *= 1.0e-8;
239 } else if (strcmp(token, "um") == 0) {
240 m_transstraggle[ntable] *= 1.0e-4;
241 } else if (strcmp(token, "mm") == 0) {
242 m_transstraggle[ntable] *= 1.0e-1;
243 } else if (strcmp(token, "m") == 0) {
244 m_transstraggle[ntable] *= 1.0e2;
245 } else if (strcmp(token, "km") == 0) {
246 m_transstraggle[ntable] *= 1.0e5;
247 } else if (strcmp(token, "cm") != 0) {
248 std::cerr << hdr << "Unknown distance unit " << token << "; aborting\n";
249 return false;
250 }
251
252 // Increment table line counter
253 ++ntable;
254 }
255
256 // Find the scaling factor and convert to MeV/cm
257 double scale = -1.;
258 while (fsrim.getline(line, size, '\n')) {
259 nread++;
260 if (strstr(line, "=============") != NULL) {
261 break;
262 } else if (strstr(line, "MeV / (mg/cm2)") != NULL ||
263 strstr(line, "MeV/(mg/cm2)") != NULL) {
264 token = strtok(line, " ");
265 scale = std::atof(token);
266 }
267 }
268 if (scale < 0) {
269 std::cerr << hdr << "Did not find stopping unit scaling; aborting.\n";
270 return false;
271 }
272 scale *= 1.e3;
273 for (unsigned int i = 0; i < ntable; ++i) {
274 m_emloss[i] *= scale;
275 m_hdloss[i] *= scale;
276 }
277
278 // Seems to have worked
279 if (m_debug) {
280 std::cout << hdr << "Successfully read " << file << "(" << nread
281 << " lines).\n";
282 }
283 return true;
284}
285
287 std::cout << "TrackSrim::Print:\n SRIM energy loss table\n\n"
288 << " Energy EM Loss HD loss Range "
289 << "l straggle t straggle\n"
290 << " [MeV] [MeV/cm] [MeV/cm] [cm] "
291 << " [cm] [cm]\n\n";
292 const unsigned int nPoints = m_emloss.size();
293 for (unsigned int i = 0; i < nPoints; ++i) {
294 printf("%10g %10g %10g %10g %10g %10g\n", m_ekin[i],
297 }
298 std::cout << "\n";
299 printf(" Work function: %g eV\n", m_work);
300 printf(" Fano factor: %g\n", m_fano);
301 printf(" Ion charge: %g\n", m_qion);
302 printf(" Mass: %g MeV\n", 1.e-6 * m_mion);
303 printf(" Density: %g g/cm3\n", m_density);
304 printf(" A, Z: %g, %g\n", m_a, m_z);
305}
306
308
309 const unsigned int nPoints = m_ekin.size();
310 std::vector<double> yE;
311 std::vector<double> yH;
312 std::vector<double> yT;
313 for (unsigned int i = 0; i < nPoints; ++i) {
314 const double em = m_emloss[i] * m_density;
315 const double hd = m_hdloss[i] * m_density;
316 yE.push_back(em);
317 yH.push_back(hd);
318 yT.push_back(em + hd);
319 }
320 const double xmin = *std::min_element(std::begin(m_ekin), std::end(m_ekin));
321 const double xmax = *std::max_element(std::begin(m_ekin), std::end(m_ekin));
322 const double ymax = *std::max_element(std::begin(yT), std::end(yT));
323 // Prepare a plot frame.
324 const std::string name = ViewBase::FindUnusedCanvasName("cSRIM");
325 TCanvas* celoss = new TCanvas(name.c_str(), "Energy loss");
326 celoss->SetLogx();
327 celoss->SetGridx();
328 celoss->SetGridy();
329 celoss->DrawFrame(xmin, 0., xmax, 1.05 * ymax, ";Ion energy [MeV];Energy loss [MeV/cm]");
330
331 // Make a graph for the 3 curves to plot.
332 TGraph gr;
333 gr.SetLineStyle(kSolid);
334 gr.SetLineWidth(2);
335 gr.SetMarkerStyle(21);
336 gr.SetLineColor(kBlue + 1);
337 gr.SetMarkerColor(kBlue + 1);
338 gr.DrawGraph(nPoints, m_ekin.data(), yE.data(), "plsame");
339
340 gr.SetLineColor(kGreen + 2);
341 gr.SetMarkerColor(kGreen + 2);
342 gr.DrawGraph(nPoints, m_ekin.data(), yH.data(), "plsame");
343
344 gr.SetLineColor(kOrange - 3);
345 gr.SetMarkerColor(kOrange - 3);
346 gr.DrawGraph(nPoints, m_ekin.data(), yT.data(), "plsame");
347
348 TLatex label;
349 double xLabel = 0.4 * xmax;
350 double yLabel = 0.9 * ymax;
351 label.SetTextColor(kBlue + 1);
352 label.SetText(xLabel, yLabel, "EM energy loss");
353 label.DrawLatex(xLabel, yLabel, "EM energy loss");
354 yLabel -= 1.5 * label.GetYsize();
355 label.SetTextColor(kGreen + 2);
356 label.DrawLatex(xLabel, yLabel, "HD energy loss");
357 yLabel -= 1.5 * label.GetYsize();
358 label.SetTextColor(kOrange - 3);
359 label.DrawLatex(xLabel, yLabel, "Total energy loss");
360 celoss->Update();
361}
362
364
365 const double xmin = *std::min_element(std::begin(m_ekin), std::end(m_ekin));
366 const double xmax = *std::max_element(std::begin(m_ekin), std::end(m_ekin));
367 const double ymax = *std::max_element(std::begin(m_range), std::end(m_range));
368
369 // Prepare a plot frame.
370 const std::string name = ViewBase::FindUnusedCanvasName("cSRIM");
371 TCanvas* crange = new TCanvas(name.c_str(), "Range");
372 crange->SetLogx();
373 crange->SetGridx();
374 crange->SetGridy();
375 crange->DrawFrame(xmin, 0., xmax, 1.05 * ymax, ";Ion energy [MeV];Projected range [cm]");
376 // Make a graph.
377 TGraph gr;
378 gr.SetLineColor(kOrange - 3);
379 gr.SetMarkerColor(kOrange - 3);
380 gr.SetLineStyle(kSolid);
381 gr.SetLineWidth(2);
382 gr.SetMarkerStyle(21);
383 gr.DrawGraph(m_ekin.size(), m_ekin.data(), m_range.data(), "plsame");
384 crange->Update();
385}
386
388
389 const double xmin = *std::min_element(std::begin(m_ekin), std::end(m_ekin));
390 const double xmax = *std::max_element(std::begin(m_ekin), std::end(m_ekin));
391 const double ymax = std::max(*std::max_element(std::begin(m_longstraggle),
392 std::end(m_longstraggle)),
393 *std::max_element(std::begin(m_transstraggle),
394 std::end(m_transstraggle)));
395 // Prepare a plot frame.
396 const std::string name = ViewBase::FindUnusedCanvasName("cSRIM");
397 TCanvas* cstraggle = new TCanvas(name.c_str(), "Straggling");
398 cstraggle->SetLogx();
399 cstraggle->SetGridx();
400 cstraggle->SetGridy();
401 cstraggle->DrawFrame(xmin, 0., xmax, 1.05 * ymax, ";Ion energy [MeV];Straggling [cm]");
402
403 // Make a graph for the 2 curves to plot.
404 const unsigned int nPoints = m_ekin.size();
405 TGraph gr;
406 gr.SetLineStyle(kSolid);
407 gr.SetLineWidth(2);
408 gr.SetMarkerStyle(21);
409
410 gr.SetLineColor(kOrange - 3);
411 gr.SetMarkerColor(kOrange - 3);
412 gr.DrawGraph(nPoints, m_ekin.data(), m_longstraggle.data(), "plsame");
413
414 gr.SetLineColor(kGreen + 2);
415 gr.SetMarkerColor(kGreen + 2);
416 gr.DrawGraph(nPoints, m_ekin.data(), m_transstraggle.data(), "plsame");
417
418 TLatex label;
419 double xLabel = 1.2 * xmin;
420 double yLabel = 0.9 * ymax;
421 label.SetTextColor(kOrange - 3);
422 label.SetText(xLabel, yLabel, "Longitudinal");
423 label.DrawLatex(xLabel, yLabel, "Longitudinal");
424 yLabel -= 1.5 * label.GetYsize();
425 label.SetTextColor(kGreen + 2);
426 label.DrawLatex(xLabel, yLabel, "Transverse");
427 cstraggle->Update();
428}
429
430double TrackSrim::DedxEM(const double e) const {
431 return Interpolate(e, m_ekin, m_emloss);
432}
433
434double TrackSrim::DedxHD(const double e) const {
435 return Interpolate(e, m_ekin, m_hdloss);
436}
437
438double TrackSrim::Xi(const double x, const double beta2) const {
439
440 constexpr double fconst = 1.e-6 * TwoPi * (
441 FineStructureConstant * FineStructureConstant * HbarC * HbarC) /
442 (ElectronMass * AtomicMassUnit);
443 return fconst * m_qion * m_qion * m_z * m_density * x / (m_a * beta2);
444}
445
446bool TrackSrim::PreciseLoss(const double step, const double estart,
447 double& deem, double& dehd) const {
448 // SRMRKS
449
450 const std::string hdr = m_className + "::PreciseLoss: ";
451 // Debugging
452 if (m_debug) {
453 std::cout << hdr << "\n"
454 << " Initial energy: " << estart << " MeV\n"
455 << " Step: " << step << " cm\n";
456 }
457 // Precision aimed for.
458 const double eps = 1.0e-2;
459 // Number of intervals.
460 unsigned int ndiv = 1;
461 // Loop until precision achieved
462 const unsigned int nMaxIter = 10;
463 bool converged = false;
464 for (unsigned int iter = 0; iter < nMaxIter; ++iter) {
465 double e4 = estart;
466 double e2 = estart;
467 deem = 0.;
468 dehd = 0.;
469 // Compute rk2 and rk4 over the number of sub-divisions
470 const double s = m_density * step / ndiv;
471 for (unsigned int i = 0; i < ndiv; i++) {
472 // rk2: initial point
473 const double de21 = s * (DedxEM(e2) + DedxHD(e2));
474 // Mid-way point
475 const double em22 = s * DedxEM(e2 - 0.5 * de21);
476 const double hd22 = s * DedxHD(e2 - 0.5 * de21);
477 // Trace the rk2 energy
478 e2 -= em22 + hd22;
479 // rk4: initial point
480 const double em41 = s * DedxEM(e4);
481 const double hd41 = s * DedxHD(e4);
482 const double de41 = em41 + hd41;
483 // Mid-way point
484 const double em42 = s * DedxEM(e4 - 0.5 * de41);
485 const double hd42 = s * DedxHD(e4 - 0.5 * de41);
486 const double de42 = em42 + hd42;
487 // Second mid-point estimate
488 const double em43 = s * DedxEM(e4 - 0.5 * de42);
489 const double hd43 = s * DedxHD(e4 - 0.5 * de42);
490 const double de43 = em43 + hd43;
491 // End point estimate
492 const double em44 = s * DedxEM(e4 - de43);
493 const double hd44 = s * DedxHD(e4 - de43);
494 const double de44 = em44 + hd44;
495 // Store the energy loss terms (according to rk4)
496 deem += (em41 + em44) / 6. + (em42 + em43) / 3.;
497 dehd += (hd41 + hd44) / 6. + (hd42 + hd43) / 3.;
498 // Store the new energy computed with rk4
499 e4 -= (de41 + de44) / 6. + (de42 + de43) / 3.;
500 }
501 if (m_debug) {
502 std::cout << hdr << "\n Iteration " << iter << " has " << ndiv
503 << " division(s). Losses:\n";
504 printf("\tde4 = %12g, de2 = %12g MeV\n", estart - e2, estart - e4);
505 printf("\tem4 = %12g, hd4 = %12g MeV\n", deem, dehd);
506 }
507 // Compare the two estimates
508 if (fabs(e2 - e4) > eps * (fabs(e2) + fabs(e4) + fabs(estart))) {
509 // Repeat with twice the number of steps.
510 ndiv *= 2;
511 } else {
512 converged = true;
513 break;
514 }
515 }
516
517 if (!converged) {
518 std::cerr << hdr << "No convergence achieved integrating energy loss.\n";
519 } else if (m_debug) {
520 std::cout << hdr << "Convergence at eps = " << eps << "\n";
521 }
522 return converged;
523}
524
525bool TrackSrim::EstimateRange(const double ekin, const double step,
526 double& stpmax) {
527 // Find distance over which the ion just does not lose all its energy
528 // ekin : Kinetic energy [MeV]
529 // step : Step length as guessed [cm]
530 // stpmax : Maximum step
531 // SRMDEZ
532
533 const std::string hdr = m_className + "::EstimateRange: ";
534 // Initial estimate
535 stpmax = step;
536
537 // Find the energy loss expected for the present step length.
538 double st1 = step;
539 double deem = 0., dehd = 0.;
540 PreciseLoss(st1, ekin, deem, dehd);
541 double de1 = deem + dehd;
542 // Do nothing if this is ok
543 if (de1 < ekin) {
544 if (m_debug) std::cout << hdr << "Initial step OK.\n";
545 return true;
546 }
547 // Find a smaller step for which the energy loss is less than EKIN.
548 double st2 = 0.5 * step;
549 double de2 = de1;
550 const unsigned int nMaxIter = 20;
551 for (unsigned int iter = 0; iter < nMaxIter; ++iter) {
552 // See where we stand
553 PreciseLoss(st2, ekin, deem, dehd);
554 de2 = deem + dehd;
555 // Below the kinetic energy: done
556 if (de2 < ekin) break;
557 // Not yet below the kinetic energy: new iteration.
558 st1 = st2;
559 de1 = de2;
560 st2 *= 0.5;
561 }
562 if (de2 >= ekin) {
563 std::cerr << hdr << "\n Did not find a smaller step in " << nMaxIter
564 << " iterations. Abandoned.\n";
565 stpmax = 0.5 * (st1 + st2);
566 return false;
567 }
568 if (m_debug)
569 printf("\tstep 1 = %g cm, de 1 = %g MeV\n\tstep 2 = %g cm, de 2 = %g MeV\n",
570 st1, de1 - ekin, st2, de2 - ekin);
571
572 // Now perform a bisection
573 for (unsigned int iter = 0; iter < nMaxIter; ++iter) {
574 // Avoid division by zero.
575 if (de2 == de1) {
576 if (m_debug) {
577 std::cerr << hdr << "Bisection failed due to equal energy loss for "
578 << "two step sizes. Abandoned.\n";
579 }
580 stpmax = 0.5 * (st1 + st2);
581 return false;
582 }
583 // Estimate step to give total energy loss.
584 double st3;
585 if ((fabs(de1 - ekin) < 0.01 * fabs(de2 - de1)) ||
586 (fabs(de1 - ekin) > 0.99 * fabs(de2 - de1))) {
587 st3 = 0.5 * (st1 + st2);
588 } else {
589 st3 = st1 - (st2 - st1) * (de1 - ekin) / (de2 - de1);
590 }
591 // See how well we are doing.
592 PreciseLoss(st3, ekin, deem, dehd);
593 const double de3 = deem + dehd;
594 if (m_debug) {
595 std::printf("\tStep 1 = %g cm, dE 1 = %g MeV\n", st1, de1 - ekin);
596 std::printf("\tStep 2 = %g cm, dE 2 = %g MeV\n", st2, de2 - ekin);
597 std::printf("\tStep 3 = %g cm, dE 3 = %g MeV\n", st3, de3 - ekin);
598 }
599 // Update the estimates above and below.
600 if (de3 > ekin) {
601 st1 = st3;
602 de1 = de3;
603 } else {
604 st2 = st3;
605 de2 = de3;
606 }
607 // See whether we've converged.
608 if (fabs(de3 - ekin) < 1e-3 * (fabs(de3) + fabs(ekin)) ||
609 fabs(st1 - st2) < 1e-3 * (fabs(st1) + fabs(st2))) {
610 stpmax = st1 - (st2 - st1) * (de1 - ekin) / (de2 - de1);
611 return true;
612 }
613 }
614 if (m_debug) {
615 std::cout << hdr << "Bisection did not converge in " << nMaxIter
616 << " steps. Abandoned.\n";
617 }
618 stpmax = st1 - (st2 - st1) * (de1 - ekin) / (de2 - de1);
619 return false;
620}
621
622bool TrackSrim::NewTrack(const double x0, const double y0, const double z0,
623 const double t0, const double dx0, const double dy0,
624 const double dz0) {
625 // Generates electrons for a SRIM track
626 // SRMGEN
627 const std::string hdr = m_className + "::NewTrack: ";
628
629 // Verify that a sensor has been set.
630 if (!m_sensor) {
631 std::cerr << hdr << "Sensor is not defined.\n";
632 return false;
633 }
634
635 // Get the bounding box.
636 double xmin = 0., ymin = 0., zmin = 0.;
637 double xmax = 0., ymax = 0., zmax = 0.;
638 if (!m_sensor->GetArea(xmin, ymin, zmin, xmax, ymax, zmax)) {
639 std::cerr << hdr << "Drift area is not set.\n";
640 return false;
641 } else if (x0 < xmin || x0 > xmax || y0 < ymin || y0 > ymax || z0 < zmin ||
642 z0 > zmax) {
643 std::cerr << hdr << "Initial position outside bounding box.\n";
644 return false;
645 }
646
647 // Make sure the initial position is inside an ionisable medium.
648 Medium* medium = nullptr;
649 if (!m_sensor->GetMedium(x0, y0, z0, medium)) {
650 std::cerr << hdr << "No medium at initial position.\n";
651 return false;
652 } else if (!medium->IsIonisable()) {
653 std::cerr << hdr << "Medium at initial position is not ionisable.\n";
654 return false;
655 }
656
657 // Normalise and store the direction.
658 const double normdir = sqrt(dx0 * dx0 + dy0 * dy0 + dz0 * dz0);
659 double xdir = dx0;
660 double ydir = dy0;
661 double zdir = dz0;
662 if (normdir < Small) {
663 if (m_debug) {
664 std::cout << hdr << "Direction vector has zero norm.\n"
665 << " Initial direction is randomized.\n";
666 }
667 // Null vector. Sample the direction isotropically.
668 RndmDirection(xdir, ydir, zdir);
669 } else {
670 // Normalise the direction vector.
671 xdir /= normdir;
672 ydir /= normdir;
673 zdir /= normdir;
674 }
675
676 // Make sure all necessary parameters have been set.
677 if (m_mion < Small) {
678 std::cerr << hdr << "Particle mass not set.\n";
679 return false;
680 } else if (!m_chargeset) {
681 std::cerr << hdr << "Particle charge not set.\n";
682 return false;
683 } else if (m_energy < Small) {
684 std::cerr << hdr << "Initial particle energy not set.\n";
685 return false;
686 } else if (m_work < Small) {
687 std::cerr << hdr << "Work function not set.\n";
688 return false;
689 } else if (m_a < Small || m_z < Small) {
690 std::cerr << hdr << "A and/or Z not set.\n";
691 return false;
692 }
693 // Check the initial energy (in MeV).
694 const double ekin0 = 1.e-6 * GetKineticEnergy();
695 if (ekin0 < 1.e-14 * m_mion || ekin0 < 1.e-6 * m_work) {
696 std::cerr << hdr << "Initial kinetic energy E = " << ekin0
697 << " MeV such that beta2 = 0 or E << W; particle stopped.\n";
698 return true;
699 }
700
701 // Get an upper limit for the track length.
702 const double tracklength = 10 * Interpolate(ekin0, m_ekin, m_range);
703
704 // Header of debugging output.
705 if (m_debug) {
706 std::cout << hdr << "Track generation with the following parameters:\n";
707 const unsigned int nTable = m_ekin.size();
708 printf(" Table size %u\n", nTable);
709 printf(" Particle kin. energy %g MeV\n", ekin0);
710 printf(" Particle mass %g MeV\n", 1.e-6 * m_mion);
711 printf(" Particle charge %g\n", m_qion);
712 printf(" Work function %g eV\n", m_work);
713 if (m_fano > 0.) {
714 printf(" Fano factor %g\n", m_fano);
715 } else {
716 std::cout << " Fano factor Not set\n";
717 }
718 printf(" Long. straggling: %d\n", m_useLongStraggle);
719 printf(" Trans. straggling: %d\n", m_useTransStraggle);
720 printf(" Cluster size %d\n", m_nsize);
721 }
722
723 // Plot.
724 if (m_viewer) PlotNewTrack(x0, y0, z0);
725
726 // Reset the cluster count.
727 m_currcluster = 0;
728 m_clusters.clear();
729
730 // Initial situation: starting position
731 double x = x0;
732 double y = y0;
733 double z = z0;
734 double t = t0;
735 // Store the energy [MeV].
736 double e = ekin0;
737 // Total distance covered
738 double dsum = 0.0;
739 // Pool of unused energy
740 double epool = 0.0;
741
742 // Loop generating clusters
743 int iter = 0;
744 while (iter < m_maxclusters || m_maxclusters < 0) {
745 // Work out what the energy loss per cm, straggling and projected range are
746 // at the start of the step.
747 const double dedxem = DedxEM(e) * m_density;
748 const double dedxhd = DedxHD(e) * m_density;
749 const double prange = Interpolate(e, m_ekin, m_range);
750 double strlon = Interpolate(e, m_ekin, m_longstraggle);
751 double strlat = Interpolate(e, m_ekin, m_transstraggle);
752
753 if (!m_useLongStraggle) strlon = 0;
754 if (!m_useTransStraggle) strlat = 0;
755
756 if (m_debug) {
757 std::cout << hdr << "\n Energy = " << e
758 << " MeV,\n dEdx em, hd = " << dedxem << ", " << dedxhd
759 << " MeV/cm,\n e-/cm = " << 1.e6 * dedxem / m_work
760 << ".\n Straggling long/lat: " << strlon << ", " << strlat
761 << " cm\n";
762 }
763 // Find the step size for which we get approximately the target # clusters.
764 double step;
765 if (m_nsize > 0) {
766 step = m_nsize * 1.e-6 * m_work / dedxem;
767 } else {
768 const double ncls = m_maxclusters > 0 ? 0.5 * m_maxclusters : 100;
769 step = ekin0 / (ncls * (dedxem + dedxhd));
770 }
771 // Truncate if this step exceeds the length.
772 bool finish = false;
773 // Make an accurate integration of the energy loss over the step.
774 double deem = 0., dehd = 0.;
775 PreciseLoss(step, e, deem, dehd);
776 // If the energy loss exceeds the particle energy, truncate step.
777 double stpmax;
778 if (deem + dehd > e) {
779 EstimateRange(e, step, stpmax);
780 step = stpmax;
781 PreciseLoss(step, e, deem, dehd);
782 deem = e * deem / (dehd + deem);
783 dehd = e - deem;
784 finish = true;
785 if (m_debug) std::cout << hdr << "Finish raised. Track length reached.\n";
786 } else {
787 stpmax = tracklength - dsum;
788 }
789 if (m_debug) {
790 std::cout << hdr << "Maximum step size set to " << stpmax << " cm.\n";
791 }
792 // Ensure that this is larger than the minimum modelable step size.
793 double stpmin;
794 if (!SmallestStep(e, deem, step, stpmin)) {
795 std::cerr << hdr << "Failure computing the minimum step size."
796 << "\n Clustering abandoned.\n";
797 return false;
798 }
799
800 double eloss;
801 if (stpmin > stpmax) {
802 // No way to find a suitable step size: use fixed energy loss.
803 if (m_debug) std::cout << hdr << "stpmin > stpmax. Deposit all energy.\n";
804 eloss = deem;
805 if (e - eloss - dehd < 0) eloss = e - dehd;
806 finish = true;
807 if (m_debug) std::cout << hdr << "Finish raised. Single deposit.\n";
808 } else if (step < stpmin) {
809 // If needed enlarge the step size
810 if (m_debug) std::cout << hdr << "Enlarging step size.\n";
811 step = stpmin;
812 PreciseLoss(step, e, deem, dehd);
813 if (deem + dehd > e) {
814 if (m_debug) std::cout << hdr << "Excess loss. Recomputing stpmax.\n";
815 EstimateRange(e, step, stpmax);
816 step = stpmax;
817 PreciseLoss(step, e, deem, dehd);
818 deem = e * deem / (dehd + deem);
819 dehd = e - deem;
820 eloss = deem;
821 } else {
822 eloss = RndmEnergyLoss(e, deem, step);
823 }
824 } else {
825 // Draw an actual energy loss for such a step.
826 if (m_debug) std::cout << hdr << "Using existing step size.\n";
827 eloss = RndmEnergyLoss(e, deem, step);
828 }
829 // Ensure we are neither below 0 nor above the total energy.
830 if (eloss < 0) {
831 if (m_debug) std::cout << hdr << "Truncating negative energy loss.\n";
832 eloss = 0;
833 } else if (eloss > e - dehd) {
834 if (m_debug) std::cout << hdr << "Excess energy loss, using mean.\n";
835 eloss = deem;
836 if (e - eloss - dehd < 0) {
837 eloss = e - dehd;
838 finish = true;
839 if (m_debug) std::cout << hdr << "Finish raised. Using mean energy.\n";
840 }
841 }
842 if (m_debug) {
843 std::cout << hdr << "Step length = " << step << " cm.\n "
844 << "Mean loss = " << deem << " MeV.\n "
845 << "Actual loss = " << eloss << " MeV.\n";
846 }
847
848 // Check that the cluster is in an ionisable medium and within bounding box
849 if (!m_sensor->GetMedium(x, y, z, medium)) {
850 if (m_debug) {
851 std::cout << hdr << "No medium at position (" << x << "," << y << ","
852 << z << ").\n";
853 }
854 break;
855 } else if (!medium->IsIonisable()) {
856 if (m_debug) {
857 std::cout << hdr << "Medium at (" << x << "," << y << "," << z
858 << ") is not ionisable.\n";
859 }
860 break;
861 } else if (!m_sensor->IsInArea(x, y, z)) {
862 if (m_debug) {
863 std::cout << hdr << "Cluster at (" << x << "," << y << "," << z
864 << ") outside bounding box.\n";
865 }
866 break;
867 }
868 // Add a cluster.
869 Cluster cluster;
870 cluster.x = x;
871 cluster.y = y;
872 cluster.z = z;
873 cluster.t = t;
874 if (m_fano < Small) {
875 // No fluctuations.
876 cluster.electrons = int((eloss + epool) / (1.e-6 * m_work));
877 cluster.ec = m_work * cluster.electrons;
878 } else {
879 double ecl = 1.e6 * (eloss + epool);
880 cluster.electrons = 0;
881 cluster.ec = 0.0;
882 while (true) {
883 // if (cluster.ec < 100) printf("ec = %g\n", cluster.ec);
884 const double ernd1 = RndmHeedWF(m_work, m_fano);
885 if (ernd1 > ecl) break;
886 cluster.electrons++;
887 cluster.ec += ernd1;
888 ecl -= ernd1;
889 }
890 if (m_debug) {
891 std::cout << hdr << "EM + pool: " << 1.e6 * (eloss + epool)
892 << " eV, W: " << m_work
893 << " eV, E/w: " << (eloss + epool) / (1.0e-6 * m_work)
894 << ", n: " << cluster.electrons << ".\n";
895 }
896 }
897 cluster.kinetic = e;
898 epool += eloss - 1.e-6 * cluster.ec;
899 if (m_debug) {
900 std::cout << hdr << "Cluster " << m_clusters.size() << "\n at ("
901 << cluster.x << ", " << cluster.y << ", " << cluster.z
902 << "),\n e = " << cluster.ec
903 << ",\n n = " << cluster.electrons
904 << ",\n pool = " << epool << " MeV.\n";
905 }
906 m_clusters.push_back(std::move(cluster));
907 if (m_viewer) PlotCluster(x, y, z);
908
909 // Keep track of the length and energy
910 dsum += step;
911 e -= eloss + dehd;
912 if (finish) {
913 // Stop if the flag is raised
914 if (m_debug) std::cout << hdr << "Finishing flag raised.\n";
915 break;
916 } else if (e < ekin0 * 1.e-9) {
917 // No energy left
918 if (m_debug) std::cout << hdr << "Energy exhausted.\n";
919 break;
920 }
921 // Draw scattering distances
922 const double scale = sqrt(step / prange);
923 const double sigt1 = RndmGaussian(0., scale * strlat);
924 const double sigt2 = RndmGaussian(0., scale * strlat);
925 const double sigl = RndmGaussian(0., scale * strlon);
926 if (m_debug)
927 std::cout << hdr << "sigma l, t1, t2: " << sigl << ", " << sigt1 << ", "
928 << sigt2 << "\n";
929 // Rotation angles to bring z-axis in line
930 double theta, phi;
931 if (xdir * xdir + zdir * zdir <= 0) {
932 if (ydir < 0) {
933 theta = -HalfPi;
934 } else if (ydir > 0) {
935 theta = +HalfPi;
936 } else {
937 std::cerr << hdr << "Zero step length; clustering abandoned.\n";
938 return false;
939 }
940 phi = 0;
941 } else {
942 phi = atan2(xdir, zdir);
943 theta = atan2(ydir, sqrt(xdir * xdir + zdir * zdir));
944 }
945
946 // Update the position.
947 const double cp = cos(phi);
948 const double ct = cos(theta);
949 const double sp = sin(phi);
950 const double st = sin(theta);
951 x += step * xdir + cp * sigt1 - sp * st * sigt2 + sp * ct * sigl;
952 y += step * ydir + ct * sigt2 + st * sigl;
953 z += step * zdir - sp * sigt1 - cp * st * sigt2 + cp * ct * sigl;
954 // Update the time.
955 const double rk = 1.e6 * e / m_mion;
956 const double gamma = 1. + rk;
957 const double beta2 = rk > 1.e-5 ? 1. - 1. / (gamma * gamma) : 2. * rk;
958 const double vmag = sqrt(beta2) * SpeedOfLight;
959 if (vmag > 0.) t += step / vmag;
960 // (Do not) update direction
961 if (false) {
962 xdir = step * xdir + cp * sigt1 - sp * st * sigt2 + sp * ct * sigl;
963 ydir = step * ydir + ct * sigt2 + st * sigl;
964 zdir = step * zdir - sp * sigt1 - cp * st * sigt2 + cp * ct * sigl;
965 double dnorm = sqrt(xdir * xdir + ydir * ydir + zdir * zdir);
966 if (dnorm <= 0) {
967 std::cerr << hdr << "Zero step length; clustering abandoned.\n";
968 return false;
969 }
970 xdir = xdir / dnorm;
971 ydir = ydir / dnorm;
972 zdir = zdir / dnorm;
973 }
974 // Next cluster
975 iter++;
976 }
977 if (iter == m_maxclusters) {
978 std::cerr << hdr << "Exceeded maximum number of clusters.\n";
979 }
980 return true;
981 // finished generating
982}
983
984bool TrackSrim::SmallestStep(const double ekin, double de, double step,
985 double& stpmin) {
986 // Determines the smallest step size for which there is little
987 // or no risk of finding negative energy fluctuations.
988 // SRMMST
989
990 const std::string hdr = m_className + "::SmallestStep: ";
991 constexpr double expmax = 30;
992
993 // By default, assume the step is right.
994 stpmin = step;
995 // Check correctness.
996 if (ekin <= 0 || de <= 0 || step <= 0) {
997 std::cerr << hdr << "Input parameters not valid.\n Ekin = " << ekin
998 << " MeV, dE = " << de << " MeV, step length = " << step
999 << " cm.\n";
1000 return false;
1001 } else if (m_mion <= 0 || fabs(m_qion) <= 0) {
1002 std::cerr << hdr
1003 << "Track parameters not valid.\n Mass = " << 1.e-6 * m_mion
1004 << " MeV, charge = " << m_qion << ".\n";
1005 return false;
1006 } else if (m_a <= 0 || m_z <= 0 || m_density <= 0) {
1007 std::cerr << hdr << "Gas parameters not valid.\n A = " << m_a
1008 << ", Z = " << m_z << " density = " << m_density << " g/cm3.\n";
1009 return false;
1010 }
1011
1012 // Basic kinematic parameters
1013 const double rkin = 1.e6 * ekin / m_mion;
1014 const double gamma = 1. + rkin;
1015 const double beta2 = rkin > 1.e-5 ? 1. - 1. / (gamma * gamma) : 2. * rkin;
1016
1017 // Compute the maximum energy transfer [MeV]
1018 const double rm = ElectronMass / m_mion;
1019 const double emax = 2 * ElectronMass * 1.e-6 * beta2 * gamma * gamma /
1020 (1. + 2 * gamma * rm + rm * rm);
1021 // Compute the Rutherford term
1022 double xi = Xi(step, beta2);
1023 // Compute the scaling parameter
1024 double rkappa = xi / emax;
1025 // Step size and energy loss
1026 double denow = de;
1027 double stpnow = step;
1028 constexpr unsigned int nMaxIter = 10;
1029 for (unsigned int iter = 0; iter < nMaxIter; ++iter) {
1030 bool retry = false;
1031 // Debugging output.
1032 if (m_debug) {
1033 PrintSettings(hdr, denow, stpnow, ekin, beta2, gamma, m_a, m_z, m_density,
1034 m_qion, m_mion, emax, xi, rkappa);
1035 }
1036 double xinew = xi;
1037 double rknew = rkappa;
1038 if (m_model <= 0 || m_model > 4) {
1039 // No fluctuations: any step is permitted
1040 stpmin = stpnow;
1041 } else if (m_model == 1) {
1042 // Landau distribution
1043 constexpr double xlmin = -3.;
1044 const double exponent = -xlmin - 1. + Gamma - beta2 - denow / xi;
1045 const double rklim = exponent < -expmax ? 0. : exp(exponent);
1046 stpmin = stpnow * (rklim / rkappa);
1047 if (m_debug) {
1048 std::cout << hdr << "Landau distribution is imposed.\n kappa_min = "
1049 << rklim << ", d_min = " << stpmin << " cm.\n";
1050 }
1051 } else if (m_model == 2) {
1052 // Vavilov distribution, ensure we're in range.
1053 const double xlmin = StepVavilov(rkappa);
1054 const double exponent = -xlmin - 1. + Gamma - beta2 - denow / xi;
1055 const double rklim = exponent < -expmax ? 0. : exp(exponent);
1056 stpmin = stpnow * (rklim / rkappa);
1057 xinew = Xi(stpmin, beta2);
1058 rknew = xinew / emax;
1059 if (m_debug) {
1060 std::cout << hdr << "Vavilov distribution is imposed.\n kappa_min = "
1061 << rklim << ", d_min = " << stpmin
1062 << " cm\n kappa_new = " << rknew << ", xi_new = " << xinew
1063 << " MeV.\n";
1064 }
1065 if (stpmin > stpnow * 1.1) {
1066 if (m_debug) std::cout << hdr << "Step size increase. New pass.\n";
1067 retry = true;
1068 }
1069 } else if (m_model == 3) {
1070 // Gaussian model
1071 const double sigma2 = xi * emax * (1 - 0.5 * beta2);
1072 stpmin = stpnow * 16 * sigma2 / (denow * denow);
1073 if (m_debug) {
1074 const double sigmaMin2 = Xi(stpmin, beta2) * emax * (1 - 0.5 * beta2);
1075 std::cout << hdr << "Gaussian distribution is imposed.\n "
1076 << "d_min = " << stpmin << " cm.\n sigma/mu (old) = "
1077 << sqrt(sigma2) / de << ",\n sigma/mu (min) = "
1078 << sqrt(sigmaMin2) / (stpmin * denow / stpnow) << "\n";
1079 }
1080 } else if (rkappa < 0.05) {
1081 // Combined model: for low kappa, use the Landau distribution.
1082 constexpr double xlmin = -3.;
1083 const double exponent = -xlmin - 1. + Gamma - beta2 - denow / xi;
1084 const double rklim = exponent < -expmax ? 0. : exp(exponent);
1085 stpmin = stpnow * (rklim / rkappa);
1086 xinew = Xi(stpmin, beta2);
1087 rknew = xinew / emax;
1088 if (m_debug) {
1089 std::cout << hdr << "Landau distribution automatic.\n kappa_min = "
1090 << rklim << ", d_min = " << stpmin << " cm.\n";
1091 }
1092 if (rknew > 0.05 || stpmin > stpnow * 1.1) {
1093 retry = true;
1094 if (m_debug) {
1095 std::cout << hdr << "Model change or step increase. New pass.\n";
1096 }
1097 }
1098 } else if (rkappa < 5) {
1099 // For medium kappa, use the Vavilov distribution
1100 const double xlmin = StepVavilov(rkappa);
1101 const double exponent = -xlmin - 1. + Gamma - beta2 - denow / xi;
1102 const double rklim = exponent < -expmax ? 0. : exp(exponent);
1103 stpmin = stpnow * (rklim / rkappa);
1104 xinew = Xi(stpmin, beta2);
1105 rknew = xinew / emax;
1106 if (m_debug) {
1107 std::cout << hdr << "Vavilov distribution automatic.\n kappa_min = "
1108 << rklim << ", d_min = " << stpmin << " cm\n kappa_new = "
1109 << rknew << ", xi_new = " << xinew << " MeV.\n";
1110 }
1111 if (rknew > 5 || stpmin > stpnow * 1.1) {
1112 retry = true;
1113 if (m_debug) {
1114 std::cout << hdr << "Model change or step increase. New pass.\n";
1115 }
1116 }
1117 } else {
1118 // And for large kappa, use the Gaussian values.
1119 const double sigma2 = xi * emax * (1 - 0.5 * beta2);
1120 stpmin = stpnow * 16 * sigma2 / (denow * denow);
1121 if (m_debug) {
1122 const double sigmaMin2 = Xi(stpmin, beta2) * emax * (1 - 0.5 * beta2);
1123 std::cout << hdr << "Gaussian distribution automatic.\n "
1124 << "d_min = " << stpmin << " cm.\n sigma/mu (old) = "
1125 << sqrt(sigma2) / de << ",\n sigma/mu (min) = "
1126 << sqrt(sigmaMin2) / (stpmin * denow / stpnow) << "\n";
1127 }
1128 }
1129 // See whether we should do another pass.
1130 if (stpnow > stpmin) {
1131 if (m_debug) {
1132 std::cout << hdr << "Step size ok, minimum: " << stpmin << " cm\n";
1133 }
1134 break;
1135 }
1136 if (!retry) {
1137 if (m_debug) {
1138 std::cerr << hdr << "Step size must be increased to " << stpmin
1139 << "cm.\n";
1140 }
1141 break;
1142 }
1143 // New iteration
1144 rkappa = rknew;
1145 xi = xinew;
1146 denow *= stpmin / stpnow;
1147 stpnow = stpmin;
1148 if (m_debug) std::cout << hdr << "Iteration " << iter << "\n";
1149 if (iter == nMaxIter - 1) {
1150 // Need interation, but ran out of tries
1151 std::cerr << hdr << "No convergence reached on step size.\n";
1152 }
1153 }
1154 return true;
1155}
1156
1157double TrackSrim::RndmEnergyLoss(const double ekin, const double de,
1158 const double step) const {
1159 // RNDDE - Generates a random energy loss.
1160 // VARIABLES : EKIN : Kinetic energy [MeV]
1161 // DE : Mean energy loss over the step [MeV]
1162 // STEP : Step length [cm]
1163 // BETA2 : Velocity-squared
1164 // GAMMA : Projectile gamma
1165 // EMAX : Maximum energy transfer per collision [MeV]
1166 // XI : Rutherford term [MeV]
1167 // FCONST : Proportionality constant
1168 // EMASS : Electron mass [MeV]
1169
1170 const std::string hdr = "TrackSrim::RndmEnergyLoss: ";
1171 // Check correctness.
1172 if (ekin <= 0 || de <= 0 || step <= 0) {
1173 std::cerr << hdr << "Input parameters not valid.\n Ekin = " << ekin
1174 << " MeV, dE = " << de << " MeV, step length = " << step
1175 << " cm.\n";
1176 return 0.;
1177 } else if (m_mion <= 0 || fabs(m_qion) <= 0) {
1178 std::cerr << hdr << "Track parameters not valid.\n Mass = "
1179 << m_mion << " MeV, charge = " << m_qion << ".\n";
1180 return 0.;
1181 } else if (m_a <= 0 || m_z <= 0 || m_density <= 0) {
1182 std::cerr << hdr << "Material parameters not valid.\n A = " << m_a
1183 << ", Z = " << m_z << ", density = " << m_density << " g/cm3.\n";
1184 return 0.;
1185 }
1186 // Basic kinematic parameters
1187 const double rkin = 1.e6 * ekin / m_mion;
1188 const double gamma = 1. + rkin;
1189 const double beta2 = rkin > 1.e-5 ? 1. - 1. / (gamma * gamma) : 2. * rkin;
1190
1191 // Compute maximum energy transfer
1192 const double rm = ElectronMass / m_mion;
1193 const double emax = 2 * ElectronMass * 1.e-6 * beta2 * gamma * gamma /
1194 (1. + 2 * gamma * rm + rm * rm);
1195 // Compute the Rutherford term
1196 const double xi = Xi(step, beta2);
1197 // Compute the scaling parameter
1198 const double rkappa = xi / emax;
1199 // Debugging output.
1200 if (m_debug) {
1201 PrintSettings(hdr, de, step, ekin, beta2, gamma, m_a, m_z, m_density,
1202 m_qion, m_mion, emax, xi, rkappa);
1203 }
1204 double rndde = de;
1205 if (m_model <= 0 || m_model > 4) {
1206 // No fluctuations.
1207 if (m_debug) std::cout << hdr << "Fixed energy loss.\n";
1208 } else if (m_model == 1) {
1209 // Landau distribution
1210 if (m_debug) std::cout << hdr << "Landau imposed.\n";
1211 const double xlmean = -(log(rkappa) + beta2 + 1. - Gamma);
1212 rndde += xi * (RndmLandau() - xlmean);
1213 } else if (m_model == 2) {
1214 // Vavilov distribution, ensure we are in range.
1215 if (m_debug) std::cout << hdr << "Vavilov imposed.\n";
1216 if (rkappa > 0.01 && rkappa < 12) {
1217 const double xvav = RndmVavilov(rkappa, beta2);
1218 rndde += xi * (xvav + log(rkappa) + beta2 + (1 - Gamma));
1219 }
1220 } else if (m_model == 3) {
1221 // Gaussian model
1222 if (m_debug) std::cout << hdr << "Gaussian imposed.\n";
1223 rndde += RndmGaussian(0., sqrt(xi * emax * (1 - 0.5 * beta2)));
1224 } else if (rkappa < 0.05) {
1225 // Combined model: for low kappa, use the landau distribution.
1226 if (m_debug) std::cout << hdr << "Landau automatic.\n";
1227 const double xlmean = -(log(rkappa) + beta2 + (1 - Gamma));
1228 const double par[] = {0.50884, 1.26116, 0.0346688, 1.46314,
1229 0.15088e-2, 1.00324, -0.13049e-3};
1230 const double xlmax = par[0] + par[1] * xlmean + par[2] * xlmean * xlmean +
1231 par[6] * xlmean * xlmean * xlmean +
1232 (par[3] + xlmean * par[4]) * exp(par[5] * xlmean);
1233 double xlan = RndmLandau();
1234 for (unsigned int iter = 0; iter < 100; ++iter) {
1235 if (xlan < xlmax) break;
1236 xlan = RndmLandau();
1237 }
1238 rndde += xi * (xlan - xlmean);
1239 } else if (rkappa < 5) {
1240 // For medium kappa, use the Vavilov distribution.
1241 if (m_debug) std::cout << hdr << "Vavilov fast automatic.\n";
1242 const double xvav = RndmVavilov(rkappa, beta2);
1243 rndde += xi * (xvav + log(rkappa) + beta2 + (1 - Gamma));
1244 } else {
1245 // And for large kappa, use the Gaussian values.
1246 if (m_debug) std::cout << hdr << "Gaussian automatic.\n";
1247 rndde = RndmGaussian(de, sqrt(xi * emax * (1 - 0.5 * beta2)));
1248 }
1249 // Debugging output
1250 if (m_debug) {
1251 std::cout << hdr << "Energy loss generated = " << rndde << " MeV.\n";
1252 }
1253 return rndde;
1254}
1255
1256bool TrackSrim::GetCluster(double& xcls, double& ycls, double& zcls,
1257 double& tcls, int& n, double& e, double& extra) {
1258 if (m_debug) {
1259 std::cout << m_className << "::GetCluster: Cluster " << m_currcluster
1260 << " of " << m_clusters.size() << "\n";
1261 }
1262 // Stop if we have exhausted the list of clusters.
1263 if (m_currcluster >= m_clusters.size()) return false;
1264
1265 const auto& cluster = m_clusters[m_currcluster];
1266 xcls = cluster.x;
1267 ycls = cluster.y;
1268 zcls = cluster.z;
1269 tcls = cluster.t;
1270
1271 n = cluster.electrons;
1272 e = cluster.ec;
1273 extra = cluster.kinetic;
1274 // Move to next cluster
1275 ++m_currcluster;
1276 return true;
1277}
1278}
Abstract base class for media.
Definition: Medium.hh:13
bool IsIonisable() const
Is charge deposition by charged particles/photon enabled in this medium?
Definition: Medium.hh:78
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
Definition: Sensor.cc:254
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
Definition: Sensor.cc:159
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
Definition: Sensor.cc:228
std::vector< double > m_longstraggle
Longitudinal straggling [cm].
Definition: TrackSrim.hh:122
double RndmEnergyLoss(const double ekin, const double de, const double step) const
Definition: TrackSrim.cc:1157
void Print()
Print the energy loss table.
Definition: TrackSrim.cc:286
double m_qion
Charge of ion.
Definition: TrackSrim.hh:101
double m_z
Effective Z of the gas.
Definition: TrackSrim.hh:107
std::vector< double > m_hdloss
Hadronic energy loss [MeV cm2/g].
Definition: TrackSrim.hh:116
void SetDensity(const double density)
Set the density [g/cm3] of the target medium.
Definition: TrackSrim.hh:27
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra) override
Definition: TrackSrim.cc:1256
std::vector< double > m_range
Projected range [cm].
Definition: TrackSrim.hh:118
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
Definition: TrackSrim.cc:622
double DedxHD(const double e) const
Definition: TrackSrim.cc:434
std::vector< Cluster > m_clusters
Definition: TrackSrim.hh:137
size_t m_currcluster
Index of the next cluster to be returned.
Definition: TrackSrim.hh:125
std::vector< double > m_emloss
EM energy loss [MeV cm2/g].
Definition: TrackSrim.hh:114
bool m_chargeset
Charge gas been defined.
Definition: TrackSrim.hh:93
bool SmallestStep(double ekin, double de, double step, double &stpmin)
Definition: TrackSrim.cc:984
double m_fano
Fano factor [-].
Definition: TrackSrim.hh:99
TrackSrim()
Constructor.
Definition: TrackSrim.cc:77
double m_mion
Mass of ion [MeV].
Definition: TrackSrim.hh:103
bool PreciseLoss(const double step, const double estart, double &deem, double &dehd) const
Definition: TrackSrim.cc:446
double m_work
Work function [eV].
Definition: TrackSrim.hh:97
unsigned int m_model
Definition: TrackSrim.hh:128
double m_density
Density [g/cm3].
Definition: TrackSrim.hh:95
bool m_useTransStraggle
Include transverse straggling.
Definition: TrackSrim.hh:88
int m_nsize
Targeted cluster size.
Definition: TrackSrim.hh:130
std::vector< double > m_ekin
Energy in energy loss table [MeV].
Definition: TrackSrim.hh:112
std::vector< double > m_transstraggle
Transverse straggling [cm].
Definition: TrackSrim.hh:120
bool EstimateRange(const double ekin, const double step, double &stpmax)
Definition: TrackSrim.cc:525
bool m_useLongStraggle
Include longitudinal straggling.
Definition: TrackSrim.hh:90
double m_a
Effective A of the gas.
Definition: TrackSrim.hh:105
void PlotRange()
Make a plot of the projected range as function of the projectile energy.
Definition: TrackSrim.cc:363
int m_maxclusters
Maximum number of clusters allowed (infinite if 0)
Definition: TrackSrim.hh:110
double DedxEM(const double e) const
Definition: TrackSrim.cc:430
bool ReadFile(const std::string &file)
Load data from a SRIM file.
Definition: TrackSrim.cc:79
double Xi(const double x, const double beta2) const
Definition: TrackSrim.cc:438
Abstract base class for track generation.
Definition: Track.hh:14
double GetKineticEnergy() const
Return the kinetic energy of the projectile.
Definition: Track.hh:60
Sensor * m_sensor
Definition: Track.hh:112
bool m_debug
Definition: Track.hh:118
ViewDrift * m_viewer
Definition: Track.hh:116
void PlotCluster(const double x0, const double y0, const double z0)
Definition: Track.cc:192
std::string m_className
Definition: Track.hh:102
double m_energy
Definition: Track.hh:107
void PlotNewTrack(const double x0, const double y0, const double z0)
Definition: Track.cc:187
static std::string FindUnusedCanvasName(const std::string &s)
Find an unused canvas name.
Definition: ViewBase.cc:310
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
Definition: Numerics.cc:1206
double RndmHeedWF(const double w, const double f)
Definition: Random.cc:699
double RndmVavilov(const double rkappa, const double beta2)
Draw a random number from a Vavilov distribution.
Definition: Random.cc:300
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
Definition: Random.hh:107
double RndmGaussian()
Draw a Gaussian random variate with mean zero and standard deviation one.
Definition: Random.hh:24
double RndmLandau()
Draw a random number from a Landau distribution.
Definition: Random.cc:104
int electrons
Number of electrons in this cluster.
Definition: TrackSrim.hh:135
double kinetic
Ion energy when cluster was created.
Definition: TrackSrim.hh:134
double ec
Energy spent to make the cluster.
Definition: TrackSrim.hh:133
double t
Cluster location and time.
Definition: TrackSrim.hh:132