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