KaliVeda
Toolkit for HIC analysis
KVRangeYanezMaterial.cpp
1 //Created by KVClassFactory on Thu Sep 27 14:50:08 2012
2 //Author: John Frankland,,,
3 
4 #include "KVRangeYanezMaterial.h"
5 #include <fstream>
6 #include "TF1.h"
7 #include "KVUnits.h"
8 #include "KVNameValueList.h"
9 #include <KVNucleus.h>
10 
11 #define NMAX 4000
12 #define NSAV 2000
13 #define NELMAX 10
14 
15 
39 
40 class range {
41  /*
42  Based on the 'range' C library of:
43  Author: Ricardo Yanez
44  Copyright (c) 2005-2011 Ricardo Yanez <ricardo.yanez@calel.org>
45 
46  Calculates energy loss of an ion in an absorber by constructing
47  range tables based on either the,
48 
49  Northcliffe-Schilling correlations (E/A < 12 MeV/A)
50  (L.C. Northcliffe, R.F. Schilling, Nucl. Data Tables A7, 233, 1970).
51 
52  or,
53 
54  Hubert-Bimbot-Gauvin correlations (2.5 < E/A < 100 MeV/A)
55  (Atomic Data and Nuclear Data Tables, 1990, 46, pp. 1-213).
56 
57  License:
58 
59  This program is free software; you can redistribute it and/or modify
60  it under the terms of the GNU General Public License as published by
61  the Free Software Foundation; either version 2 of the License, or
62  (at your option) any later version.
63 
64  This program is distributed in the hope that it will be useful,
65  but WITHOUT ANY WARRANTY; without even the implied warranty of
66  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
67  GNU General Public License for more details.
68 
69  You should have received a copy of the GNU General Public License
70  along with this program; if not, write to the Free Software
71  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
72 
73  */
74 private:
75  const double FMT = 0.005;
76 
77  int isw1 = 0;
78  int isw2 = 0;
79  int isw3 = 0;
80  int isw4 = 0;
81  int isw5 = 0;
82 
83  int numel;
84  struct elem {
85  double w;
86  int z;
87  int a;
88  };
89 
90  std::vector<elem> cmpnd{NELMAX};
91 
92 public:
93  int nelem;
94  int is_gas = 0;
95  std::vector<elem> absorb{NELMAX};
96 
100  double passage(int icorr, int zp, int ap, int iabso, int zt, int at,
101  double ein, double t, double* err);
102 
107  double egassap(int icorr, int zp, int ap, int iabso, int zt, int at,
108  double t, double eut, double* err);
112  double thickn(int icorr, int zp, int ap, int iabso, int zt, int at,
113  double ein, double delen);
114 
118  double rangen(int icorr, int zp, int ap, int iabso, int zt, int at,
119  double ein);
120 
126  void rangetab(int icorr, int zp, int ap, int iabso, int zt, int at,
127  double* em, double* r, int* n);
128 
129  void dedxtab(int icorr, int zp, int ap, int iabso, int zt, int at,
130  double e, double* tdedxe, double* tdedxn);
131 
132  void def_absorber(int zt, int at, int iabso);
133 
136  double SMOOTHERSTEP(double x)
137  {
138  return (pow(x, 3) * (3 * x * (2 * x - 5) + 10));
139  }
140 
144  double dedx(int icorr, double e, int zp, int ap, int zt, int at);
145 
149  double ededx(double e, int zp, int zt);
150 
151  /*
152  Computes current value of nuclear stopping power for a
153  Projectile of energy E (MeV/A), mass ap and charge zp in an
154  absorber with mass at and charge zt using a function fit to
155  the LSS universal nuclear stopping power curve.
156  */
157  double ndedx(double e, int zp, int ap, int zt, int at);
158 
159  /*
160  Compute the "electrical" energy loss rate in any material
161  (-dE/dx) above 2.5 MeV/A according to Hubert et al.
162  */
163  double ededxh(double e, int zp, int zt);
164 
165  /*
166  Given E/A (MeV/A) and the material atomic number Z, A
167  logarithm F is calculated that should be used to convert
168  stopping power in Aluminium to stopping power of the material
169  with atomic number Z by adding it to the log of -(1/Z2)(dE/dx)
170  This routine is for solids only.
171  */
172  void mpyers(double* le, int zt, double* f);
173 
174  /*
175  Given E/A (MeV/A) and the gas material atomic number Z, A
176  logarithm F is calculated that should be used to convert
177  stopping power in Aluminium to stopping power of the gas
178  with atomic number Z by adding it to the log of -(1/Z2)(dE/dx)
179  This routine is for gases only.
180  */
181  void gfact(double* le, int zt, double* f);
182 
183  /*
184  Get corresponding vectors of log(E/A) and log((dE/dx)/Z2) for a
185  specific ion with charge Z in Aluminium
186  */
187  void alion(int zp, double* dedxz2);
188 
189  /*
190  Data pool for aluminium -dE/dx/Z2.
191  Call with energy E/A and get a series of log(Z) and log((-dE/dx)/Z2)
192  */
193  void alref(double le, double* lz, double* dedx);
194 
195  /*
196  Interpolate for current value of s(2,a) function as given by
197  F.Hubert, R.Rimbot and H.Gauvin, Atomic Data and Nuclear Data
198  Tables, 1990, 46, pp. 1-213.
199  */
200  double s2az(double e, int zt);
201 
202  double polint(double* xa, double* ya, int n, double x, double* dy)
203  {
204  double y, c[n], d[n];
205  int i, m, ns = 1;
206  double den, dif, dift, ho, hp, w;
207 
208  dif = fabs(x - *xa);
209 
210  for (i = 0 ; i < n ; i++) {
211  dift = fabs(x - * (xa + i));
212  if (dift < dif) {
213  ns = i + 1;
214  dif = dift;
215  }
216  c[i] = d[i] = *(ya + i);
217  }
218  y = *(ya + ns - 1);
219  ns--;
220  for (m = 0 ; m < n - 1 ; m++) {
221  for (i = 0 ; i < n - m - 1 ; i++) {
222  ho = *(xa + i) - x;
223  hp = *(xa + i + m + 1) - x;
224  w = c[i + 1] - d[i];
225  den = ho - hp;
226  if (den == 0.0)
227  exit(0);
228  den = w / den;
229  d[i] = hp * den;
230  c[i] = ho * den;
231  }
232  if (2 * ns < n - m - 1)
233  *dy = c[ns];
234  else {
235  *dy = d[ns - 1];
236  ns--;
237  }
238  y += *dy;
239  }
240  *dy = fabs(*dy);
241  return y;
242  }
243 
244  unsigned int locate(double y[], int n, double x)
245  {
246 
247  unsigned int jl, jm, ju;
248 
249  jl = 0;
250  ju = n;
251 
252  while ((ju - jl) > 1) {
253  jm = (ju + jl) / 2;
254  if ((y[n - 1] > y[1]) && (x > y[jm]))
255  jl = jm;
256  else
257  ju = jm;
258  }
259  return jl;
260  }
261 
262  void splie2(double TOTO[], double x2a[], double* ya[],
263  int m, int n, double* y2a[])
264  {
265  double ytmp[n], y2tmp[n];
266  int i, j;
267 
268  for (j = 0 ; j < m ; j++) {
269  for (i = 0 ; i < n ; i++) {
270  ytmp[i] = *(ya[i] + j);
271  }
272 
273  spline(x2a, ytmp, n, 1.0e30, 1.0e30, &y2tmp[0]);
274  for (i = 0 ; i < n ; i++) {
275  *(y2a[i] + j) = y2tmp[i];
276  }
277  }
278  }
279 
280 
281  void splin2(double x1a[], double x2a[], double* ya[], double* y2a[],
282  int m, int n, double x1, double x2, double* y)
283  {
284  double ytmp[n], y2tmp[n], yytmp[m];
285  int i, j;
286 
287  for (j = 0 ; j < m ; j++) {
288  for (i = 0 ; i < n ; i++) {
289  ytmp[i] = *(ya[i] + j);
290  y2tmp[i] = *(y2a[i] + j);
291  }
292  splint(x2a, ytmp, y2tmp, n, x2, &yytmp[j]);
293  }
294  spline(x1a, yytmp, m, 1.0e30, 1.0e30, &y2tmp[0]);
295  splint(x1a, yytmp, y2tmp, m, x1, y);
296  }
297 
298  void spline(double x[], double y[], int n, double yp1, double ypn, double* y2)
299  {
300 
301  int k;
302  double p, qn, sig, un, u[200];
303  int m = n - 1;
304 
305  if (yp1 > 0.99e30)
306  *y2 = u[0] = 0.0;
307  else {
308  *y2 = -0.5;
309  u[0] = (3.0 / (x[1] - x[0])) * ((y[1] - y[0]) / (x[1] - x[0]) - yp1);
310  }
311  for (k = 1 ; k < n - 1 ; k++) {
312  sig = (x[k] - x[k - 1]) / (x[k + 1] - x[k - 1]);
313  p = *(y2 + k - 1) * sig + 2.0;
314  *(y2 + k) = (sig - 1.0) / p;
315  u[k] = (6.0 * ((y[k + 1] - y[k]) / (x[k + 1] - x[k]) -
316  (y[k] - y[k - 1]) / (x[k] - x[k - 1])) /
317  (x[k + 1] - x[k - 1]) - sig * u[k - 1]) / p;
318  }
319  if (ypn > 0.99e30)
320  qn = un = 0.0;
321  else {
322  qn = 0.5;
323  un = (3.0 / (x[m] - x[m - 1])) *
324  (ypn - (y[m] - y[m - 1]) / (x[m] - x[m - 1]));
325  }
326  *(y2 + m) = (un - qn * u[m - 1]) / (*(y2 + m - 1) * qn + 1.0);
327  for (k = m - 1 ; k >= 0 ; k--) {
328  *(y2 + k) *= *(y2 + k + 1);
329  *(y2 + k) += u[k];
330  }
331  }
332 
333  void splint(double xa[], double ya[], double y2a[],
334  int n, double x, double* y)
335  {
336 
337  int k, klo, khi;
338  double a, b, h;
339 
340  klo = 1;
341  khi = n;
342 
343  while (khi - klo > 1) {
344  k = (khi + klo) >> 1;
345  if (xa[k - 1] > x)
346  khi = k;
347  else
348  klo = k;
349  }
350  klo--;
351  khi--;
352  h = xa[khi] - xa[klo];
353  if (h == 0.0) {
354  printf("Bad xa input in splint.\n");
355  exit(0);
356  }
357  a = (xa[khi] - x) / h;
358  b = (x - xa[klo]) / h;
359  *y = a * ya[klo] + b * ya[khi] +
360  ((pow(a, 3) - a) * y2a[klo] + (pow(b, 3) - b) * y2a[khi]) * (h * h) / 6.0;
361  }
362 
363 };
364 
365 
366 
368 
369 void range::def_absorber(int zt, int at, int iabso)
370 {
371 
372  int i;
373 
374  switch (iabso) {
375 
376  /* User defined */
377  case -1:
378  for (i = 0 ; i < nelem ; i++) {
379  cmpnd[i].z = absorb[i].z;
380  cmpnd[i].a = absorb[i].a;
381  cmpnd[i].w = absorb[i].w;
382  }
383  numel = nelem;
384  break;
385 
386  /* Single element */
387 
388  case 0:
389  numel = 1;
390  cmpnd[0].z = zt;
391  cmpnd[0].a = at;
392  cmpnd[0].w = at * 1;
393  break;
394 
395  /* Solids */
396 
397  case 1: /* Mylar (C-10, H-8, O-4) */
398  numel = 3;
399  cmpnd[0].z = 6;
400  cmpnd[0].a = 12;
401  cmpnd[0].w = cmpnd[0].a * 10.;
402  cmpnd[1].z = 1;
403  cmpnd[1].a = 1;
404  cmpnd[1].w = cmpnd[1].a * 8.;
405  cmpnd[2].z = 8;
406  cmpnd[2].a = 16;
407  cmpnd[2].w = cmpnd[2].a * 4.;
408  break;
409  case 2: /* Polyethylene (C-2, H-4) */
410  numel = 2;
411  cmpnd[0].z = 6;
412  cmpnd[0].a = 12;
413  cmpnd[0].w = cmpnd[0].a * 2.;
414  cmpnd[1].z = 1;
415  cmpnd[1].a = 1;
416  cmpnd[1].w = cmpnd[1].a * 4.;
417  break;
418  case 3: /* Polypropylene (C-3, H-6) */
419  numel = 2;
420  cmpnd[0].z = 6;
421  cmpnd[0].a = 12;
422  cmpnd[0].w = cmpnd[0].a * 3.;
423  cmpnd[1].z = 1;
424  cmpnd[1].a = 1;
425  cmpnd[1].w = cmpnd[1].a * 6.;
426  break;
427  case 4: /* Kapton (H-10, C-22, N-2, O-5) */
428  numel = 4;
429  cmpnd[0].z = 1;
430  cmpnd[0].a = 1;
431  cmpnd[0].w = cmpnd[0].a * 10.;
432  cmpnd[1].z = 6;
433  cmpnd[1].a = 12;
434  cmpnd[1].w = cmpnd[1].a * 22.;
435  cmpnd[2].z = 7;
436  cmpnd[2].a = 14;
437  cmpnd[2].w = cmpnd[2].a * 2.;
438  cmpnd[3].z = 8;
439  cmpnd[3].a = 16;
440  cmpnd[3].w = cmpnd[3].a * 5.;
441  break;
442  case 5: /* Cesium Iodine (CsI) */
443  numel = 2;
444  cmpnd[0].z = 55;
445  cmpnd[0].a = 133;
446  cmpnd[0].w = cmpnd[0].a * 1.;
447  cmpnd[1].z = 53;
448  cmpnd[1].a = 127;
449  cmpnd[1].w = cmpnd[1].a * 1.;
450  break;
451  case 6: /* Sodium Iodine (NaI) */
452  numel = 2;
453  cmpnd[0].z = 11;
454  cmpnd[0].a = 23;
455  cmpnd[0].w = cmpnd[0].a * 1.;
456  cmpnd[1].z = 53;
457  cmpnd[1].a = 127;
458  cmpnd[1].w = cmpnd[1].a * 1.;
459  break;
460  case 7: /* Aluminum Oxide (Al2O3) */
461  numel = 2;
462  cmpnd[0].z = 13;
463  cmpnd[0].a = 27;
464  cmpnd[0].w = cmpnd[0].a * 2.;
465  cmpnd[1].z = 8;
466  cmpnd[1].a = 16;
467  cmpnd[1].w = cmpnd[1].a * 3.;
468  break;
469  case 8: /* Tin-Lead (Sn60/Pb40) */
470  numel = 4;
471  cmpnd[0].z = 50;
472  cmpnd[0].a = 116;
473  cmpnd[0].w = 0.206 * 60.;
474  cmpnd[1].z = 50;
475  cmpnd[1].a = 118;
476  cmpnd[1].w = 0.340 * 60.;
477  cmpnd[2].z = 50;
478  cmpnd[2].a = 120;
479  cmpnd[2].w = 0.454 * 60.;
480  cmpnd[3].z = 82;
481  cmpnd[3].a = 208;
482  cmpnd[3].w = 40.;
483  break;
484  case 9: /* Natural Ni (Ni-nat) */
485  numel = 5;
486  cmpnd[0].z = 28;
487  cmpnd[0].a = 58;
488  cmpnd[0].w = cmpnd[0].a * 0.680769;
489  cmpnd[1].z = 28;
490  cmpnd[1].a = 60;
491  cmpnd[1].w = cmpnd[1].a * 0.262231;
492  cmpnd[2].z = 28;
493  cmpnd[2].a = 61;
494  cmpnd[2].w = cmpnd[2].a * 0.011399;
495  cmpnd[3].z = 28;
496  cmpnd[3].a = 62;
497  cmpnd[3].w = cmpnd[3].a * 0.036345;
498  cmpnd[4].z = 28;
499  cmpnd[4].a = 64;
500  cmpnd[4].w = cmpnd[4].a * 0.009256;
501  break;
502 
503  /* Gases */
504 
505  case 100: /* Carbon Tetrafluoride (CF4) */
506  numel = 2;
507  cmpnd[0].z = 6;
508  cmpnd[0].a = 12;
509  cmpnd[0].w = cmpnd[0].a * 1.;
510  cmpnd[1].z = 9;
511  cmpnd[1].a = 19;
512  cmpnd[1].w = cmpnd[1].a * 4.;
513  is_gas = 1;
514  break;
515  case 101: /* Propane (H-8, C-3) */
516  numel = 2;
517  cmpnd[0].z = 1;
518  cmpnd[0].a = 1;
519  cmpnd[0].w = cmpnd[0].a * 8.;
520  cmpnd[1].z = 6;
521  cmpnd[1].a = 12;
522  cmpnd[1].w = cmpnd[1].a * 3.;
523  is_gas = 1;
524  break;
525  case 102: /* Butane (H-10, C-4) */
526  numel = 2;
527  cmpnd[0].z = 1;
528  cmpnd[0].a = 1;
529  cmpnd[0].w = cmpnd[0].a * 10.;
530  cmpnd[1].z = 6;
531  cmpnd[1].a = 12;
532  cmpnd[1].w = cmpnd[1].a * 4.;
533  is_gas = 1;
534  break;
535  case 103: /* Octane (H-18, C-8) */
536  numel = 2;
537  cmpnd[0].z = 1;
538  cmpnd[0].a = 1;
539  cmpnd[0].w = cmpnd[0].a * 18.;
540  cmpnd[1].z = 6;
541  cmpnd[1].a = 12;
542  cmpnd[1].w = cmpnd[1].a * 8.;
543  is_gas = 1;
544  break;
545  default:
546  printf("No valid absorber data.\n");
547  exit(0);
548  }
549 }
550 
551 
552 
554 
555 double range::passage(int icorr, int zp, int ap, int iabso, int zt, int at,
556  double ein, double t, double* err)
557 {
558  double em[NMAX], r[NMAX];
559 
560  double elin, elut, rin, rut, lerr;
561  int n;
562 
563  int jj, jjj;
564 
565  /* check correlation */
566  if (icorr == 0 && ein / ap > 12.0)
567  icorr = 1; /* switch to H-B-G */
568  if (icorr == 1 && ein / ap <= 2.5)
569  icorr = 0; /* switch to N-S */
570 
571  rangetab(icorr, zp, ap, iabso, zt, at, &em[0], &r[0], &n);
572  jjj = n - 3;
573 
574  elin = log10(ein / (double)ap);
575  jj = locate(em, n, elin);
576  if (jj > jjj) jj = jjj;
577  rin = polint(&em[jj], &r[jj], 3, elin, &lerr);
578  rut = rin - t;
579  if (rut <= 0.0) {
580  *err = 0.0;
581  return 0.0;
582  }
583  else {
584  jj = locate(r, n, rut);
585  if (jj > jjj) jj = jjj;
586  elut = polint(&r[jj], &em[jj], 3, rut, &lerr);
587  *err = fabs(pow(10.0, elut - lerr * 3) - pow(10.0, elut + lerr * 3)) / pow(10.0, elut);
588  return pow(10.0, elut) * ap;
589  }
590 }
591 
592 
593 
595 
596 double range::egassap(int icorr, int zp, int ap, int iabso, int zt, int at,
597  double t, double eut, double* err)
598 {
599  double em[NMAX], r[NMAX];
600 
601  double elut, elin, eaut, rut, rin, lerr;
602  int n;
603 
604  int jj, jjj;
605 
606  if (eut / (double)ap != 0.0) {
607  if (icorr == 0 && eut / ap > 12.0)
608  icorr = 1; /* switch to H-B-G */
609  }
610 
611  rangetab(icorr, zp, ap, iabso, zt, at, &em[0], &r[0], &n);
612  jjj = n - 3;
613 
614  if (eut / (double)ap != 0.0) {
615  elut = log10(eut / (double)ap);
616  jj = locate(em, n, elut);
617  if (jj > jjj) jj = jjj;
618  rut = polint(&em[jj], &r[jj], 3, elut, &lerr);
619  }
620  else
621  rut = 0.0;
622 
623  rin = rut + t;
624  jj = locate(r, n, rin);
625  if (jj > jjj) jj = jjj;
626  elin = polint(&r[jj], &em[jj], 3, rin, &lerr);
627  *err = fabs(pow(10.0, elin - lerr * 3) - pow(10.0, elin + lerr * 3)) / pow(10.0, elin);
628  eaut = pow(10.0, elin);
629 
630  if (icorr == 0 && eaut > 12.0)
631  printf("warning: Hubert-Bimbot-Gauvin correlations should be used in this case.\n");
632  if (icorr == 1 && eaut <= 2.5)
633  printf("Warning: Northcliffe-Schilling correlations should be used in this case.\n");
634 
635  return eaut * ap;
636 }
637 
638 
639 
641 
642 double range::thickn(int icorr, int zp, int ap, int iabso, int zt, int at,
643  double ein, double delen)
644 {
645  double em[NMAX], r[NMAX];
646 
647  double elin, elut, rin, rut, rerr;
648  int n;
649 
650  int jj, jjj;
651 
652  if (icorr == 0 && ein / (double)ap > 12.0)
653  icorr = 1; /* switch to H-B-G */
654  if (icorr == 1 && ein / (double)ap <= 2.5)
655  icorr = 0; /* switch to N-S */
656 
657  rangetab(icorr, zp, ap, iabso, zt, at, &em[0], &r[0], &n);
658  jjj = n - 3;
659 
660  elin = log10(ein / (double)ap);
661  jj = locate(em, n, elin);
662  if (jj > jjj) jj = jjj;
663  rin = polint(&em[jj], &r[jj], 3, elin, &rerr);
664  if (ein - delen <= 0.0)
665  rut = 0.0;
666  else {
667  elut = log10((ein - delen) / (double)ap);
668  jj = locate(em, n, elut);
669  if (jj > jjj) jj = jjj;
670  rut = polint(&em[jj], &r[jj], 3, elut, &rerr);
671  }
672  return rin - rut;
673 }
674 
675 
676 
678 
679 double range::rangen(int icorr, int zp, int ap, int iabso, int zt, int at,
680  double ein)
681 {
682  double em[NMAX], r[NMAX];
683 
684  double elin, rerr;
685  int n;
686 
687  int jj, jjj;
688 
689  if (icorr == 0 && ein / (double)ap > 12.0)
690  icorr = 1; /* switch to H-B-G */
691  if (icorr == 1 && ein / (double)ap <= 2.5)
692  icorr = 0; /* switch to N-S */
693 
694  rangetab(icorr, zp, ap, iabso, zt, at, &em[0], &r[0], &n);
695  jjj = n - 3;
696 
697  elin = log10(ein / (double)ap);
698  jj = locate(em, n, elin);
699  if (jj > jjj) jj = jjj;
700  return polint(&em[jj], &r[jj], 3, elin, &rerr);
701 
702 }
703 
704 
705 
707 
708 void range::rangetab(int icorr, int zp, int ap, int iabso, int zt, int at,
709  double* em, double* r, int* n)
710 {
711  double elog[62] = {
712  -2.0000000000, -1.9030899870, -1.7958800173, -1.6989700043, -1.6020599913,
713  -1.4948500217, -1.3979400087, -1.3010299957, -1.2218487496, -1.1549019600,
714  -1.0969100130, -1.0457574906, -1.0000000000, -0.9030899870, -0.7958800173,
715  -0.6989700043, -0.6020599913, -0.4948500217, -0.3979400087, -0.3010299957,
716  -0.2218487496, -0.1549019600, -0.0969100130, -0.0457574906, 0.0000000000,
717  0.0969100130, 0.2041199827, 0.3010299957, 0.3979400087, 0.5051499783,
718  0.6020599913, 0.6989700043, 0.7781512504, 0.8450980400, 0.9030899870,
719  0.9542425094, 1.0000000000, 1.0413926852, 1.0791812460, 1.1760912591,
720  1.3010299957, 1.3979400087, 1.4771212547, 1.5440680444, 1.6020599913,
721  1.6532125138, 1.6989700043, 1.7781512504, 1.8450980400, 1.9030899870,
722  1.9542425094, 2.0000000000, 2.0413926852, 2.0791812460, 2.1760912591,
723  2.3010299957, 2.3979400087, 2.4771212547, 2.5440680444, 2.6020599913,
724  2.6532125138, 2.6989700043
725  };
726 
727  int ntalel;
728  double dedxt[NELMAX][NMAX];
729  double wtot;
730 
731  int i, j, k;
732  double est, elg, e;
733  double rng, rval, rnow, rold;
734  double etot, eold;
735  double dedxnow;
736 
737  int jsav;
738  static int isav = 0;
739  static int icc[NSAV], ijj[NSAV], iab[NSAV];
740  static int izp[NSAV], iap[NSAV], izt[NSAV], iat[NSAV];
741  static double esav[NSAV][NMAX], rsav[NSAV][NMAX];
742 
743 
744  /* check if table saved */
745  if (isav > 0) {
746  for (i = 0 ; i < isav ; i++) {
747  jsav = 0;
748  if (icorr == icc[i]) jsav++;
749  if (iabso == iab[i]) jsav++;
750  if (zp == izp[i]) jsav++;
751  if (ap == iap[i]) jsav++;
752  if (zt == izt[i]) jsav++;
753  if (at == iat[i]) jsav++;
754  if (jsav == 6) {
755  k = ijj[i];
756  for (j = 0 ; j < k ; j++) {
757  *(em + j) = esav[i][j];
758  *(r + j) = rsav[i][j];
759  }
760  *n = k;
761  return;
762  }
763  }
764  }
765 
766  if (isav > NSAV) {
767  printf("WARNING in rangetab: NSAV too small\n");
768  exit(0);
769  }
770  /* Save this call */
771  icc[isav] = icorr;
772  iab[isav] = iabso;
773  izp[isav] = zp;
774  iap[isav] = ap;
775  izt[isav] = zt;
776  iat[isav] = at;
777 
778  switch (icorr) {
779  case 0:
780  ntalel = 38;
781  break;
782  case 1:
783  ntalel = 61;
784  break;
785  case 2:
786  ntalel = 61;
787  break;
788  default:
789  printf("No valid range correlation.\n");
790  exit(0);
791  }
792 
793  /* define absorber */
794  def_absorber(zt, at, iabso);
795 
796  /* Compute a range table */
797  est = 3 * elog[0];
798  wtot = 0.0;
799  for (i = 0 ; i < numel ; i++) {
800  isw1 = isw2 = 0;
801  zt = cmpnd[i].z;
802  at = cmpnd[i].a;
803  wtot += cmpnd[i].w;
804  *n = 0;
805  for (j = 1 ; j <= 7000 ; j++) {
806  elg = FMT * (j - 1200);
807  e = pow(exp(elg), log(10.0));
808  if (elg >= est) {
809  if (elg <= elog[ntalel]) {
810  dedxt[i][*n] = dedx(icorr, e, zp, ap, zt, at);
811  *(em + (*n)) = elg;
812  (*n)++;
813  }
814  else
815  break;
816  if (*n > 3999) break;
817  }
818  }
819  }
820 
821  ijj[isav] = *n;
822 
823  rng = 0.0;
824  rold = 0.0;
825  eold = 0.0;
826  for (j = 0 ; j < *n ; j++) {
827  elg = *(em + j);
828  e = pow(exp(elg), log(10.0));
829  etot = e * ap;
830  dedxnow = 0.0;
831  for (i = 0 ; i < numel ; i++)
832  dedxnow += dedxt[i][j] * cmpnd[i].w;
833  dedxnow /= wtot;
834  rnow = 1.0 / dedxnow;
835  rval = 0.5 * (rold + rnow) * (etot - eold);
836  rng += rval;
837  *(r + j) = rng;
838  eold = etot;
839  rold = rnow;
840  esav[isav][j] = *(em + j);
841  rsav[isav][j] = *(r + j);
842  }
843  isav++;
844 }
845 
846 
847 
849 
850 void range::dedxtab(int icorr, int zp, int ap, int iabso, int zt, int at,
851  double e, double* tdedxe, double* tdedxn)
852 {
853  double dedxn[NELMAX], dedxe[NELMAX];
854  double tw;
855  double e_int;
856  int i;
857 
858  def_absorber(zt, at, iabso);
859 
860  for (i = 0 ; i < numel ; i++) {
861  isw1 = isw2 = 0;
862  zt = cmpnd[i].z;
863  at = cmpnd[i].a;
864  if (e < 2.5) {
865  dedxn[i] = ndedx(e, zp, ap, zt, at) * pow(zp, 2);
866  dedxe[i] = ededx(e, zp, zt) * pow(zp, 2);
867  }
868  else {
869  if (icorr == 1) {
870  dedxn[i] = 0.0;
871  dedxe[i] = ededxh(e, zp, zt);
872  }
873  else if (icorr == 0) {
874  dedxn[i] = ndedx(e, zp, ap, zt, at) * pow(zp, 2);
875  dedxe[i] = ededx(e, zp, zt) * pow(zp, 2);
876  }
877  else if (icorr == 2) {
878  if (e < 12) {
879  /* interpolate */
880  e_int = (e - 2.5) / 9.5;
881  dedxn[i] = (1. - e_int) * ndedx(e, zp, ap, zt, at) * pow(zp, 2);
882  dedxe[i] = e_int * ededxh(e, zp, zt) + (1. - e_int) * ededx(e, zp, zt) * pow(zp, 2);
883  }
884  else {
885  dedxn[i] = 0.0;
886  dedxe[i] = ededxh(e, zp, zt);
887  }
888  }
889  }
890  }
891  *tdedxn = *tdedxe = 0.0;
892  tw = 0.0;
893  for (i = 0 ; i < numel ; i++) {
894  *tdedxn += dedxn[i] * cmpnd[i].w;
895  *tdedxe += dedxe[i] * cmpnd[i].w;
896  tw += cmpnd[i].w;
897  }
898  *tdedxn /= tw;
899  *tdedxe /= tw;
900 }
901 
902 
903 
905 
906 double range::dedx(int icorr, double e, int zp, int ap, int zt, int at)
907 {
908  double dedxn, dedxe, e_int, smoothe;
909  static double emin = 0.1, emax = 2.5;
910 
911  if (icorr == 1 && e < 2.5) icorr = 0;
912  if (icorr == 2) {
913  if (e < emin) icorr = 0;
914  else if (e > emax) icorr = 1;
915  }
916  switch (icorr) {
917  case 0:
918  dedxn = ndedx(e, zp, ap, zt, at);
919  dedxe = ededx(e, zp, zt);
920  return (dedxn + dedxe) * pow(zp, 2);
921  break;
922  case 1:
923  return ededxh(e, zp, zt);
924  break;
925  case 2: /* interpolate in range emin<e<emax */
926  e_int = (e - emin) / (emax - emin);
927  smoothe = SMOOTHERSTEP(e_int);
928  return (smoothe * ededxh(e, zp, zt) +
929  (1. - smoothe) * (ndedx(e, zp, ap, zt, at) + ededx(e, zp, zt)) * pow(zp, 2));
930  break;
931  default:
932  exit(0);
933  }
934 }
935 
936 
937 
939 
940 double range::ededx(double e, int zp, int zt)
941 {
942 
943  double elog[42] = {-1.903090, -1.795880, -1.698970, -1.602060, -1.494850,
944  -1.397940, -1.301030, -1.221849, -1.154902, -1.096910,
945  -1.045757, -1.000000, -0.903090, -0.795880, -0.698970,
946  -0.602060, -0.494850, -0.397940, -0.301030, -0.221849,
947  -0.154902, -0.096910, -0.045757, 0.000000, +0.096910,
948  +0.204120, +0.301030, +0.397940, +0.505150, +0.602060,
949  +0.698970, +0.778151, +0.845098, +0.903090, +0.954243,
950  +1.000000, +1.041393, +1.079181, +1.301030, +1.602060,
951  +1.845098, +2.000000
952  };
953 
954  int zgases[11] = {1, 2, 7, 8, 9, 10, 17, 18, 36, 54, 86};
955 
956  int i, j;
957  unsigned int jj;
958  static double dedxz2[42];
959  static double ak, a, ftargl;
960  double b, ftarg, el, err;
961  int gas;
962 
963  if (!isw1) {
964  alion(zp, &dedxz2[0]);
965  ak = (dedxz2[2] - dedxz2[0]) / (elog[2] - elog[0]);
966  a = dedxz2[0] - ak * elog[0];
967  //ftarg = 1.0;
968 
969  /* Is it a gas? */
970  gas = 0;
971  for (i = 0 ; i < 11 ; i++)
972  if (zt == zgases[i]) gas = 1;
973 
974  /* Special case for gases */
975  if (gas) {
976  gfact(&elog[0], zt, &ftargl);
977  dedxz2[0] += ftargl;
978  for (j = 1 ; j < 42 ; j++) {
979  gfact(&elog[j], zt, &ftarg);
980  dedxz2[j] += ftarg;
981  }
982  }
983  /* It is a solid */
984  else {
985  if (zt != 13) {
986  mpyers(&elog[0], zt, &ftargl);
987  dedxz2[0] += ftargl;
988  for (j = 1 ; j < 42 ; j++) {
989  mpyers(&elog[j], zt, &ftarg);
990  dedxz2[j] += ftarg;
991  }
992  }
993  }
994  isw1 = 1;
995  }
996  el = log10(e);
997  if (el < elog[0])
998  b = a + ak * el + ftargl;
999  else {
1000  jj = locate(elog, 42, el);
1001  if (jj > 39) jj = 39;
1002  b = polint(&elog[jj], &dedxz2[jj], 3, el, &err);
1003  }
1004  return pow(10.0, b);
1005 }
1006 
1007 
1008 
1010 
1011 double range::ndedx(double e, int zp, int ap, int zt, int at)
1012 {
1013 
1014  double zexpo, dedr;
1015  double x, y, xl, yl;
1016 
1017  /* Compute current x-coordinate from projectile and target data */
1018  zexpo = sqrt(pow(zp, 2.0 / 3.0) + pow(zt, 2.0 / 3.0));
1019  x = sqrt((e * ap * at / (double)(ap + at)) / (zp * zt * zexpo));
1020 
1021  /* Get log for polynomial fit */
1022  xl = log10(x);
1023  yl = -6.80959 + xl *
1024  (-6.60315 + xl * (-1.73474 + xl * xl * (0.04937 + xl * 0.00486)));
1025  y = pow(10.0, yl);
1026 
1027  /* Convert y to -dEPS/dRHO */
1028  dedr = y * zt * ap / (zp * at * (ap + at) * zexpo);
1029  return dedr;
1030 }
1031 
1032 
1033 
1035 
1036 double range::ededxh(double e, int zp, int zt)
1037 {
1038  static double b, c, d;
1039  static double x1, x2, x3, x4;
1040  double xg1;
1041 
1042  /* Special case for He */
1043  if (zp == 2) {
1044  isw2 = 1;
1045  return s2az(e, zt) * 4.0;
1046  }
1047 
1048  if (!isw2) {
1049  if (is_gas) {
1050  /* If the element or compound is a gas
1051  we use the (unpublished) effective charge
1052  parametrisation for gaseous absorbers
1053  (M.F. Rivet, IPN Orsay, private communication).
1054  */
1055  b = 2.6878;
1056  c = 0.07452;
1057  d = 0.51216 + 0.5693 * exp(-0.002158 * zt);
1058  x2 = 11.109 + 0.1057 * log(zt);
1059  x3 = 0.5622 - 0.07901 * log(zt);
1060  x4 = 0.8276 - 0.03059 * log(zt);
1061  }
1062  else {
1063  if (zt == 4) {
1064  b = 2.0;
1065  c = 0.04369;
1066  d = 2.045;
1067  x2 = 7.0;
1068  x3 = 0.2643;
1069  x4 = 0.4171;
1070  }
1071  else {
1072  if (zt == 6) {
1073  b = 1.91;
1074  c = 0.03958;
1075  d = 2.584;
1076  x2 = 6.933;
1077  x3 = 0.2433;
1078  x4 = 0.3969;
1079  }
1080  else {
1081  b = 1.658;
1082  c = 0.0517;
1083  d = 1.164 + 0.2319 * exp(-0.004302 * zt);
1084  x2 = 8.144 + 0.09876 * log(zt);
1085  x3 = 0.314 + 0.01072 * log(zt);
1086  x4 = 0.5218 + 0.02521 * log(zt);
1087  }
1088  }
1089  }
1090  x1 = d + b * exp(-c * zp);
1091  isw2 = 1;
1092  }
1093  xg1 = 1.0 - x1 * exp(-x2 * pow(e, x3) / pow(zp, x4));
1094  return s2az(e, zt) * pow((xg1 * zp), 2);
1095 }
1096 
1097 
1098 
1100 
1101 void range::mpyers(double* le, int zt, double* f)
1102 {
1103  double za[12] = {4., 6., 13., 22., 28., 32., 40., 47., 63., 73., 79., 92.};
1104 
1105  double ea[38] = {0.0125, .016, .02, .025, .032, .04, .05, .06, .07, .08, .09, .1,
1106  .125, .16, .2, .25, .32, .4, .5, .6, .7, .8, .9, 1., 1.25, 1.6, 2.,
1107  2.5, 3.2, 4., 5., 6., 7., 8., 9.0, 10., 11., 12.0
1108  };
1109 
1110  double fb[38][12] = {
1111  {1.6499, 1.3651, 1., .665, .5395, .4901, .454, .42, .2669, .2284, .21, .1799},
1112  {1.6501, 1.3652, 1., .6651, .54, .4901, .4541, .4201, .2671, .2285, .2099, .1801},
1113  {1.65, 1.365, 1., .6649, .54, .49, .454, .4199, .267, .2285, .2101, .1801},
1114  {1.65, 1.365, 1., .6649, .5401, .49, .4541, .42, .2669, .2286, .2105, .1801},
1115  {1.6482, 1.3662, 1., .6651, .541, .492, .454, .42, .267, .229, .2111, .1811},
1116  {1.6441, 1.3691, 1., .667, .544, .496, .457, .4231, .269, .2316, .2129, .1831},
1117  {1.638, 1.3729, 1., .67, .549, .4999, .461, .427, .271, .2354, .217, .1865},
1118  {1.632, 1.381, 1., .674, .554, .503, .465, .431, .276, .24, .2215, .1905},
1119  {1.625, 1.3901, 1., .678, .56, .508, .47, .436, .28, .245, .227, .195},
1120  {1.6181, 1.4021, 1., .683, .565, .5121, .473, .44, .285, .25, .231, .199},
1121  {1.611, 1.414, 1., .688, .57, .518, .477, .444, .291, .2545, .236, .2035},
1122  {1.6049, 1.427, 1., .693, .576, .522, .48, .449, .296, .259, .24, .207},
1123  {1.59, 1.467, 1., .704, .588, .533, .49, .459, .3075, .269, .2515, .2175},
1124  {1.569, 1.5251, 1., .716, .602, .546, .501, .47, .32, .281, .2635, .23},
1125  {1.548, 1.571, 1., .727, .615, .557, .511, .479, .332, .293, .275, .24},
1126  {1.523, 1.6021, 1., .739, .628, .57, .522, .489, .345, .305, .286, .25},
1127  {1.492, 1.621, 1., .753, .643, .587, .535, .5, .358, .318, .2985, .2625},
1128  {1.462, 1.623, 1., .762, .656, .6, .547, .511, .372, .331, .311, .273},
1129  {1.43, 1.607, 1., .773, .67, .615, .56, .522, .386, .344, .3235, .285},
1130  {1.402, 1.579, 1., .783, .682, .628, .57, .531, .397, .355, .3345, .297},
1131  {1.379, 1.548, 1., .79, .691, .637, .577, .538, .406, .363, .343, .304},
1132  {1.358, 1.514, 1., .797, .699, .647, .585, .545, .414, .372, .351, .312},
1133  {1.342, 1.483, 1., .802, .707, .656, .592, .551, .423, .38, .359, .32},
1134  {1.327, 1.453, 1., .809, .714, .663, .599, .557, .43, .387, .366, .327},
1135  {1.297, 1.389, 1., .82, .73, .68, .613, .569, .447, .403, .381, .341},
1136  {1.266, 1.327, 1., .835, .746, .697, .628, .58, .465, .421, .398, .357},
1137  {1.239, 1.293, 1., .847, .761, .714, .64, .593, .48, .438, .414, .373},
1138  {1.213, 1.271, 1., .859, .777, .73, .653, .606, .498, .453, .43, .39},
1139  {1.189, 1.256, 1., .87, .792, .745, .668, .618, .516, .472, .448, .407},
1140  {1.17, 1.246, 1., .879, .805, .758, .683, .63, .533, .488, .465, .423},
1141  {1.154, 1.237, 1., .886, .816, .77, .694, .642, .55, .502, .481, .439},
1142  {1.144, 1.23, 1., .892, .823, .78, .704, .651, .558, .515, .493, .453},
1143  {1.134, 1.225, 1., .896, .828, .787, .712, .659, .57, .524, .503, .463},
1144  {1.127, 1.22, 1., .9, .833, .792, .72, .665, .578, .533, .512, .472},
1145  {1.123, 1.215, 1., .903, .837, .797, .725, .671, .585, .54, .519, .48},
1146  {1.118, 1.21, 1., .907, .84, .801, .729, .677, .59, .547, .527, .487},
1147  {1.114, 1.207, 1., .91, .844, .805, .733, .683, .596, .554, .535, .493},
1148  {1.111, 1.203, 1., .912, .847, .809, .737, .688, .601, .56, .54, .498}
1149  };
1150 
1151  int i, j;
1152  static double lfb[38][12], y2b[38][12];
1153  static double* pfb[38], *py2b[38];
1154  static double lza[12], lea[38];
1155  double zl;
1156 
1157  if (!isw4) {
1158  for (i = 0 ; i < 38 ; i++) {
1159  lea[i] = log10(ea[i]);
1160  for (j = 0 ; j < 12 ; j++)
1161  lfb[i][j] = log10(fb[i][j]);
1162  }
1163  for (j = 0 ; j < 12 ; j++)
1164  lza[j] = log10(za[j]);
1165  for (i = 0 ; i < 38 ; i++) {
1166  pfb[i] = &(lfb[i])[0];
1167  py2b[i] = &(y2b[i])[0];
1168  }
1169  splie2(lza, lea, &pfb[0], 12, 38, &py2b[0]);
1170  isw4 = 1;
1171  }
1172 
1173  zl = log10(zt);
1174  if (*le < lea[0])
1175  *le = lea[0];
1176  splin2(lza, lea, &pfb[0], &py2b[0], 12, 38, zl, *le, f);
1177 }
1178 
1179 
1180 
1182 
1183 void range::gfact(double* le, int zt, double* f)
1184 {
1185  double za[9] = {1., 2., 7., 8., 10., 18., 36., 54., 86.};
1186 
1187  double ea[38] = {0.0125, .016, .02, .025, .032, .04, .05, .06, .07, .08, .09, .1,
1188  .125, .16, .2, .25, .32, .4, .5, .6, .7, .8, .9, 1., 1.25, 1.6, 2.,
1189  2.5, 3.2, 4., 5., 6., 7., 8., 9.0, 10., 11., 12.0
1190  };
1191 
1192  double far[38] = {0.584, .5891, .595, .607, .6251, .6461, .671, .6970, .722, .747,
1193  .771, .795, .851, .923, .985, 1.034, 1.0550, 1.051, 1.031, 1.001,
1194  .964, .928, .894, .871, .843, .8350, .842, .862, .885, .9, .909,
1195  .914, .918, .919, .9200, .92, .92, .9200
1196  };
1197 
1198  double fa[38][9] = {
1199  {6.4213, 2.7396, 1.8149, 1.7055, 1.5153, 1., .6132, .4349, .2868},
1200  {6.2480, 2.6351, 1.7709, 1.6654, 1.4804, 1., .6195, .4415, .2931},
1201  {6.0166, 2.4705, 1.7091, 1.6102, 1.4403, 1., .6235, .4536, .3036},
1202  {5.7334, 2.2751, 1.6345, 1.5404, 1.3939, 1., .6393, .4661, .3378},
1203  {5.4557, 2.0702, 1.5344, 1.4608, 1.3343, 1., .6495, .4895, .3408},
1204  {5.1854, 1.8853, 1.4441, 1.3900, 1.2800, 1., .6625, .5078, .3607},
1205  {4.9329, 1.7348, 1.3695, 1.3323, 1.2400, 1., .6736, .5246, .3800},
1206  {4.7917, 1.6397, 1.3299, 1.2854, 1.2152, 1., .6800, .5351, .3931},
1207  {4.6953, 1.5900, 1.3020, 1.2645, 1.1994, 1., .6815, .5415, .4003},
1208  {4.6454, 1.5596, 1.2866, 1.2504, 1.1901, 1., .6855, .5462, .4056},
1209  {4.6042, 1.5446, 1.2775, 1.2451, 1.1854, 1., .6873, .5499, .4072},
1210  {4.6038, 1.5447, 1.2805, 1.2403, 1.1799, 1., .6906, .5497, .4076},
1211  {4.6652, 1.5571, 1.2868, 1.2397, 1.1845, 1., .6886, .5488, .4078},
1212  {4.7995, 1.6002, 1.3131, 1.2654, 1.2004, 1., .6858, .5460, .4052},
1213  {4.9848, 1.6548, 1.3553, 1.3015, 1.2305, 1., .6802, .5401, .3990},
1214  {5.2418, 1.7351, 1.4072, 1.3549, 1.2650, 1., .6760, .5300, .3897},
1215  {5.6304, 1.8379, 1.4948, 1.4265, 1.3052, 1., .6673, .5204, .3801},
1216  {5.9373, 1.9601, 1.5566, 1.4796, 1.3397, 1., .6613, .5148, .3749},
1217  {6.1881, 2.0563, 1.5975, 1.5199, 1.3598, 1., .6595, .5150, .3754},
1218  {6.3436, 2.1279, 1.6104, 1.5305, 1.3676, 1., .6633, .5175, .3786},
1219  {6.3693, 2.1888, 1.6100, 1.5301, 1.3693, 1., .6691, .5207, .3848},
1220  {6.3254, 2.2198, 1.6045, 1.5301, 1.3674, 1., .6767, .5280, .3922},
1221  {6.2639, 2.2326, 1.5995, 1.5246, 1.3557, 1., .6823, .5391, .4005},
1222  {6.1883, 2.2377, 1.5970, 1.5155, 1.3502, 1., .6889, .5488, .4076},
1223  {5.8955, 2.2076, 1.5718, 1.5006, 1.3333, 1., .7022, .5670, .4282},
1224  {5.4133, 2.1198, 1.5306, 1.4635, 1.3102, 1., .7174, .5832, .4515},
1225  {4.9170, 2.0155, 1.4941, 1.4228, 1.2874, 1., .7304, .5998, .4715},
1226  {4.3851, 1.8851, 1.4397, 1.3701, 1.2645, 1., .7424, .6183, .4919},
1227  {3.9549, 1.7571, 1.3933, 1.3299, 1.2350, 1., .7605, .6328, .5073},
1228  {3.6667, 1.6578, 1.3555, 1.3000, 1.2222, 1., .7667, .6456, .5222},
1229  {3.4983, 1.5896, 1.3300, 1.2805, 1.2068, 1., .7778, .6567, .5390},
1230  {3.4027, 1.5471, 1.3151, 1.2648, 1.1980, 1., .7812, .6641, .5471},
1231  {3.3442, 1.5185, 1.3050, 1.2582, 1.1895, 1., .7865, .6710, .5545},
1232  {3.3080, 1.4962, 1.3003, 1.2546, 1.1850, 1., .7867, .6746, .5604},
1233  {3.2717, 1.4870, 1.2946, 1.2500, 1.1826, 1., .7881, .6772, .5652},
1234  {3.2500, 1.4837, 1.2946, 1.2500, 1.1826, 1., .7902, .6793, .5685},
1235  {3.2282, 1.4804, 1.2945, 1.2500, 1.1826, 1., .7924, .6815, .5707},
1236  {3.2066, 1.4772, 1.2946, 1.2500, 1.1826, 1., .7946, .6837, .5728}
1237  };
1238 
1239  int i, j;
1240  unsigned int jj;
1241  static double lfa[38][9], y2a[38][9];
1242  static double* pfa[38], *py2a[38];
1243  static double lza[9], lea[38], lfar[38];
1244  double zl, fgl, fgal, err;
1245 
1246  if (!isw3) {
1247  for (i = 0 ; i < 38 ; i++) {
1248  lea[i] = log10(ea[i]);
1249  lfar[i] = log10(far[i]);
1250  for (j = 0 ; j < 9 ; j++)
1251  lfa[i][j] = log10(fa[i][j]);
1252  }
1253  for (j = 0 ; j < 9 ; j++)
1254  lza[j] = log10(za[j]);
1255  for (i = 0 ; i < 38 ; i++) {
1256  pfa[i] = &(lfa[i])[0];
1257  py2a[i] = &(y2a[i])[0];
1258  }
1259  splie2(lza, lea, &pfa[0], 9, 38, &py2a[0]);
1260  isw3 = 1;
1261  }
1262  zl = log10(zt);
1263  if (*le < lea[0])
1264  *le = lea[0];
1265  splin2(lza, lea, &pfa[0], &py2a[0], 9, 38, zl, *le, &fgl);
1266  jj = locate(lea, 38, *le);
1267  if (jj > 35) jj = 35;
1268  fgal = polint(&lea[jj], &lfar[jj], 3, *le, &err);
1269  *f = fgl + fgal;
1270 }
1271 
1272 
1273 
1275 
1276 void range::alion(int zp, double* dedxz2)
1277 {
1278  double elog[42] = {-1.903090, -1.795880, -1.698970, -1.602060, -1.494850,
1279  -1.397940, -1.301030, -1.221849, -1.154902, -1.096910,
1280  -1.045757, -1.000000, -0.903090, -0.795880, -0.698970,
1281  -0.602060, -0.494850, -0.397940, -0.301030, -0.221849,
1282  -0.154902, -0.096910, -0.045757, 0.000000, +0.096910,
1283  +0.204120, +0.301030, +0.397940, +0.505150, +0.602060,
1284  +0.698970, +0.778151, +0.845098, +0.903090, +0.954243,
1285  +1.000000, +1.041393, +1.079181, +1.301030, +1.602060,
1286  +1.845098, +2.000000
1287  };
1288 
1289  int j;
1290  unsigned int jj;
1291  double dedx[22], lz[22];
1292  double zlog, err;
1293 
1294  zlog = log10(zp);
1295 
1296  for (j = 0; j < 42 ; j++) {
1297  alref(elog[j], &lz[0], &dedx[0]);
1298  jj = locate(lz, 22, zlog);
1299  if (jj > 19) jj = 19;
1300  *(dedxz2 + j) = polint(&lz[jj], &dedx[jj], 3, zlog, &err);
1301  }
1302 }
1303 
1304 
1305 
1307 
1308 void range::alref(double le, double* lz, double* dedx)
1309 {
1310  double zval[22] = {1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 16.,
1311  20., 25., 32., 40., 50., 61., 79., 100.
1312  };
1313 
1314  double elog[42] = {-1.903090, -1.795880, -1.698970, -1.602060, -1.494850,
1315  -1.397940, -1.301030, -1.221849, -1.154902, -1.096910,
1316  -1.045757, -1.000000, -0.903090, -0.795880, -0.698970,
1317  -0.602060, -0.494850, -0.397940, -0.301030, -0.221849,
1318  -0.154902, -0.096910, -0.045757, 0.000000, +0.096910,
1319  +0.204120, +0.301030, +0.397940, +0.505150, +0.602060,
1320  +0.698970, +0.778151, +0.845098, +0.903090, +0.954243,
1321  +1.000000, +1.041393, +1.079181, +1.301030, +1.602060,
1322  +1.845098, +2.000000
1323  };
1324 
1325  /* H; Z=1 */
1326  double h[42] = {-0.67572, -0.62160, -0.57349, -0.52143, -0.46980, -0.43063,
1327  -0.39794, -0.37779, -0.36653, -0.35952, -0.35655, -0.35655,
1328  -0.36351, -0.38195, -0.40782, -0.44129, -0.48545, -0.53018,
1329  -0.58004, -0.62709, -0.66555, -0.70115, -0.73049, -0.75945,
1330  -0.82102, -0.88941, -0.95468, -1.02228, -1.10237, -1.17393,
1331  -1.24413, -1.30103, -1.35655, -1.39794, -1.43180, -1.46852,
1332  -1.49485, -1.52288, -1.72584, -1.95861, -2.14874, -2.25181
1333  };
1334 
1335  /* He; Z=2 */
1336  double he[42] = {-0.87615, -0.82246, -0.77404, -0.72584, -0.67162, -0.62663,
1337  -0.58503, -0.55479, -0.53276, -0.51606, -0.50376, -0.49485,
1338  -0.48247, -0.48050, -0.48845, -0.50585, -0.53387, -0.56623,
1339  -0.60380, -0.64589, -0.68194, -0.71388, -0.74232, -0.76828,
1340  -0.82536, -0.89279, -0.95664, -1.02342, -1.09963, -1.17070,
1341  -1.24222, -1.30103, -1.35164, -1.39523, -1.43474, -1.46852,
1342  -1.49826, -1.53018, -1.70997, -1.94310, -2.13077, -2.23657
1343  };
1344 
1345  /* Li; Z=3 */
1346  double li[42] = {-1.03990, -0.98623, -0.93763, -0.88904, -0.83532, -0.78870,
1347  -0.74473, -0.71145, -0.68566, -0.66577, -0.65018, -0.63827,
1348  -0.61939, -0.60985, -0.61103, -0.61979, -0.63618, -0.65561,
1349  -0.67778, -0.70358, -0.72816, -0.75175, -0.77440, -0.79588,
1350  -0.84568, -0.90658, -0.96658, -1.02996, -1.10360, -1.17249,
1351  -1.24328, -1.30200, -1.35218, -1.39674, -1.43573, -1.47137,
1352  -1.50246, -1.53264, -1.70333, -1.93554, -2.11919, -2.22915
1353  };
1354 
1355  /* Be; Z=4 */
1356  double be[42] = {-1.16630, -1.11280, -1.06424, -1.01604, -0.96232, -0.91498,
1357  -0.86967, -0.83472, -0.80740, -0.78615, -0.76955, -0.75696,
1358  -0.73785, -0.72830, -0.72801, -0.73283, -0.74172, -0.75187,
1359  -0.76321, -0.77875, -0.79385, -0.80879, -0.82373, -0.83863,
1360  -0.87533, -0.92445, -0.97675, -1.03533, -1.10582, -1.17352,
1361  -1.24365, -1.30266, -1.35347, -1.39794, -1.43771, -1.47334,
1362  -1.50515, -1.53480, -1.70115, -1.92812, -2.11919, -2.22915
1363  };
1364 
1365  /* B; Z=5 */
1366  double b[42] = {-1.26857, -1.21496, -1.16673, -1.11827, -1.06449, -1.01827,
1367  -0.97469, -0.94157, -0.91578, -0.89551, -0.87943, -0.86646,
1368  -0.84430, -0.82786, -0.81976, -0.81667, -0.81793, -0.82206,
1369  -0.82822, -0.83791, -0.84783, -0.85811, -0.86864, -0.87943,
1370  -0.90714, -0.94692, -0.99140, -1.04383, -1.10969, -1.17496,
1371  -1.24413, -1.30277, -1.35379, -1.39881, -1.43890, -1.47470,
1372  -1.50808, -1.53820, -1.69897, -1.92812, -2.11351, -2.22915
1373  };
1374 
1375  /* C; Z=6 */
1376  double c[42] = {-1.36597, -1.31227, -1.26382, -1.21546, -1.16168, -1.11335,
1377  -1.06803, -1.03230, -1.00327, -0.97939, -0.95960, -0.94310,
1378  -0.91315, -0.88941, -0.87642, -0.86967, -0.86708, -0.86753,
1379  -0.86967, -0.87751, -0.88587, -0.89477, -0.90406, -0.91364,
1380  -0.93825, -0.97356, -1.01323, -1.06048, -1.12094, -1.18192,
1381  -1.24795, -1.30491, -1.35491, -1.39915, -1.43903, -1.47496,
1382  -1.50786, -1.53843, -1.69897, -1.92812, -2.10791, -2.23657
1383  };
1384 
1385  /* N; Z=7 */
1386  double n[42] = {-1.43468, -1.38099, -1.33264, -1.28417, -1.23065, -1.18207,
1387  -1.13789, -1.10316, -1.07498, -1.05171, -1.03209, -1.01543,
1388  -0.98331, -0.95348, -0.93242, -0.91721, -0.90694, -0.90295,
1389  -0.90309, -0.90873, -0.91546, -0.92289, -0.93091, -0.93930,
1390  -0.96160, -0.99419, -1.03152, -1.07625, -1.13389, -1.19230,
1391  -1.25579, -1.31053, -1.35877, -1.40150, -1.44002, -1.47509,
1392  -1.50693, -1.53638, -1.69897, -1.92812, -2.11351, -2.23657
1393  };
1394 
1395  /* O; Z=8 */
1396  double o[42] = {-1.51481, -1.46120, -1.41260, -1.36417, -1.31064, -1.26211,
1397  -1.21679, -1.18066, -1.15095, -1.12594, -1.10461, -1.08619,
1398  -1.04977, -1.01456, -0.98855, -0.96859, -0.95321, -0.94441,
1399  -0.93930, -0.94119, -0.94554, -0.95157, -0.95873, -0.96658,
1400  -0.98809, -1.01971, -1.05552, -1.09793, -1.15220, -1.20728,
1401  -1.26761, -1.32032, -1.36701, -1.40876, -1.44672, -1.48149,
1402  -1.51348, -1.54302, -1.70115, -1.92996, -2.11919, -2.23657
1403  };
1404 
1405  /* F; Z=9 */
1406  double f[42] = {-1.59335, -1.53983, -1.49135, -1.44295, -1.38931, -1.34087,
1407  -1.29243, -1.25489, -1.22382, -1.19744, -1.17473, -1.15490,
1408  -1.11504, -1.07534, -1.04469, -1.02002, -0.99984, -0.98727,
1409  -0.97881, -0.97881, -0.98183, -0.98680, -0.99298, -1.00000,
1410  -1.01946, -1.04827, -1.08092, -1.11968, -1.16939, -1.22033,
1411  -1.27653, -1.32619, -1.37067, -1.41086, -1.44759, -1.48149,
1412  -1.51281, -1.54206, -1.70333, -1.93554, -2.11919, -2.23657
1413  };
1414 
1415  /* Ne; Z=10 */
1416  double ne[42] = {-1.667562, -1.614036, -1.565431, -1.516984, -1.463442,
1417  -1.414991, -1.366532, -1.327440, -1.294821, -1.267044,
1418  -1.242984, -1.221849, -1.178814, -1.134837, -1.099851,
1419  -1.070581, -1.045661, -1.029374, -1.017729, -1.016554,
1420  -1.018544, -1.022505, -1.027751, -1.033858, -1.051294,
1421  -1.077638, -1.107905, -1.144118, -1.190912, -1.239050,
1422  -1.292430, -1.339704, -1.382161, -1.420559, -1.455684,
1423  -1.488117, -1.518128, -1.546223, -1.709965, -1.939302,
1424  -2.124939, -2.244125
1425  };
1426 
1427  /* Na; Z=11 */
1428  double na[42] = {-1.721058, -1.667311, -1.618892, -1.570501, -1.516820,
1429  -1.468416, -1.419933, -1.380786, -1.348226, -1.320407,
1430  -1.296389, -1.275318, -1.232446, -1.188746, -1.153929,
1431  -1.124556, -1.098528, -1.080058, -1.065126, -1.059368,
1432  -1.057930, -1.059286, -1.062507, -1.067048, -1.081744,
1433  -1.106069, -1.135107, -1.170563, -1.216675, -1.264098,
1434  -1.316521, -1.362626, -1.403721, -1.440816, -1.474580,
1435  -1.505523, -1.534028, -1.560602, -1.712198, -1.939302,
1436  -2.130768, -2.244125
1437  };
1438 
1439  /* Mg; Z=13 */
1440  double mg[42] = {-1.772578, -1.719030, -1.670517, -1.622183, -1.568525,
1441  -1.520073, -1.471637, -1.432043, -1.398770, -1.370336,
1442  -1.345650, -1.324005, -1.280013, -1.235689, -1.201234,
1443  -1.172622, -1.146496, -1.126147, -1.107635, -1.098269,
1444  -1.093830, -1.092708, -1.093905, -1.096722, -1.108137,
1445  -1.129466, -1.156326, -1.189926, -1.234445, -1.280588,
1446  -1.331705, -1.376751, -1.416896, -1.453012, -1.485895,
1447  -1.515997, -1.543782, -1.569643, -1.716699, -1.943095,
1448  -2.130768, -2.244125
1449  };
1450 
1451  /* Al; Z=13 */
1452  double al[42] = {-1.819477, -1.765788, -1.717342, -1.668938, -1.615315,
1453  -1.566832, -1.518362, -1.478769, -1.445056, -1.416178,
1454  -1.391120, -1.368989, -1.324092, -1.279083, -1.244712,
1455  -1.216274, -1.189346, -1.167302, -1.146395, -1.134325,
1456  -1.127585, -1.124494, -1.123946, -1.125247, -1.133765,
1457  -1.152303, -1.177082, -1.208978, -1.251858, -1.296734,
1458  -1.346616, -1.390614, -1.429789, -1.465058, -1.497104,
1459  -1.526405, -1.553485, -1.578552, -1.721246, -1.946922,
1460  -2.130768, -2.251812
1461  };
1462 
1463  /* S; Z=16 */
1464  double s[42] = {-1.940188, -1.886579, -1.838047, -1.789669, -1.736050,
1465  -1.687585, -1.639158, -1.599556, -1.566068, -1.536375,
1466  -1.510448, -1.487543, -1.440552, -1.393635, -1.358479,
1467  -1.328625, -1.298392, -1.272186, -1.246453, -1.229580,
1468  -1.218301, -1.211042, -1.206734, -1.204636, -1.206133,
1469  -1.217544, -1.236718, -1.263853, -1.302185, -1.343333,
1470  -1.389623, -1.430608, -1.467176, -1.500023, -1.529833,
1471  -1.557043, -1.582165, -1.605329, -1.739929, -1.954677,
1472  -2.142668, -2.251812
1473  };
1474 
1475  /* Ca; Z=20 */
1476  double ca[42] = {-2.075979, -2.021819, -1.972854, -1.923997, -1.869827,
1477  -1.820951, -1.771985, -1.732066, -1.698265, -1.668978,
1478  -1.643162, -1.619608, -1.570934, -1.521578, -1.483795,
1479  -1.450506, -1.415782, -1.385129, -1.354774, -1.334630,
1480  -1.320095, -1.309649, -1.302291, -1.297333, -1.292494,
1481  -1.296881, -1.310092, -1.331777, -1.364667, -1.401237,
1482  -1.443065, -1.480402, -1.513818, -1.543900, -1.571217,
1483  -1.596151, -1.619111, -1.640402, -1.761954, -1.968592,
1484  -2.148742, -2.259637
1485  };
1486 
1487  /* Mn; Z=25 */
1488  double mn[42] = {-2.225571, -2.170053, -2.119827, -2.069642, -2.014125,
1489  -1.963946, -1.913754, -1.872740, -1.838033, -1.807990,
1490  -1.781528, -1.757816, -1.707638, -1.655498, -1.613580,
1491  -1.575563, -1.535642, -1.500335, -1.465385, -1.442445,
1492  -1.425311, -1.412424, -1.402779, -1.395653, -1.385879,
1493  -1.384429, -1.392159, -1.408383, -1.435381, -1.466787,
1494  -1.503646, -1.537003, -1.567095, -1.594346, -1.619152,
1495  -1.641882, -1.662852, -1.682271, -1.787812, -1.982967,
1496  -2.161151, -2.267606
1497  };
1498 
1499  /* Ge; Z=32 */
1500  double ge[42] = {-2.400492, -2.343408, -2.291715, -2.240111, -2.182995,
1501  -2.131319, -2.079657, -2.037496, -2.001785, -1.970886,
1502  -1.943639, -1.919231, -1.867598, -1.810463, -1.763088,
1503  -1.718798, -1.672302, -1.631539, -1.591419, -1.564711,
1504  -1.544545, -1.529158, -1.517344, -1.508296, -1.494333,
1505  -1.487817, -1.490406, -1.500990, -1.521326, -1.546631,
1506  -1.577491, -1.606185, -1.632484, -1.656595, -1.678751,
1507  -1.699203, -1.718177, -1.735865, -1.823909, -2.002177,
1508  -2.187087, -2.283997
1509  };
1510 
1511  /* Zr; Z=40 */
1512  double zr[42] = {-2.561853, -2.503243, -2.450231, -2.397194, -2.338601,
1513  -2.285565, -2.232566, -2.189264, -2.152620, -2.120904,
1514  -2.092925, -2.067907, -2.014882, -1.956245, -1.903242,
1515  -1.852826, -1.799560, -1.753563, -1.709214, -1.677858,
1516  -1.654381, -1.636459, -1.622625, -1.611899, -1.594536,
1517  -1.584078, -1.582591, -1.588464, -1.602995, -1.622728,
1518  -1.648011, -1.672309, -1.695106, -1.716345, -1.736142,
1519  -1.754611, -1.771888, -1.788112, -1.861697, -2.026872,
1520  -2.200659, -2.301030
1521  };
1522 
1523  /* Sn; Z=50 */
1524  double sn[42] = {-2.721064, -2.660906, -2.606460, -2.552036, -2.491821,
1525  -2.437422, -2.383000, -2.338528, -2.300926, -2.268347,
1526  -2.239608, -2.213930, -2.159492, -2.099283, -2.044871,
1527  -1.990430, -1.930228, -1.878243, -1.830878, -1.794081,
1528  -1.766699, -1.745819, -1.729629, -1.716934, -1.695613,
1529  -1.680761, -1.675002, -1.676195, -1.685072, -1.699405,
1530  -1.719194, -1.739042, -1.758185, -1.776379, -1.793595,
1531  -1.809859, -1.825243, -1.839808, -1.913640, -2.060481,
1532  -2.229148, -2.318759
1533  };
1534 
1535  /* Pm; Z=61 */
1536  double pm[42] = {-2.859781, -2.798118, -2.742451, -2.686714, -2.625043,
1537  -2.569315, -2.513602, -2.468054, -2.429555, -2.396193,
1538  -2.366784, -2.340466, -2.284742, -2.223095, -2.167368,
1539  -2.111644, -2.049993, -1.993513, -1.942557, -1.901296,
1540  -1.870337, -1.846490, -1.827754, -1.812816, -1.786842,
1541  -1.766889, -1.756499, -1.752995, -1.756612, -1.766145,
1542  -1.781107, -1.797027, -1.812861, -1.828225, -1.842956,
1543  -1.857026, -1.870441, -1.883229, -1.958607, -2.096910,
1544  -2.259637, -2.337242
1545  };
1546 
1547  /* Au; Z=79 */
1548  double au[42] = {-3.042131, -2.978549, -2.921062, -2.863644, -2.800058,
1549  -2.742599, -2.685136, -2.638190, -2.598525, -2.564142,
1550  -2.533801, -2.506670, -2.449214, -2.385659, -2.328194,
1551  -2.270741, -2.207173, -2.149714, -2.092264, -2.047594,
1552  -2.012868, -1.985183, -1.962713, -1.944210, -1.910215,
1553  -1.881059, -1.862441, -1.851053, -1.846528, -1.849244,
1554  -1.857867, -1.868886, -1.880736, -1.892718, -1.904509,
1555  -1.915963, -1.927015, -1.937638, -2.017729, -2.154902,
1556  -2.288193, -2.361511
1557  };
1558 
1559  /* Fm; Z=100 */
1560  double fm[42] = {-3.208520, -3.143211, -3.084231, -3.025212, -2.959912,
1561  -2.900872, -2.841879, -2.793633, -2.752862, -2.717559,
1562  -2.686407, -2.658546, -2.599514, -2.534231, -2.475215,
1563  -2.416189, -2.350899, -2.291885, -2.232866, -2.188767,
1564  -2.152835, -2.122876, -2.097497, -2.075726, -2.033033,
1565  -1.992295, -1.962577, -1.940611, -1.925985, -1.921442,
1566  -1.924464, -1.931844, -1.941035, -1.950879, -1.960828,
1567  -1.970612, -1.980107, -1.989259, -2.070581, -2.193820,
1568  -2.309804, -2.371611
1569  };
1570 
1571  double err;
1572  int j;
1573  unsigned int jj;
1574 
1575  for (j = 0; j < 22 ; j++)
1576  *(lz + j) = log10(zval[j]);
1577 
1578  jj = locate(elog, 42, le);
1579  if (jj > 39) jj = 39;
1580 
1581  /* Interpolate data for Al to given energy for all standard ions */
1582 
1583  /* First hydrogen ions */
1584  *(dedx++) = polint(&elog[jj], &h[jj], 3, le, &err);
1585 
1586  /* Then Helium ions */
1587  *(dedx++) = polint(&elog[jj], &he[jj], 3, le, &err);
1588 
1589  /* Then Lithium ions */
1590  *(dedx++) = polint(&elog[jj], &li[jj], 3, le, &err);
1591 
1592  /* Then Beryllium ions */
1593  *(dedx++) = polint(&elog[jj], &be[jj], 3, le, &err);
1594 
1595  /* Then Bor ions */
1596  *(dedx++) = polint(&elog[jj], &b[jj], 3, le, &err);
1597 
1598  /* Then Carbon ions */
1599  *(dedx++) = polint(&elog[jj], &c[jj], 3, le, &err);
1600 
1601  /* Then Nitrogen ions */
1602  *(dedx++) = polint(&elog[jj], &n[jj], 3, le, &err);
1603 
1604  /* Then Oxygen ions */
1605  *(dedx++) = polint(&elog[jj], &o[jj], 3, le, &err);
1606 
1607  /* Then Fluorine ions */
1608  *(dedx++) = polint(&elog[jj], &f[jj], 3, le, &err);
1609 
1610  /* Then Neon ions */
1611  *(dedx++) = polint(&elog[jj], &ne[jj], 3, le, &err);
1612 
1613  /* Then Sodium ions */
1614  *(dedx++) = polint(&elog[jj], &na[jj], 3, le, &err);
1615 
1616  /* Then Magnesium ions */
1617  *(dedx++) = polint(&elog[jj], &mg[jj], 3, le, &err);
1618 
1619  /* Then Aluminium ions */
1620  *(dedx++) = polint(&elog[jj], &al[jj], 3, le, &err);
1621 
1622  /* Then Sulfur ions */
1623  *(dedx++) = polint(&elog[jj], &s[jj], 3, le, &err);
1624 
1625  /* Then Calcium ions */
1626  *(dedx++) = polint(&elog[jj], &ca[jj], 3, le, &err);
1627 
1628  /* Then Manganese ions */
1629  *(dedx++) = polint(&elog[jj], &mn[jj], 3, le, &err);
1630 
1631  /* Then Germanium ions */
1632  *(dedx++) = polint(&elog[jj], &ge[jj], 3, le, &err);
1633 
1634  /* Then Zirconium ions */
1635  *(dedx++) = polint(&elog[jj], &zr[jj], 3, le, &err);
1636 
1637  /* Then Tin ions */
1638  *(dedx++) = polint(&elog[jj], &sn[jj], 3, le, &err);
1639 
1640  /* Then Prometium ions */
1641  *(dedx++) = polint(&elog[jj], &pm[jj], 3, le, &err);
1642 
1643  /* Then Gold ions */
1644  *(dedx++) = polint(&elog[jj], &au[jj], 3, le, &err);
1645 
1646  /* Then Fermium ions */
1647  *(dedx++) = polint(&elog[jj], &fm[jj], 3, le, &err);
1648 }
1649 
1650 
1651 
1653 
1654 double range::s2az(double e, int zt)
1655 {
1656  double za[18] = {4.0, 6.0, 13.0, 14.0, 22.0, 26.0, 28.0, 29.0, 32.0, 34.0, 40.0,
1657  47.0, 50.0, 64.0, 73.0, 79.0, 82.0, 92.0
1658  };
1659 
1660  double el[38] = {0.916291, 1.098612, 1.252763, 1.386294, 1.504077, 1.609438,
1661  1.704748, 1.791759, 1.871802, 1.945910, 2.079442, 2.197225,
1662  2.302585, 2.397895, 2.484907, 2.708050, 2.995732, 3.218876,
1663  3.401197, 3.555348, 3.688879, 3.806662, 3.912023, 4.007333,
1664  4.094345, 4.174387, 4.248495, 4.382027, 4.499810, 4.605170,
1665  5.010635, 5.298317, 5.521461, 5.703782, 5.857933, 5.991465,
1666  6.109248, 6.214608
1667  };
1668 
1669  static double sa2[18][38] = {
1670  {
1671  2.19060, 2.32662, 2.44357, 2.54625, 2.63806, 2.72038, 2.79565, 2.86470, 2.92901,
1672  2.98826, 3.09555, 3.19114, 3.27677, 3.35455, 3.42575, 3.60822, 3.84436, 4.02575,
1673  4.17501, 4.29953, 4.40693, 4.50104, 4.58488, 4.66015, 4.72819, 4.79060, 4.84820,
1674  4.95048, 5.04019, 5.11892, 5.41429, 5.61371, 5.75009, 5.85432, 5.93698, 6.00455,
1675  6.06082, 6.10800
1676  },
1677  {
1678  2.12549, 2.25929, 2.37462, 2.47575, 2.56623, 2.64754, 2.72190, 2.79035, 2.85337,
1679  2.91231, 3.01849, 3.11283, 3.19785, 3.27479, 3.34529, 3.52676, 3.76038, 3.94119,
1680  4.08936, 4.21313, 4.31999, 4.41372, 4.49699, 4.57174, 4.63976, 4.70168, 4.75890,
1681  4.86103, 4.94978, 5.02829, 5.32210, 5.51089, 5.64575, 5.75009, 5.83275, 5.89980,
1682  5.95610, 6.00253
1683  },
1684  {
1685  2.35942, 2.48561, 2.59461, 2.69082, 2.77700, 2.85467, 2.92574, 2.99124, 3.05177,
1686  3.10834, 3.21078, 3.30158, 3.38360, 3.45777, 3.52591, 3.70095, 3.92841, 4.10440,
1687  4.24925, 4.37010, 4.47392, 4.56547, 4.64677, 4.72002, 4.78639, 4.84724, 4.90290,
1688  5.00267, 5.09008, 5.16685, 5.45439, 5.64292, 5.77635, 5.87814, 5.95900, 6.02606,
1689  6.08139, 6.12728
1690  },
1691  {
1692  2.33666, 2.46187, 2.57047, 2.66607, 2.75161, 2.82895, 2.89951, 2.96472, 3.02516,
1693  3.08129, 3.18327, 3.27346, 3.35527, 3.42883, 3.49661, 3.67202, 3.89837, 4.07454,
1694  4.21821, 4.33897, 4.44241, 4.53378, 4.61522, 4.68828, 4.75454, 4.81527, 4.87109,
1695  4.97081, 5.05812, 5.13492, 5.42275, 5.61783, 5.75324, 5.85519, 5.93603, 6.00152,
1696  6.05654, 6.10240
1697  },
1698  {
1699  2.52604, 2.64649, 2.75161, 2.84387, 2.92667, 3.00175, 3.07046, 3.13385, 3.19236,
1700  3.24740, 3.34671, 3.43501, 3.51409, 3.58723, 3.65351, 3.82470, 4.04698, 4.21991,
1701  4.36027, 4.47854, 4.58194, 4.67184, 4.75222, 4.82426, 4.88952, 4.94942, 5.00490,
1702  5.10316, 5.18946, 5.26537, 5.54999, 5.73837, 5.87013, 5.97166, 6.05122, 6.11703,
1703  6.17179, 6.21711
1704  },
1705  {
1706  2.58296, 2.70008, 2.80222, 2.89273, 2.97348, 3.04703, 3.11452, 3.17666, 3.23399,
1707  3.28809, 3.38582, 3.47216, 3.55086, 3.62216, 3.68788, 3.85730, 4.07601, 4.24750,
1708  4.38805, 4.50533, 4.60667, 4.69592, 4.77537, 4.84692, 4.91204, 4.97154, 5.02600,
1709  5.12394, 5.20939, 5.28491, 5.56816, 5.75324, 5.88351, 5.98449, 6.06404, 6.12958,
1710  6.18384, 6.22972
1711  },
1712  {
1713  2.59931, 2.71394, 2.81383, 2.90270, 2.98232, 3.05442, 3.12073, 3.18206, 3.23844,
1714  3.29145, 3.38803, 3.47377, 3.55173, 3.62216, 3.68688, 3.85493, 4.07307, 4.24227,
1715  4.38203, 4.49856, 4.60018, 4.68910, 4.76828, 4.83931, 4.90425, 4.96328, 5.01804,
1716  5.11558, 5.20028, 5.27557, 5.55709, 5.74071, 5.87102, 5.97166, 6.05122, 6.11590,
1717  6.17059, 6.21586
1718  },
1719  {
1720  2.65962, 2.77379, 2.87307, 2.96133, 3.04073, 3.11227, 3.17845, 3.23972, 3.29616,
1721  3.34884, 3.44515, 3.53102, 3.60822, 3.67794, 3.74334, 3.91077, 4.12738, 4.29769,
1722  4.43543, 4.55400, 4.65384, 4.74271, 4.82177, 4.89285, 4.95721, 5.01653, 5.07078,
1723  5.16817, 5.25334, 5.32826, 5.60961, 5.79425, 5.92474, 6.02399, 6.10351, 6.16820,
1724  6.22214, 6.26722
1725  },
1726  {
1727  2.71018, 2.82304, 2.92155, 3.00932, 3.08785, 3.15943, 3.22515, 3.28542, 3.34175,
1728  3.39397, 3.49003, 3.57466, 3.65158, 3.72140, 3.78649, 3.95285, 4.16853, 4.33897,
1729  4.47634, 4.59275, 4.69373, 4.78221, 4.86103, 4.93194, 4.99636, 5.05537, 5.10977,
1730  5.20665, 5.29134, 5.36660, 5.64717, 5.80914, 5.93982, 6.03961, 6.12044, 6.18505,
1731  6.23993, 6.28584
1732  },
1733  {
1734  2.75436,
1735  2.86558, 2.96278, 3.04913, 3.12641, 3.19724, 3.26231, 3.32216, 3.37773, 3.42960,
1736  3.52421, 3.60914, 3.68489, 3.75502, 3.81899, 3.98459, 4.19971, 4.36812, 4.50533,
1737  4.62283, 4.72283, 4.81097, 4.88986, 4.96042, 5.02486, 5.08401, 5.13790, 5.23534,
1738  5.32005, 5.39483, 5.67520, 5.84305, 5.97264, 6.07268, 6.15281, 6.21837, 6.27118,
1739  6.31720
1740  },
1741  {
1742  2.77780, 2.88822, 2.98479, 3.07046, 3.14714, 3.21763, 3.28208, 3.34104, 3.39621,
1743  3.44750, 3.54132, 3.62497, 3.70095, 3.77009, 3.83275, 3.99676, 4.20976, 4.37605,
1744  4.51442, 4.62793, 4.72819, 4.81558, 4.89352, 4.96363, 5.02753, 5.08603, 5.14003,
1745  5.23628, 5.32056, 5.39483, 5.67374, 5.85956, 5.98847, 6.08798, 6.16582, 6.23099,
1746  6.28449, 6.32974
1747  },
1748  {
1749  2.86734, 2.97446, 3.06884, 3.15239, 3.22766, 3.29616, 3.35886, 3.41733, 3.47135,
1750  3.52167, 3.61377, 3.69590, 3.77009, 3.83854, 3.90084, 4.06139, 4.27228, 4.43543,
1751  4.57077, 4.68313, 4.78161, 4.86783, 4.94485, 5.01389, 5.07678, 5.13450, 5.18767,
1752  5.28244, 5.36553, 5.43873, 5.71383, 5.89797, 6.02502, 6.12385, 6.20219, 6.26590,
1753  6.31858, 6.36398
1754  },
1755  {
1756  2.92667, 3.03240, 3.12527, 3.20769, 3.28208, 3.34955, 3.41201, 3.46974, 3.52337,
1757  3.57288, 3.66419, 3.74651, 3.82013, 3.88733, 3.94895, 4.10895, 4.31811, 4.48141,
1758  4.61522, 4.72847, 4.82644, 4.91238, 4.98900, 5.05812, 5.12101, 5.17876, 5.23159,
1759  5.32620, 5.40925, 5.48284, 5.75718, 5.94077, 6.06835, 6.16701, 6.24507, 6.30892,
1760  6.36253, 6.40698
1761  },
1762  {
1763  3.03812, 3.14018, 3.22956, 3.30907, 3.38140, 3.44672, 3.50739, 3.56313, 3.61563,
1764  3.66419, 3.75289, 3.83275, 3.90455, 3.97124, 4.03137, 4.18827, 4.39369, 4.55448,
1765  4.68638, 4.79815, 4.89485, 4.97986, 5.05537, 5.12394, 5.18588, 5.24288, 5.29532,
1766  5.38934, 5.47089, 5.54358, 5.81583, 6.00051, 6.12958, 6.22719, 6.30481, 6.36834,
1767  6.42071, 6.46467
1768  },
1769  {
1770  3.15063, 3.24612, 3.32981, 3.40521, 3.47377, 3.53702, 3.59448, 3.64870, 3.69893,
1771  3.74545, 3.83160, 3.90828, 3.97790, 4.04270, 4.10137, 4.25416, 4.45524, 4.61295,
1772  4.74271, 4.85235, 4.94766, 5.03135, 5.10605, 5.17345, 5.23487, 5.29134, 5.34279,
1773  5.43586, 5.51710, 5.58867, 5.85781, 6.04066, 6.17299, 6.26986, 6.34671, 6.40850,
1774  6.46147, 6.50563
1775  },
1776  {
1777  3.21016, 3.30294, 3.38508, 3.45856, 3.52591, 3.58723, 3.64391, 3.69590, 3.74545,
1778  3.79091, 3.87521, 3.95154, 4.02017, 4.08341, 4.14144, 4.29182, 4.49006, 4.64573,
1779  4.77389, 4.88257, 4.97660, 5.05969, 5.13365, 5.20028, 5.26150, 5.31699, 5.36820,
1780  5.46025, 5.54103, 5.61234, 5.88082, 6.06189, 6.19358, 6.28987, 6.36543, 6.42842,
1781  6.47923, 6.52420
1782  },
1783  {
1784  3.21763, 3.31113, 3.39323, 3.46734, 3.53530, 3.59630, 3.65351, 3.70603, 3.75609,
1785  3.80205, 3.88733, 3.96332, 4.03137, 4.09535, 4.15250, 4.30451, 4.50329, 4.65963,
1786  4.78819, 4.89720, 4.99157, 5.07477, 5.14904, 5.21582, 5.27656, 5.33239, 5.38388,
1787  5.47625, 5.55644, 5.62821, 5.89615, 6.07811, 6.20962, 6.30481, 6.38155, 6.44402,
1788  6.49565, 6.53965
1789  },
1790  {
1791  3.26035, 3.35455, 3.43812, 3.51325, 3.58092, 3.64391, 3.70095, 3.75395, 3.80429,
1792  3.85022, 3.93478, 4.01184, 4.08044, 4.14459, 4.20304, 4.35363, 4.55234, 4.70859,
1793  4.83679, 4.94555, 5.03942, 5.12227, 5.19666, 5.26343, 5.32415, 5.37953, 5.43071,
1794  5.52271, 5.60349, 5.67447, 5.94172, 6.12044, 6.24507, 6.34102, 6.41611, 6.47923,
1795  6.53103, 6.57307
1796  }
1797  };
1798 
1799  int i;
1800  static double y2a[18][38];
1801  static double* psa2[18], *py2a[18];
1802  double sa2ln;
1803  double le;
1804 
1805  if (!isw5) {
1806  for (i = 0 ; i < 18 ; i++) {
1807  psa2[i] = &(sa2[i])[0];
1808  py2a[i] = &(y2a[i])[0];
1809  }
1810  splie2(el, za, &psa2[0], 38, 18, &py2a[0]);
1811  isw5 = 1;
1812  }
1813  le = log(e);
1814  splin2(el, za, &psa2[0], &py2a[0], 38, 18, le, zt, &sa2ln);
1815  return exp(-sa2ln);
1816 }
1817 
1818 
1820 
1821 
1822 
1826  : range_lib(new range)
1827 {
1828  // Default constructor
1829 }
1830 
1831 
1832 
1835 
1837  const KVIonRangeTable* t, const Char_t* name, const Char_t* symbol, const Char_t* state, Double_t density,
1838  Int_t Z, Int_t A)
1839  : KVIonRangeTableMaterial(t, name, symbol, state, density, Z, A),
1840  range_lib(new range)
1841 {
1842  // Create material (single-element absorber)
1843  fTableType = 1; // Hubert for E/A>2.5MeV, switches to Northcliffe for <2.5AMeV.
1844  fNelem = 1;
1845  iabso = 0;
1846  fAbsorb[0].z = Z;
1847  fAbsorb[0].a = A;
1848  fAbsorb[0].w = A;
1849 
1851 }
1852 
1853 
1854 
1855 
1862 
1864  : KVIonRangeTableMaterial(), range_lib(new range)
1865 {
1866  // Copy constructor
1867  // This ctor is used to make a copy of an existing object (for example
1868  // when a method returns an object), and it is always a good idea to
1869  // implement it.
1870  // If your class allocates memory in its constructor(s) then it is ESSENTIAL :-)
1871 
1872  obj.Copy(*this);
1873 }
1874 
1875 
1876 
1877 
1885 
1887 {
1888  // This method copies the current state of 'this' object into 'obj'
1889  // You should add here any member variables, for example:
1890  // (supposing a member variable KVRangeYanezMaterial::fToto)
1891  // CastedObj.fToto = fToto;
1892  // or
1893  // CastedObj.SetToto( GetToto() );
1894 
1896  //KVRangeYanezMaterial& CastedObj = (KVRangeYanezMaterial&)obj;
1897 }
1898 
1899 
1900 
1902 
1904 {
1905  fDeltaE = new TF1(Form("KVRangeYanezMaterial:%s:EnergyLoss", GetType()), this,
1907  0., 1.e+03, 0, "KVRangeYanezMaterial", "DeltaEFunc");
1908  fDeltaE->SetNpx(100);
1909  fRange = new TF1(Form("KVRangeYanezMaterial:%s:Range", GetType()), this,
1911  0., 1.e+03, 0, "KVRangeYanezMaterial", "RangeFunc");
1912  fRange->SetNpx(100);
1913  fEres = new TF1(Form("KVRangeYanezMaterial:%s:ResidualEnergy", GetType()), this,
1915  0., 1.e+03, 0, "KVRangeYanezMaterial", "EResFunc");
1916  fEres->SetNpx(100);
1917 }
1918 
1919 
1920 KVRangeYanezMaterial::~KVRangeYanezMaterial() = default; // this for PIMPL to work : https://www.cppstories.com/2018/01/pimpl/
1921 
1922 
1930 
1932 {
1933  // Function parameterising the energy loss of charged particles in this material.
1934  // This is simply an interface to the passage() function of the Range C library.
1935  // The incident energy E[0] is given in MeV.
1936  // The energy loss is calculated in MeV.
1937  // To avoid divergences as E->0, for any incident energy E<=1.e-3MeV (i.e. 1keV)
1938  // we return dE=E i.e. all particles with E<=1keV are stopped.
1939 
1940  return (E[0] - EResFunc(E, p));
1941 }
1942 
1943 
1944 
1952 
1954 {
1955  // Function parameterising the residual energy of charged particles after traversing this material.
1956  // This is simply an interface to the passage() function of the Range C library.
1957  // The incident energy E[0] is given in MeV.
1958  // The residual energy is calculated in MeV.
1959  // To avoid divergences as E->0, for any incident energy E<=1.e-3MeV (i.e. 1keV)
1960  // we return Eres=0 i.e. all particles with E<=1keV are stopped.
1961 
1962  if (E[0] < 1.e-3) return 0.0;
1963  return range_lib->passage(fTableType, Zp, Ap, iabso, fAbsorb[0].z, fAbsorb[0].a, E[0], thickness / KVUnits::mg, &error);
1964 }
1965 
1966 
1967 
1975 
1977 {
1978  // Function parameterising the range of charged particles in this material.
1979  // This is simply an interface to the rangen() function of the Range C library.
1980  // The incident energy E[0] is given in MeV.
1981  // The range is calculated in g/cm**2.
1982  // To avoid divergences as E->0, for any incident energy E<=1.e-3MeV (i.e. 1keV)
1983  // we return range=0 i.e. all particles with E<=1keV are stopped.
1984 
1985  if (E[0] < 1.e-3) return 0.;
1986  Double_t R = range_lib->rangen(fTableType, Zp, Ap, iabso, fAbsorb[0].z, fAbsorb[0].a, E[0]);
1987  return (R * KVUnits::mg);
1988 }
1989 
1990 
1991 
1993 
1995 {
1996  range_lib->nelem = fNelem;
1997  range_lib->is_gas = (int)IsGas(); // special treatment for effective charge in gases (M.F. Rivet, R. Bimbot et al)
1998  if (iabso < 0) {
1999  //cout << "nelem="<<nelem<<endl;
2000  for (int k = 0; k < fNelem; k++) {
2001  range_lib->absorb[k].z = fAbsorb[k].z;
2002  range_lib->absorb[k].a = fAbsorb[k].a;
2003  range_lib->absorb[k].w = fAbsorb[k].w;
2004  //cout << k << " " << absorb[k].z << " " << absorb[k].a << " " << absorb[k].w << endl;
2005  }
2006  }
2007  Zp = Z;
2008  Ap = A;
2009 }
2010 
2011 
2012 
2017 
2019 {
2020  // Return function giving energy loss (in MeV) as a function of incident energy (in MeV) for
2021  // charged particles (Z,A) traversing (or not) the thickness e (in g/cm**2) of this material.
2022  // isotopic mass isoAmat argument is not used.
2023 
2025  thickness = e;
2026  fDeltaE->SetRange(0, 400 * Ap);
2027  return fDeltaE;
2028 }
2029 
2030 
2031 
2036 
2038 {
2039  // Return function giving residual energy (in MeV) as a function of incident energy (in MeV) for
2040  // charged particles (Z,A) traversing (or not) the thickness e (in g/cm**2) of this material.
2041  // isotopic mass isoAmat argument is not used.
2042 
2044  thickness = e;
2045  fEres->SetRange(0, 400 * Ap);
2046  return fEres;
2047 }
2048 
2049 
2050 
2055 
2057 {
2058  // Return function giving range (in g/cm**2) as a function of incident energy (in MeV) for
2059  // charged particles (Z,A) traversing this material.
2060  // isotopic mass isoAmat argument is not used.
2061 
2063  // Yanez' range functions tend to have a (negative) minimum at very low energy,
2064  // below which the range diverges towards +infty at E=0.
2065  // Therefore we limit the range of the function to E>=Emin where Emin
2066  // is the energy for which the minimum occurs
2067  Double_t minX = fRange->GetMinimumX(1.e-6, 400 * Ap);
2068  fRange->SetRange(TMath::Max(1.e-6, minX), 400 * Ap);
2069  return fRange;
2070 }
2071 
2072 
2073 
2075 
2077 {
2079  if (IsCompound()) {
2080  fNelem = 0;
2081  iabso = -1;
2082  TIter next(fComposition);
2083  KVNameValueList* nvl;
2084  while ((nvl = (KVNameValueList*)next())) {
2085  fAbsorb[fNelem].z = nvl->GetIntValue("Z");
2086  fAbsorb[fNelem].a = nvl->GetIntValue("A");
2087  fAbsorb[fNelem].w = nvl->GetDoubleValue("Ar*Weight");
2088  fNelem++;
2089  }
2090  }
2091  else if (IsMixture()) {
2092  fNelem = 0;
2093  iabso = -1;
2094  TIter next(fComposition);
2095  KVNameValueList* nvl;
2096  while ((nvl = (KVNameValueList*)next())) {
2097  fAbsorb[fNelem].z = nvl->GetIntValue("Z");
2098  fAbsorb[fNelem].a = nvl->GetIntValue("A");
2099  fAbsorb[fNelem].w = nvl->GetDoubleValue("Ar*Weight");
2100  fNelem++;
2101  }
2102  }
2103 }
2104 
2105 
2106 
2111 
2113 {
2114  // Overrides KVIonRangeTableMaterial method to use the egassap() function of the Range C library.
2115  // Calculates incident energy (in MeV) of an ion (Z,A) with residual energy Eres (MeV) after thickness e (in g/cm**2).
2116  // isotopic mass isoAmat argument is not used.
2118  return range_lib->egassap(fTableType, Zp, Ap, iabso, fAbsorb[0].z, fAbsorb[0].a, e / KVUnits::mg, Eres, &error);
2119 }
2120 
2121 
2122 
2129 
2130 void KVRangeYanezMaterial::SaveMaterial(std::ofstream& matfile)
2131 {
2132  // Write definition of material in a file in the directory
2133  //
2134  // $(WORKING_DIR)/range_yanez/[name].dat
2135  //
2136  // All files in this directory are read when the table is initialised
2137 
2138  matfile << "// " << GetName() << " (" << GetSymbol() << ") generated ";
2139  TDatime now;
2140  matfile << now.AsString() << std::endl << std::endl;
2141  if (IsCompound()) matfile << "COMPOUND" << std::endl;
2142  else if (IsMixture()) matfile << "MIXTURE" << std::endl;
2143  else matfile << "ELEMENT" << std::endl;
2144  matfile << "name=" << GetName() << std::endl;
2145  matfile << "symbol=" << GetSymbol() << std::endl;
2146  matfile << "state=";
2147  if (IsGas()) matfile << "gas" << std::endl;
2148  else {
2149  matfile << "solid" << std::endl;
2150  matfile << "density=" << GetDensity() << std::endl;
2151  }
2152  if (IsCompound() || IsMixture()) {
2153  matfile << "nelem=" << GetNElem() << std::endl;
2154  TString listname = (IsCompound() ? "Compound element %d" : "Mixture element %d");
2155 
2156  for (int nel = 1; nel <= GetNElem(); ++nel) {
2157  KVNameValueList* compo = (KVNameValueList*)fComposition->FindObject(Form(listname.Data(), nel));
2158  KVNucleus nuc(compo->GetIntValue("Z"), compo->GetIntValue("A"));
2159  matfile << nuc.GetSymbol() << " " << compo->GetIntValue("Natoms");
2160  if (IsMixture()) matfile << " " << compo->GetDoubleValue("Proportion");
2161  matfile << std::endl;
2162  }
2163  }
2164  matfile << std::endl;
2165 }
2166 
2167 
int Int_t
ROOT::R::TRInterface & r
#define d(i)
#define f(i)
#define c(i)
#define e(i)
char Char_t
double Double_t
winID w
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t b
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
char name[80]
char * Form(const char *fmt,...)
virtual const Char_t * GetType() const
Definition: KVBase.h:176
Material for use in energy loss & range calculations.
TF1 * fDeltaE
function parameterising energy loss in material
const Char_t * GetSymbol() const
TF1 * fEres
function parameterising residual energy after crossing material
TF1 * fRange
function parameterising range of charged particles in material
void Copy(TObject &) const override
KVList * fComposition
composition of compound/mixture
Abstract base class for calculation of range & energy loss of charged particles in matter.
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Int_t GetIntValue(const Char_t *name) const
Double_t GetDoubleValue(const Char_t *name) const
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
const Char_t * GetSymbol(Option_t *opt="") const
Definition: KVNucleus.cpp:71
Description of absorber for the Range dE/dx and range library.
Int_t Ap
Z,A of incident projectile ion.
Double_t thickness
in g/cm**2
Int_t fNelem
number of elements in material
Int_t fTableType
=0 for Northcliffe-Schilling (<12 MeV/u), =1 for Hubert et al (2.5<E/A<500 MeV), =2 for interpolated ...
void PrepareRangeLibVariables(Int_t Z, Int_t A)
TF1 * GetRangeFunction(Int_t Z, Int_t A, Double_t isoAmat=0) override
TF1 * GetDeltaEFunction(Double_t e, Int_t Z, Int_t A, Double_t isoAmat=0) override
KVRangeYanezMaterial()
Default constructor.
Int_t iabso
value of iabso argument for function calls
Double_t RangeFunc(Double_t *, Double_t *)
void Copy(TObject &) const override
Double_t GetEIncFromEResOfIon(Int_t Z, Int_t A, Double_t Eres, Double_t e, Double_t isoAmat=0.) override
Double_t EResFunc(Double_t *, Double_t *)
std::unique_ptr< range > range_lib
Double_t DeltaEFunc(Double_t *, Double_t *)
Double_t error
calculated error in MeV
struct KVRangeYanezMaterial::Elem fAbsorb[10]
list of elements
TF1 * GetEResFunction(Double_t e, Int_t Z, Int_t A, Double_t isoAmat=0) override
void SaveMaterial(std::ofstream &matfile)
TObject * FindObject(const char *name) const override
const char * AsString() const
virtual Double_t GetMinimumX(Double_t xmin=0, Double_t xmax=0, Double_t epsilon=1.E-10, Int_t maxiter=100, Bool_t logx=false) const
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual void SetNpx(Int_t npx=100)
const char * GetName() const override
const char * Data() const
std::vector< elem > cmpnd
double polint(double *xa, double *ya, int n, double x, double *dy)
double ededxh(double e, int zp, int zt)
void splie2(double TOTO[], double x2a[], double *ya[], int m, int n, double *y2a[])
std::vector< elem > absorb
void gfact(double *le, int zt, double *f)
double thickn(int icorr, int zp, int ap, int iabso, int zt, int at, double ein, double delen)
void alref(double le, double *lz, double *dedx)
void splint(double xa[], double ya[], double y2a[], int n, double x, double *y)
void def_absorber(int zt, int at, int iabso)
void splin2(double x1a[], double x2a[], double *ya[], double *y2a[], int m, int n, double x1, double x2, double *y)
double rangen(int icorr, int zp, int ap, int iabso, int zt, int at, double ein)
double SMOOTHERSTEP(double x)
double dedx(int icorr, double e, int zp, int ap, int zt, int at)
const double FMT
double egassap(int icorr, int zp, int ap, int iabso, int zt, int at, double t, double eut, double *err)
void alion(int zp, double *dedxz2)
double ededx(double e, int zp, int zt)
void rangetab(int icorr, int zp, int ap, int iabso, int zt, int at, double *em, double *r, int *n)
void mpyers(double *le, int zt, double *f)
void dedxtab(int icorr, int zp, int ap, int iabso, int zt, int at, double e, double *tdedxe, double *tdedxn)
double ndedx(double e, int zp, int ap, int zt, int at)
double passage(int icorr, int zp, int ap, int iabso, int zt, int at, double ein, double t, double *err)
double s2az(double e, int zt)
void spline(double x[], double y[], int n, double yp1, double ypn, double *y2)
unsigned int locate(double y[], int n, double x)
Expr< UnaryOp< Sqrt< T >, SMatrix< T, D, D2, R >, T >, T, D, D2, R > sqrt(const SMatrix< T, D, D2, R > &rhs)
Expr< UnaryOp< Fabs< T >, SMatrix< T, D, D2, R >, T >, T, D, D2, R > fabs(const SMatrix< T, D, D2, R > &rhs)
RVec< PromoteTypes< T0, T1 > > pow(const T0 &x, const RVec< T1 > &v)
RVec< PromoteType< T > > log(const RVec< T > &v)
RVec< PromoteType< T > > log10(const RVec< T > &v)
RVec< PromoteType< T > > exp(const RVec< T > &v)
Double_t y[n]
Double_t x[n]
const Int_t n
TH1 * h
const long double mg
Definition: KVUnits.h:74
constexpr Double_t E()
constexpr Double_t R()
Double_t Max(Double_t a, Double_t b)
TMarker m
TArc a
ClassImp(TPyArg)