20double StepVavilov(
const double rkappa) {
24 }
else if (rkappa < 1) {
26 }
else if (rkappa < 2) {
28 }
else if (rkappa < 3) {
30 }
else if (rkappa < 4) {
32 }
else if (rkappa < 5) {
34 }
else if (rkappa < 6) {
36 }
else if (rkappa < 7) {
38 }
else if (rkappa < 8) {
44double Interpolate(
const double x,
const std::vector<double>& xtab,
45 const std::vector<double>& ytab) {
48 }
else if (x > xtab.back()) {
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";
81 const std::string hdr =
m_className +
"::ReadFile:\n ";
84 fsrim.open(file.c_str(), std::ios::in);
86 std::cerr << hdr <<
"Could not open SRIM file " << file
87 <<
" for reading.\n The file perhaps does not exist.\n";
90 unsigned int nread = 0;
94 std::cout << hdr <<
"SRIM header records from file " << file <<
"\n";
98 while (fsrim.getline(line, 100,
'\n')) {
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) {
111 token = strtok(line,
" []=");
112 token = strtok(NULL,
" []=");
113 token = strtok(NULL,
" []=");
115 m_q = std::atof(token);
117 token = strtok(NULL,
" []=");
118 token = strtok(NULL,
" []=");
119 token = strtok(NULL,
" []=");
121 m_mass = std::atof(token) * AtomicMassUnitElectronVolt;
124 if (!fsrim.getline(line, 100,
'\n')) {
125 std::cerr << hdr <<
"Premature EOF looking for target density (line "
130 if (!fsrim.getline(line, 100,
'\n')) {
131 std::cerr << hdr <<
"Premature EOF looking for target density (line "
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,
" ");
144 while (fsrim.getline(line, 100,
'\n')) {
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) {
150 std::cout << hdr <<
"Stopping units: MeV / (mg/cm2) as expected.\n";
154 std::cerr << hdr <<
"Unknown stopping units. Aborting (line " << nread
160 while (fsrim.getline(line, 100,
'\n')) {
162 if (strstr(line,
"-----------") != NULL)
break;
172 unsigned int ntable = 0;
173 while (fsrim.getline(line, 100,
'\n')) {
175 if (strstr(line,
"-----------") != NULL)
break;
177 token = strtok(line,
" ");
178 m_ekin.push_back(atof(token));
179 token = strtok(NULL,
" ");
180 if (strcmp(token,
"eV") == 0) {
182 }
else if (strcmp(token,
"keV") == 0) {
184 }
else if (strcmp(token,
"GeV") == 0) {
186 }
else if (strcmp(token,
"MeV") != 0) {
187 std::cerr << hdr <<
"Unknown energy unit " << token <<
"; aborting\n";
191 token = strtok(NULL,
" ");
194 token = strtok(NULL,
" ");
197 token = strtok(NULL,
" ");
198 m_range.push_back(atof(token));
199 token = strtok(NULL,
" ");
200 if (strcmp(token,
"A") == 0) {
202 }
else if (strcmp(token,
"um") == 0) {
204 }
else if (strcmp(token,
"mm") == 0) {
206 }
else if (strcmp(token,
"m") == 0) {
208 }
else if (strcmp(token,
"km") == 0) {
210 }
else if (strcmp(token,
"cm") != 0) {
211 std::cerr << hdr <<
"Unknown distance unit " << token <<
"; aborting\n";
215 token = strtok(NULL,
" ");
217 token = strtok(NULL,
" ");
218 if (strcmp(token,
"A") == 0) {
220 }
else if (strcmp(token,
"um") == 0) {
222 }
else if (strcmp(token,
"mm") == 0) {
224 }
else if (strcmp(token,
"m") == 0) {
226 }
else if (strcmp(token,
"km") == 0) {
228 }
else if (strcmp(token,
"cm") != 0) {
229 std::cerr << hdr <<
"Unknown distance unit " << token <<
"; aborting\n";
233 token = strtok(NULL,
" ");
235 token = strtok(NULL,
" ");
236 if (strcmp(token,
"A") == 0) {
238 }
else if (strcmp(token,
"um") == 0) {
240 }
else if (strcmp(token,
"mm") == 0) {
242 }
else if (strcmp(token,
"m") == 0) {
244 }
else if (strcmp(token,
"km") == 0) {
246 }
else if (strcmp(token,
"cm") != 0) {
247 std::cerr << hdr <<
"Unknown distance unit " << token <<
"; aborting\n";
257 while (fsrim.getline(line, 100,
'\n')) {
259 if (strstr(line,
"=============") != NULL) {
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);
268 std::cerr << hdr <<
"Did not find stopping unit scaling; aborting.\n";
272 for (
unsigned int i = 0; i < ntable; ++i) {
279 std::cout << hdr <<
"Successfully read " << file <<
"(" << nread
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] "
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],
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);
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) {
317 yT.push_back(em + hd);
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));
323 TCanvas* celoss =
new TCanvas();
327 celoss->DrawFrame(xmin, 0., xmax, 1.05 * ymax,
";Ion energy [MeV];Energy loss [MeV/cm]");
331 gr.SetLineStyle(kSolid);
333 gr.SetMarkerStyle(21);
334 gr.SetLineColor(kBlue + 1);
335 gr.SetMarkerColor(kBlue + 1);
336 gr.DrawGraph(nPoints,
m_ekin.data(), yE.data(),
"plsame");
338 gr.SetLineColor(kGreen + 2);
339 gr.SetMarkerColor(kGreen + 2);
340 gr.DrawGraph(nPoints,
m_ekin.data(), yH.data(),
"plsame");
342 gr.SetLineColor(kOrange - 3);
343 gr.SetMarkerColor(kOrange - 3);
344 gr.DrawGraph(nPoints,
m_ekin.data(), yT.data(),
"plsame");
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");
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));
368 TCanvas* crange =
new TCanvas();
372 crange->DrawFrame(xmin, 0., xmax, 1.05 * ymax,
";Ion energy [MeV];Projected range [cm]");
374 const unsigned int nPoints =
m_ekin.size();
376 gr.SetLineColor(kOrange - 3);
377 gr.SetMarkerColor(kOrange - 3);
378 gr.SetLineStyle(kSolid);
380 gr.SetMarkerStyle(21);
381 gr.DrawGraph(nPoints,
m_ekin.data(),
m_range.data(),
"plsame");
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),
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]");
401 const unsigned int nPoints =
m_ekin.size();
403 gr.SetLineStyle(kSolid);
405 gr.SetMarkerStyle(21);
407 gr.SetLineColor(kOrange - 3);
408 gr.SetMarkerColor(kOrange - 3);
411 gr.SetLineColor(kGreen + 2);
412 gr.SetMarkerColor(kGreen + 2);
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");
437 constexpr double fconst = 1.e-6 * TwoPi * (
438 FineStructureConstant * FineStructureConstant * HbarC * HbarC) /
439 (ElectronMass * AtomicMassUnit);
444 double& deem,
double& dehd)
const {
447 const std::string hdr =
m_className +
"::PreciseLoss: ";
450 std::cout << hdr <<
"\n"
451 <<
" Initial energy: " << estart <<
" MeV\n"
452 <<
" Step: " << step <<
" cm\n";
455 const double eps = 1.0e-2;
457 unsigned int ndiv = 1;
459 const unsigned int nMaxIter = 10;
460 bool converged =
false;
461 for (
unsigned int iter = 0; iter < nMaxIter; ++iter) {
467 const double s =
m_density * step / ndiv;
468 for (
unsigned int i = 0; i < ndiv; i++) {
472 const double em22 = s *
DedxEM(e2 - 0.5 * de21);
473 const double hd22 = s *
DedxHD(e2 - 0.5 * de21);
477 const double em41 = s *
DedxEM(e4);
478 const double hd41 = s *
DedxHD(e4);
479 const double de41 = em41 + hd41;
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;
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;
489 const double em44 = s *
DedxEM(e4 - de43);
490 const double hd44 = s *
DedxHD(e4 - de43);
491 const double de44 = em44 + hd44;
493 deem += (em41 + em44) / 6. + (em42 + em43) / 3.;
494 dehd += (hd41 + hd44) / 6. + (hd42 + hd43) / 3.;
496 e4 -= (de41 + de44) / 6. + (de42 + de43) / 3.;
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);
505 if (fabs(e2 - e4) > eps * (fabs(e2) + fabs(e4) + fabs(estart))) {
515 std::cerr << hdr <<
"No convergence achieved integrating energy loss.\n";
517 std::cout << hdr <<
"Convergence at eps = " << eps <<
"\n";
530 const std::string hdr =
m_className +
"::EstimateRange: ";
536 double deem = 0., dehd = 0.;
538 double de1 = deem + dehd;
541 if (
m_debug) std::cout << hdr <<
"Initial step OK.\n";
545 double st2 = 0.5 * step;
547 const unsigned int nMaxIter = 20;
548 for (
unsigned int iter = 0; iter < nMaxIter; ++iter) {
553 if (de2 < ekin)
break;
560 std::cerr << hdr <<
"\n Did not find a smaller step in " << nMaxIter
561 <<
" iterations. Abandoned.\n";
562 stpmax = 0.5 * (st1 + st2);
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);
570 for (
unsigned int iter = 0; iter < nMaxIter; ++iter) {
574 std::cerr << hdr <<
"Bisection failed due to equal energy loss for "
575 <<
"two step sizes. Abandoned.\n";
577 stpmax = 0.5 * (st1 + st2);
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);
586 st3 = st1 - (st2 - st1) * (de1 - ekin) / (de2 - de1);
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);
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);
610 std::cout << hdr <<
"Bisection did not converge in " << nMaxIter
611 <<
" steps. Abandoned.\n";
613 stpmax = st1 - (st2 - st1) * (de1 - ekin) / (de2 - de1);
618 const double t0,
const double dx0,
const double dy0,
622 const std::string hdr =
m_className +
"::NewTrack: ";
626 std::cerr << hdr <<
"Sensor is not defined.\n";
631 double xmin = 0., ymin = 0., zmin = 0.;
632 double xmax = 0., ymax = 0., zmax = 0.;
634 std::cerr << hdr <<
"Drift area is not set.\n";
636 }
else if (x0 < xmin || x0 > xmax || y0 < ymin || y0 > ymax || z0 < zmin ||
638 std::cerr << hdr <<
"Initial position outside bounding box.\n";
645 std::cerr << hdr <<
"No medium at initial position.\n";
648 std::cerr << hdr <<
"Medium at initial position is not ionisable.\n";
653 const double normdir = sqrt(dx0 * dx0 + dy0 * dy0 + dz0 * dz0);
657 if (normdir < Small) {
659 std::cout << hdr <<
"Direction vector has zero norm.\n"
660 <<
" Initial direction is randomized.\n";
673 std::cerr << hdr <<
"Particle mass not set.\n";
676 std::cerr << hdr <<
"Particle charge not set.\n";
679 std::cerr << hdr <<
"Initial particle energy not set.\n";
681 }
else if (
m_work < Small) {
682 std::cerr << hdr <<
"Work function not set.\n";
684 }
else if (
m_a < Small ||
m_z < Small) {
685 std::cerr << hdr <<
"A and/or Z not set.\n";
690 if (ekin0 < 1.e-14 *
m_mass || ekin0 < 1.e-3 *
m_work) {
692 std::cout << hdr <<
"Initial kinetic energy E = " << ekin0
693 <<
" MeV such that beta2 = 0 or E << W; particle stopped.\n";
699 const double tracklength = 10 * Interpolate(ekin0,
m_ekin,
m_range);
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);
711 printf(
" Fano factor %g\n",
m_fano);
713 std::cout <<
" Fano factor Not set\n";
717 printf(
" Cluster size %d\n",
m_nsize);
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
763 step = ekin0 / (ncls * (dedxem + dedxhd));
768 double deem = 0., dehd = 0.;
772 if (deem + dehd > e) {
776 deem = e * deem / (dehd + deem);
779 if (
m_debug) std::cout << hdr <<
"Finish raised. Track length reached.\n";
781 stpmax = tracklength - dsum;
784 std::cout << hdr <<
"Maximum step size set to " << stpmax <<
" cm.\n";
789 std::cerr << hdr <<
"Failure computing the minimum step size."
790 <<
"\n Clustering abandoned.\n";
795 if (stpmin > stpmax) {
797 if (
m_debug) std::cout << hdr <<
"stpmin > stpmax. Deposit all energy.\n";
799 if (e - eloss - dehd < 0) eloss = e - dehd;
801 if (
m_debug) std::cout << hdr <<
"Finish raised. Single deposit.\n";
802 }
else if (step < stpmin) {
804 if (
m_debug) std::cout << hdr <<
"Enlarging step size.\n";
807 if (deem + dehd > e) {
808 if (
m_debug) std::cout << hdr <<
"Excess loss. Recomputing stpmax.\n";
812 deem = e * deem / (dehd + deem);
820 if (
m_debug) std::cout << hdr <<
"Using existing step size.\n";
825 if (
m_debug) std::cout << hdr <<
"Truncating negative energy loss.\n";
827 }
else if (eloss > e - dehd) {
828 if (
m_debug) std::cout << hdr <<
"Excess energy loss, using mean.\n";
830 if (e - eloss - dehd < 0) {
833 if (
m_debug) std::cout << hdr <<
"Finish raised. Using mean energy.\n";
837 std::cout << hdr <<
"Step length = " << step <<
" cm.\n "
838 <<
"Mean loss = " << deem <<
" MeV.\n "
839 <<
"Actual loss = " << eloss <<
" MeV.\n";
845 std::cout << hdr <<
"No medium at position (" << x <<
"," << y <<
","
851 std::cout << hdr <<
"Medium at (" << x <<
"," << y <<
"," << z
852 <<
") is not ionisable.\n";
857 std::cout << hdr <<
"Cluster at (" << x <<
"," << y <<
"," << z
858 <<
") outside bounding box.\n";
873 double ecl = 1.e6 * (eloss + epool);
879 if (ernd1 > ecl)
break;
885 std::cout << hdr <<
"EM + pool: " << 1.e6 * (eloss + epool)
887 <<
" eV, E/w: " << (eloss + epool) / (1.0e-6 *
m_work)
888 <<
", n: " << cluster.
electrons <<
".\n";
891 epool += eloss - 1.e-6 * cluster.
ec;
893 std::cout << hdr <<
"Cluster " <<
m_clusters.size() <<
"\n at ("
894 << cluster.
x <<
", " << cluster.
y <<
", " << cluster.
z
895 <<
"),\n e = " << cluster.
ec
897 <<
",\n pool = " << epool <<
" MeV.\n";
906 if (
m_debug) std::cout << hdr <<
"Finishing flag raised.\n";
908 }
else if (e < ekin0 * 1.e-9) {
910 if (
m_debug) std::cout << hdr <<
"Energy exhausted.\n";
914 const double scale = sqrt(step / prange);
919 std::cout << hdr <<
"sigma l, t1, t2: " << sigl <<
", " << sigt1 <<
", "
923 if (xdir * xdir + zdir * zdir <= 0) {
926 }
else if (ydir > 0) {
929 std::cerr << hdr <<
"Zero step length; clustering abandoned.\n";
934 phi = atan2(xdir, zdir);
935 theta = atan2(ydir, sqrt(xdir * xdir + zdir * zdir));
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;
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);
954 std::cerr << hdr <<
"Zero step length; clustering abandoned.\n";
965 std::cerr << hdr <<
"Exceeded maximum number of clusters.\n";
977 const std::string hdr =
m_className +
"::SmallestStep: ";
978 constexpr double expmax = 30;
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
988 }
else if (
m_mass <= 0 || fabs(
m_q) <= 0) {
990 <<
"Track parameters not valid.\n Mass = " << 1.e-6 *
m_mass
991 <<
" MeV, charge = " <<
m_q <<
".\n";
994 std::cerr << hdr <<
"Gas parameters not valid.\n A = " <<
m_a
995 <<
", Z = " <<
m_z <<
" density = " <<
m_density <<
" g/cm3.\n";
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;
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);
1009 double xi =
Xi(step, beta2);
1011 double rkappa = xi / emax;
1014 double stpnow = step;
1015 constexpr unsigned int nMaxIter = 10;
1016 for (
unsigned int iter = 0; iter < nMaxIter; ++iter) {
1020 PrintSettings(hdr, denow, stpnow, ekin, beta2, gamma,
m_a,
m_z,
m_density,
1024 double rknew = rkappa;
1025 if (m_model <= 0 || m_model > 4) {
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);
1035 std::cout << hdr <<
"Landau distribution is imposed.\n kappa_min = "
1036 << rklim <<
", d_min = " << stpmin <<
" cm.\n";
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;
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
1052 if (stpmin > stpnow * 1.1) {
1053 if (
m_debug) std::cout << hdr <<
"Step size increase. New pass.\n";
1058 const double sigma2 = xi * emax * (1 - 0.5 * beta2);
1059 stpmin = stpnow * 16 * sigma2 / (denow * denow);
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";
1067 }
else if (rkappa < 0.05) {
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;
1076 std::cout << hdr <<
"Landau distribution automatic.\n kappa_min = "
1077 << rklim <<
", d_min = " << stpmin <<
" cm.\n";
1079 if (rknew > 0.05 || stpmin > stpnow * 1.1) {
1082 std::cout << hdr <<
"Model change or step increase. New pass.\n";
1085 }
else if (rkappa < 5) {
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;
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";
1098 if (rknew > 5 || stpmin > stpnow * 1.1) {
1101 std::cout << hdr <<
"Model change or step increase. New pass.\n";
1106 const double sigma2 = xi * emax * (1 - 0.5 * beta2);
1107 stpmin = stpnow * 16 * sigma2 / (denow * denow);
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";
1117 if (stpnow > stpmin) {
1119 std::cout << hdr <<
"Step size ok, minimum: " << stpmin <<
" cm\n";
1125 std::cerr << hdr <<
"Step size must be increased to " << stpmin
1133 denow *= stpmin / stpnow;
1135 if (
m_debug) std::cout << hdr <<
"Iteration " << iter <<
"\n";
1136 if (iter == nMaxIter - 1) {
1138 std::cerr << hdr <<
"No convergence reached on step size.\n";
1145 const double step)
const {
1157 const std::string hdr =
"TrackSrim::RndmEnergyLoss: ";
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
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";
1169 std::cerr << hdr <<
"Material parameters not valid.\n A = " <<
m_a
1170 <<
", Z = " <<
m_z <<
", density = " <<
m_density <<
" g/cm3.\n";
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;
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);
1183 const double xi =
Xi(step, beta2);
1185 const double rkappa = xi / emax;
1188 PrintSettings(hdr, de, step, ekin, beta2, gamma,
m_a,
m_z,
m_density,
m_q,
1189 m_mass, emax, xi, rkappa);
1192 if (m_model <= 0 || m_model > 4) {
1194 if (
m_debug) std::cout << hdr <<
"Fixed energy loss.\n";
1197 if (
m_debug) std::cout << hdr <<
"Landau imposed.\n";
1198 const double xlmean = -(log(rkappa) + beta2 + 1. - Gamma);
1202 if (
m_debug) std::cout << hdr <<
"Vavilov imposed.\n";
1203 if (rkappa > 0.01 && rkappa < 12) {
1205 rndde += xi * (xvav + log(rkappa) + beta2 + (1 - Gamma));
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) {
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);
1221 for (
unsigned int iter = 0; iter < 100; ++iter) {
1222 if (xlan < xlmax)
break;
1225 rndde += xi * (xlan - xlmean);
1226 }
else if (rkappa < 5) {
1228 if (
m_debug) std::cout << hdr <<
"Vavilov fast automatic.\n";
1230 rndde += xi * (xvav + log(rkappa) + beta2 + (1 - Gamma));
1233 if (
m_debug) std::cout << hdr <<
"Gaussian automatic.\n";
1234 rndde =
RndmGaussian(de, sqrt(xi * emax * (1 - 0.5 * beta2)));
1238 std::cout << hdr <<
"Energy loss generated = " << rndde <<
" MeV.\n";
1244 double& tcls,
int& n,
double& e,
double& extra) {
1258 n = cluster.electrons;
1260 extra = cluster.kinetic;
Abstract base class for media.
bool IsIonisable() const
Is charge deposition by charged particles/photon enabled in this medium?
bool IsInArea(const double x, const double y, const double z)
Check if a point is inside the user area.
bool GetMedium(const double x, const double y, const double z, Medium *&medium)
Get the medium at (x, y, z).
bool GetArea(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Return the current user area.
std::vector< double > m_longstraggle
Longitudinal straggling [cm].
double RndmEnergyLoss(const double ekin, const double de, const double step) const
unsigned int m_currcluster
Index of the next cluster to be returned.
std::vector< double > m_hdloss
Hadronic energy loss [MeV cm2/g].
void SetDensity(const double density)
Set the density [g/cm3] of the target medium.
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra) override
std::vector< double > m_range
Projected range [cm].
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
double DedxHD(const double e) const
std::vector< Cluster > m_clusters
std::vector< double > m_emloss
EM energy loss [MeV cm2/g].
double m_mass
Mass of ion [MeV].
bool m_chargeset
Charge gas been defined.
bool SmallestStep(double ekin, double de, double step, double &stpmin)
double m_fano
Fano factor [-].
bool PreciseLoss(const double step, const double estart, double &deem, double &dehd) const
double m_work
Work function [eV].
double m_density
Density [g/cm3].
bool m_useTransStraggle
Include transverse straggling.
int m_nsize
Targeted cluster size.
std::vector< double > m_ekin
Energy in energy loss table [MeV].
std::vector< double > m_transstraggle
Transverse straggling [cm].
bool EstimateRange(const double ekin, const double step, double &stpmax)
bool m_useLongStraggle
Include longitudinal straggling.
double m_a
A and Z of the gas.
void PlotRange()
Make a plot of the projected range as function of the projectile energy.
int m_maxclusters
Maximum number of clusters allowed (infinite if 0)
double DedxEM(const double e) const
bool ReadFile(const std::string &file)
Load data from a SRIM file.
double Xi(const double x, const double beta2) const
Abstract base class for track generation.
double GetKineticEnergy() const
Return the kinetic energy of the projectile.
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
double RndmHeedWF(const double w, const double f)
double RndmVavilov(const double rkappa, const double beta2)
Draw a random number from a Vavilov distribution.
void RndmDirection(double &dx, double &dy, double &dz, const double length=1.)
Draw a random (isotropic) direction vector.
double RndmGaussian()
Draw a Gaussian random variate with mean zero and standard deviation one.
double RndmLandau()
Draw a random number from a Landau distribution.
int electrons
Number of electrons in this cluster.
double kinetic
Ion energy when cluster was created.
double ec
Energy spent to make the cluster.
double t
Cluster location and time.