Version: SMASH-1.5
parametrizations.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2013-2018
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #include "smash/parametrizations.h"
12 
13 #include <cmath>
14 #include <initializer_list>
15 #include <iostream>
16 #include <memory>
17 #include <vector>
18 
19 #include "smash/average.h"
20 #include "smash/clebschgordan.h"
21 #include "smash/constants.h"
22 #include "smash/cxx14compat.h"
23 #include "smash/interpolation.h"
24 #include "smash/kinematics.h"
25 #include "smash/lowess.h"
26 #include "smash/pow.h"
27 
28 namespace smash {
29 
30 double xs_high_energy(double mandelstam_s, bool is_opposite_charge, double ma,
31  double mb, double P, double R1, double R2) {
32  const double M = 2.1206;
33  const double H = 0.272;
34  const double eta1 = 0.4473;
35  const double eta2 = 0.5486;
36  const double s_sab = mandelstam_s / (ma + mb + M) / (ma + mb + M);
37  double xs =
38  H * std::log(s_sab) * std::log(s_sab) + P + R1 * std::pow(s_sab, -eta1);
39  xs = is_opposite_charge ? xs + R2 * std::pow(s_sab, -eta2)
40  : xs - R2 * std::pow(s_sab, -eta2);
41  return xs;
42 }
43 
44 double pp_high_energy(double mandelstam_s) {
45  return xs_high_energy(mandelstam_s, false, 0.939, 0.939, 34.41, 13.07, 7.394);
46 }
47 
48 double ppbar_high_energy(double mandelstam_s) {
49  return xs_high_energy(mandelstam_s, true, 0.939, 0.939, 34.41, 13.07, 7.394);
50 }
51 
52 double np_high_energy(double mandelstam_s) {
53  return xs_high_energy(mandelstam_s, false, 0.939, 0.939, 34.41, 12.52, 6.66);
54 }
55 
56 double npbar_high_energy(double mandelstam_s) {
57  return xs_high_energy(mandelstam_s, true, 0.939, 0.939, 34.41, 12.52, 6.66);
58 }
59 
60 double piplusp_high_energy(double mandelstam_s) {
61  return xs_high_energy(mandelstam_s, false, 0.939, 0.138, 18.75, 9.56, 1.767);
62 }
63 
64 double piminusp_high_energy(double mandelstam_s) {
65  return xs_high_energy(mandelstam_s, true, 0.939, 0.138, 18.75, 9.56, 1.767);
66 }
67 
68 double xs_ppbar_annihilation(double mandelstam_s) {
69  const double xs_ref = 120.;
70  const double s_ref = 4. * nucleon_mass * nucleon_mass;
71  const double constant_a = 0.05;
72  const double constant_b = 0.6;
73  const double factor = constant_a * constant_a * s_ref /
74  ((mandelstam_s - s_ref) * (mandelstam_s - s_ref) +
75  constant_a * constant_a * s_ref) +
76  constant_b;
77  return xs_ref * (s_ref / mandelstam_s) * factor;
78 }
79 
80 double xs_string_hard(double mandelstam_s, double xs_0, double e_0,
81  double lambda_pow) {
82  const double sqrts = std::sqrt(mandelstam_s);
83  if (sqrts < e_0) {
84  return 0.;
85  } else {
86  double xs = xs_0 * std::pow(std::log(sqrts / e_0), lambda_pow);
87  return xs;
88  }
89 }
90 
91 double NN_string_hard(double mandelstam_s) {
92  return xs_string_hard(mandelstam_s, 0.087, 4.1, 3.8);
93 }
94 
95 double Npi_string_hard(double mandelstam_s) {
96  return xs_string_hard(mandelstam_s, 0.042, 3.5, 4.2);
97 }
98 
99 double pipi_string_hard(double mandelstam_s) {
100  return xs_string_hard(mandelstam_s, 0.013, 2.3, 4.7);
101 }
102 
103 /* pi+ p elastic cross section parametrization, PDG data.
104  *
105  * The PDG data is smoothed using the LOWESS algorithm. If more than one
106  * cross section was given for one p_lab value, the corresponding cross sections
107  * are averaged. */
108 static double piplusp_elastic_pdg(double mandelstam_s) {
109  if (piplusp_elastic_interpolation == nullptr) {
110  std::vector<double> x = PIPLUSP_ELASTIC_P_LAB;
111  std::vector<double> y = PIPLUSP_ELASTIC_SIG;
112  std::vector<double> dedup_x;
113  std::vector<double> dedup_y;
114  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
115  dedup_y = smooth(dedup_x, dedup_y, 0.1, 5);
117  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
118  }
119  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
120  return (*piplusp_elastic_interpolation)(p_lab);
121 }
122 
123 double piplusp_elastic_high_energy(double mandelstam_s, double m1, double m2) {
124  const double p_lab = (m1 > m2) ? plab_from_s(mandelstam_s, m2, m1)
125  : plab_from_s(mandelstam_s, m1, m2);
126  const auto logp = std::log(p_lab);
127  return 11.4 * std::pow(p_lab, -0.4) + 0.079 * logp * logp;
128 }
129 
130 double piplusp_elastic_AQM(double mandelstam_s, double m1, double m2) {
131  const double p_lab = (m1 > m2) ? plab_from_s(mandelstam_s, m2, m1)
132  : plab_from_s(mandelstam_s, m1, m2);
133  if (p_lab < 3.05) { // the plab from which the param starts to explode
134  return 7.5; // this will be scaled down by 2/3 for meson-meson
135  } else {
136  return piplusp_elastic_high_energy(mandelstam_s, m1, m2);
137  }
138 }
139 
140 double piplusp_elastic(double mandelstam_s) {
141  double sigma;
142  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
143  if (mandelstam_s < 2.25) {
144  sigma = really_small;
145  } else if (mandelstam_s > 4.84) {
146  const auto logp = std::log(p_lab);
147  sigma = 11.4 * std::pow(p_lab, -0.4) + 0.079 * logp * logp;
148  } else {
149  sigma = piplusp_elastic_pdg(mandelstam_s);
150  }
151 
152  // The elastic contributions from decays still need to be subtracted.
153  if (piplusp_elastic_res_interpolation == nullptr) {
154  std::vector<double> x = PIPLUSP_RES_SQRTS;
155  for (auto& i : x) {
156  i = i * i;
157  }
158  std::vector<double> y = PIPLUSP_RES_SIG;
160  make_unique<InterpolateDataSpline>(x, y);
161  }
162  // Interpolation may not be good if sqrts is beyond 2.2 GeV
163  if (mandelstam_s < 4.84) {
164  sigma -= (*piplusp_elastic_res_interpolation)(mandelstam_s);
165  }
166  if (sigma < 0) {
167  sigma = really_small;
168  }
169  return sigma;
170 }
171 
172 /* pi+ p to Sigma+ K+ cross section parametrization, PDG data.
173  *
174  * The PDG data is smoothed using the LOWESS algorithm. If more than one
175  * cross section was given for one p_lab value, the corresponding cross sections
176  * are averaged. */
177 double piplusp_sigmapluskplus_pdg(double mandelstam_s) {
178  if (piplusp_sigmapluskplus_interpolation == nullptr) {
179  std::vector<double> x = PIPLUSP_SIGMAPLUSKPLUS_P_LAB;
180  std::vector<double> y = PIPLUSP_SIGMAPLUSKPLUS_SIG;
181  std::vector<double> dedup_x;
182  std::vector<double> dedup_y;
183  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
184  dedup_y = smooth(dedup_x, dedup_y, 0.2, 5);
186  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
187  }
188  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
189  return (*piplusp_sigmapluskplus_interpolation)(p_lab);
190 }
191 
192 /* pi- p elastic cross section parametrization, PDG data.
193  *
194  * The PDG data is smoothed using the LOWESS algorithm. If more than one
195  * cross section was given for one p_lab value, the corresponding cross sections
196  * are averaged. */
197 static double piminusp_elastic_pdg(double mandelstam_s) {
198  if (piminusp_elastic_interpolation == nullptr) {
199  std::vector<double> x = PIMINUSP_ELASTIC_P_LAB;
200  std::vector<double> y = PIMINUSP_ELASTIC_SIG;
201  std::vector<double> dedup_x;
202  std::vector<double> dedup_y;
203  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
204  dedup_y = smooth(dedup_x, dedup_y, 0.03, 6);
206  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
207  }
208  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
209  return (*piminusp_elastic_interpolation)(p_lab);
210 }
211 
212 double piminusp_elastic(double mandelstam_s) {
213  double sigma;
214  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
215  const auto logp = std::log(p_lab);
216  if (mandelstam_s < 1.69) {
217  sigma = really_small;
218  } else if (mandelstam_s > 4.84) {
219  sigma = 1.76 + 11.2 * std::pow(p_lab, -0.64) + 0.043 * logp * logp;
220  } else {
221  sigma = piminusp_elastic_pdg(mandelstam_s);
222  }
223  /* Tune down the elastic cross section when sqrt s is between 1.8 GeV
224  * and 1.97 GeV so that the total cross section can fit the data. The
225  * scaling factor is chosen so that the it's equal to one and its
226  * derivate vanishes at the both ends. The minimum scaling factor in this
227  * region is 0.88-0.12=0.76. */
228  if (mandelstam_s > 3.24 && mandelstam_s < 3.8809) {
229  sigma *= (0.12 * cos(2 * M_PI * (sqrt(mandelstam_s) - 1.8) / (1.97 - 1.8)) +
230  0.88);
231  }
232  // The elastic contributions from decays still need to be subtracted.
233  if (piminusp_elastic_res_interpolation == nullptr) {
234  std::vector<double> x = PIMINUSP_RES_SQRTS;
235  for (auto& i : x) {
236  i = i * i;
237  }
238  std::vector<double> y = PIMINUSP_RES_SIG;
240  make_unique<InterpolateDataSpline>(x, y);
241  }
242  // Interpolation may not be good if sqrts is beyond 2.2 GeV
243  if (mandelstam_s < 4.84) {
244  sigma -= (*piminusp_elastic_res_interpolation)(mandelstam_s);
245  }
246  if (sigma < 0) {
247  sigma = really_small;
248  }
249  return sigma;
250 }
251 
252 /* pi- p -> Lambda K0 cross section parametrization, PDG data.
253  *
254  * The PDG data is smoothed using the LOWESS algorithm. If more than one
255  * cross section was given for one p_lab value, the corresponding cross sections
256  * are averaged. */
257 double piminusp_lambdak0_pdg(double mandelstam_s) {
258  if (piminusp_lambdak0_interpolation == nullptr) {
259  std::vector<double> x = PIMINUSP_LAMBDAK0_P_LAB;
260  std::vector<double> y = PIMINUSP_LAMBDAK0_SIG;
261  std::vector<double> dedup_x;
262  std::vector<double> dedup_y;
263  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
264  dedup_y = smooth(dedup_x, dedup_y, 0.2, 6);
266  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
267  }
268  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
269  return (*piminusp_lambdak0_interpolation)(p_lab);
270 }
271 
272 /* pi- p -> Sigma- K+ cross section parametrization, PDG data.
273  *
274  * The PDG data is smoothed using the LOWESS algorithm. If more than one
275  * cross section was given for one p_lab value, the corresponding cross sections
276  * are averaged. */
277 double piminusp_sigmaminuskplus_pdg(double mandelstam_s) {
279  std::vector<double> x = PIMINUSP_SIGMAMINUSKPLUS_P_LAB;
280  std::vector<double> y = PIMINUSP_SIGMAMINUSKPLUS_SIG;
281  std::vector<double> dedup_x;
282  std::vector<double> dedup_y;
283  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
284  dedup_y = smooth(dedup_x, dedup_y, 0.2, 6);
286  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
287  }
288  const double p_lab = plab_from_s(mandelstam_s, pion_mass, nucleon_mass);
290 }
291 
292 /* pi- p -> Sigma0 K0 cross section parametrization, resonance contribution.
293  *
294  * The data is smoothed using the LOWESS algorithm. If more than one
295  * cross section was given for one sqrts value, the corresponding cross sections
296  * are averaged. */
297 double piminusp_sigma0k0_res(double mandelstam_s) {
298  if (piminusp_sigma0k0_interpolation == nullptr) {
299  std::vector<double> x = PIMINUSP_SIGMA0K0_RES_SQRTS;
300  std::vector<double> y = PIMINUSP_SIGMA0K0_RES_SIG;
301  std::vector<double> dedup_x;
302  std::vector<double> dedup_y;
303  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
304  dedup_y = smooth(dedup_x, dedup_y, 0.03, 6);
306  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
307  }
308  const double sqrts = sqrt(mandelstam_s);
309  return (*piminusp_sigma0k0_interpolation)(sqrts);
310 }
311 
312 double pp_elastic(double mandelstam_s) {
313  const double p_lab = plab_from_s(mandelstam_s);
314  if (p_lab < 0.435) {
315  return 5.12 * nucleon_mass /
316  (mandelstam_s - 4 * nucleon_mass * nucleon_mass) +
317  1.67;
318  } else if (p_lab < 0.8) {
319  return 23.5 + 1000 * pow_int(p_lab - 0.7, 4);
320  } else if (p_lab < 2.0) {
321  return 1250 / (p_lab + 50) - 4 * (p_lab - 1.3) * (p_lab - 1.3);
322  } else if (p_lab < 2.776) {
323  return 77 / (p_lab + 1.5);
324  } else {
325  const auto logp = std::log(p_lab);
326  return 11.9 + 26.9 * std::pow(p_lab, -1.21) + 0.169 * logp * logp -
327  1.85 * logp;
328  }
329 }
330 
331 double pp_elastic_high_energy(double mandelstam_s, double m1, double m2) {
332  const double p_lab = (m1 > m2) ? plab_from_s(mandelstam_s, m2, m1)
333  : plab_from_s(mandelstam_s, m1, m2);
334  const auto logp = std::log(p_lab);
335  return 11.9 + 26.9 * std::pow(p_lab, -1.21) + 0.169 * logp * logp -
336  1.85 * logp;
337 }
338 
339 double pp_total(double mandelstam_s) {
340  const double p_lab = plab_from_s(mandelstam_s);
341  if (p_lab < 0.4) {
342  return 34 * std::pow(p_lab / 0.4, -2.104);
343  } else if (p_lab < 0.8) {
344  return 23.5 + 1000 * pow_int(p_lab - 0.7, 4);
345  } else if (p_lab < 1.5) {
346  return 23.5 + 24.6 / (1 + std::exp(-(p_lab - 1.2) / 0.1));
347  } else if (p_lab < 5.0) {
348  return 41 + 60 * (p_lab - 0.9) * std::exp(-1.2 * p_lab);
349  } else {
350  const auto logp = std::log(p_lab);
351  return 48.0 + 0.522 * logp * logp - 4.51 * logp;
352  }
353 }
354 
355 double np_elastic(double mandelstam_s) {
356  const double p_lab = plab_from_s(mandelstam_s);
357  if (p_lab < 0.525) {
358  return 17.05 * nucleon_mass /
359  (mandelstam_s - 4 * nucleon_mass * nucleon_mass) -
360  6.83;
361  } else if (p_lab < 0.8) {
362  return 33 + 196 * std::pow(std::abs(p_lab - 0.95), 2.5);
363  } else if (p_lab < 2.0) {
364  return 31 / std::sqrt(p_lab);
365  } else if (p_lab < 2.776) {
366  return 77 / (p_lab + 1.5);
367  } else {
368  const auto logp = std::log(p_lab);
369  return 11.9 + 26.9 * std::pow(p_lab, -1.21) + 0.169 * logp * logp -
370  1.85 * logp;
371  }
372 }
373 
374 double np_total(double mandelstam_s) {
375  const double p_lab = plab_from_s(mandelstam_s);
376  const auto logp = std::log(p_lab);
377  if (p_lab < 0.4) {
378  return 6.3555 * std::pow(p_lab, -3.2481) * std::exp(-0.377 * logp * logp);
379  } else if (p_lab < 1.0) {
380  return 33 + 196 * std::pow(std::abs(p_lab - 0.95), 2.5);
381  } else if (p_lab < 2.0) {
382  return 24.2 + 8.9 * p_lab;
383  } else if (p_lab < 5.0) {
384  return 42;
385  } else {
386  return 48.0 + 0.522 * logp * logp - 4.51 * logp;
387  }
388 }
389 
390 double ppbar_elastic(double mandelstam_s) {
391  const double p_lab = plab_from_s(mandelstam_s);
392  if (p_lab < 0.3) {
393  return 78.6;
394  } else if (p_lab < 5.0) {
395  return 31.6 + 18.3 / p_lab - 1.1 / (p_lab * p_lab) - 3.8 * p_lab;
396  } else {
397  const auto logp = std::log(p_lab);
398  return 10.2 + 52.7 * std::pow(p_lab, -1.16) + 0.125 * logp * logp -
399  1.28 * logp;
400  }
401 }
402 
403 double ppbar_total(double mandelstam_s) {
404  const double p_lab = plab_from_s(mandelstam_s);
405  if (p_lab < 0.3) {
406  return 271.6 * std::exp(-1.1 * p_lab * p_lab);
407  } else if (p_lab < 5.0) {
408  return 75.0 + 43.1 / p_lab + 2.6 / (p_lab * p_lab) - 3.9 * p_lab;
409  } else {
410  const auto logp = std::log(p_lab);
411  return 38.4 + 77.6 * std::pow(p_lab, -0.64) + 0.26 * logp * logp -
412  1.2 * logp;
413  }
414 }
415 
416 double deuteron_pion_elastic(double mandelstam_s) {
417  const double tmp = std::sqrt(mandelstam_s) - 2.172;
418  return 4.0 + 0.27 / (tmp * tmp + 0.065 * 0.065);
419 }
420 
421 double deuteron_nucleon_elastic(double mandelstam_s) {
422  const double s = mandelstam_s;
423  return 2500.0 * std::exp(-smash::square(s - 7.93) / 0.003) +
424  600.0 * std::exp(-smash::square(s - 7.93) / 0.1) + 10.0;
425 }
426 
427 double kplusp_elastic_background(double mandelstam_s) {
428  constexpr double a0 = 10.508; // mb
429  constexpr double a1 = -3.716; // mb/GeV
430  constexpr double a2 = 1.845; // mb/GeV^2
431  constexpr double a3 = -0.764; // GeV^-1
432  constexpr double a4 = 0.508; // GeV^-2
433 
434  const double p_lab = plab_from_s(mandelstam_s, kaon_mass, nucleon_mass);
435  const double p_lab2 = p_lab * p_lab;
436 
437  return (a0 + a1 * p_lab + a2 * p_lab2) / (1 + a3 * p_lab + a4 * p_lab2);
438 }
439 
440 double kplusn_elastic_background(double mandelstam_s) {
441  return 0.25 * kplusp_elastic_background(mandelstam_s);
442 }
443 
444 double kplusn_k0p(double mandelstam_s) {
445  return 0.25 * kplusp_elastic_background(mandelstam_s);
446 }
447 
448 /* K- p elastic cross section parametrization, PDG data.
449  *
450  * The PDG data is smoothed using the LOWESS algorithm. If more than one
451  * cross section was given for one p_lab value, the corresponding cross sections
452  * are averaged. */
453 static double kminusp_elastic_pdg(double mandelstam_s) {
454  if (kminusp_elastic_interpolation == nullptr) {
455  std::vector<double> x = KMINUSP_ELASTIC_P_LAB;
456  std::vector<double> y = KMINUSP_ELASTIC_SIG;
457  std::vector<double> dedup_x;
458  std::vector<double> dedup_y;
459  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
460  dedup_y = smooth(dedup_x, dedup_y, 0.1, 5);
462  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
463  }
464  const double p_lab = plab_from_s(mandelstam_s, kaon_mass, nucleon_mass);
465  return (*kminusp_elastic_interpolation)(p_lab);
466 }
467 
468 double kminusp_elastic_background(double mandelstam_s) {
469  const double p_lab = plab_from_s(mandelstam_s, kaon_mass, nucleon_mass);
470  double sigma;
471  if (std::sqrt(mandelstam_s) < 1.68) {
472  /* The parametrization here also works for anti-K0 n, Lambda pi0,
473  * Sigma+ pi-, Sigma- pi+, Sigma0 pi0 with different parameters a0, a1, a2.
474  *
475  * The values of the parameters are *not* taken from the source above,
476  * they come from a fit to PDG data. */
477  constexpr double a0 = 186.03567644; // mb GeV^2
478  constexpr double a1 = 0.22002795; // Gev
479  constexpr double a2 = 0.64907116;
480 
481  const double p_i = p_lab;
482  const double p_f = p_lab;
483 
484  const double ratio = a1 * a1 / (a1 * a1 + p_f * p_f);
485  sigma = a0 * p_f / (p_i * mandelstam_s) * std::pow(ratio, a2);
486  } else {
487  sigma = kminusp_elastic_pdg(mandelstam_s);
488  }
489  // The elastic contributions from decays still need to be subtracted.
490  if (kminusp_elastic_res_interpolation == nullptr) {
491  std::vector<double> x = KMINUSP_RES_SQRTS;
492  for (auto& i : x) {
493  i = plab_from_s(i * i, kaon_mass, nucleon_mass);
494  }
495  std::vector<double> y = KMINUSP_RES_SIG;
497  make_unique<InterpolateDataSpline>(x, y);
498  }
499  const auto old_sigma = sigma;
500  sigma -= (*kminusp_elastic_res_interpolation)(p_lab);
501  if (sigma < 0) {
502  std::cout << "NEGATIVE SIGMA: sigma=" << sigma
503  << ", sqrt(s)=" << std::sqrt(mandelstam_s)
504  << ", sig_el_exp=" << old_sigma
505  << ", sig_el_res=" << (*kminusp_elastic_res_interpolation)(p_lab)
506  << std::endl;
507  }
508  assert(sigma >= 0);
509  return sigma;
510 }
511 
512 double kminusn_elastic_background(double) { return 4.0; }
513 
514 double k0p_elastic_background(double mandelstam_s) {
515  // by isospin symmetry
516  return kplusn_elastic_background(mandelstam_s);
517 }
518 
519 double k0n_elastic_background(double mandelstam_s) {
520  // by isospin symmetry
521  return kplusp_elastic_background(mandelstam_s);
522 }
523 
524 double kbar0p_elastic_background(double mandelstam_s) {
525  // by isospin symmetry
526  return kminusn_elastic_background(mandelstam_s);
527 }
528 
529 double kbar0n_elastic_background(double mandelstam_s) {
530  // by isospin symmetry
531  return kminusp_elastic_background(mandelstam_s);
532 }
533 
534 double kplusp_inelastic_background(double mandelstam_s) {
535  if (kplusp_total_interpolation == nullptr) {
536  std::vector<double> x = KPLUSP_TOT_PLAB;
537  std::vector<double> y = KPLUSP_TOT_SIG;
538  std::vector<double> dedup_x;
539  std::vector<double> dedup_y;
540  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
541  dedup_y = smooth(dedup_x, dedup_y, 0.1, 5);
543  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
544  }
545  const double p_lab = plab_from_s(mandelstam_s, kaon_mass, nucleon_mass);
547  mandelstam_s);
548 }
549 
550 double kplusn_inelastic_background(double mandelstam_s) {
551  if (kplusn_total_interpolation == nullptr) {
552  std::vector<double> x = KPLUSN_TOT_PLAB;
553  std::vector<double> y = KPLUSN_TOT_SIG;
554  std::vector<double> dedup_x;
555  std::vector<double> dedup_y;
556  std::tie(dedup_x, dedup_y) = dedup_avg(x, y);
557  dedup_y = smooth(dedup_x, dedup_y, 0.05, 5);
559  make_unique<InterpolateDataLinear<double>>(dedup_x, dedup_y);
560  }
561  const double p_lab = plab_from_s(mandelstam_s, kaon_mass, nucleon_mass);
563  mandelstam_s) -
564  kplusn_k0p(mandelstam_s);
565 }
566 
575 static void initialize(std::unordered_map<std::pair<uint64_t, uint64_t>, double,
576  pair_hash>& ratios) {
577  const auto& type_p = ParticleType::find(pdg::p);
578  const auto& type_n = ParticleType::find(pdg::n);
579  const auto& type_K_p = ParticleType::find(pdg::K_p);
580  const auto& type_K_z = ParticleType::find(pdg::K_z);
581  const auto& type_Delta_pp = ParticleType::find(pdg::Delta_pp);
582  const auto& type_Delta_p = ParticleType::find(pdg::Delta_p);
583  const auto& type_Delta_z = ParticleType::find(pdg::Delta_z);
584  const auto& type_Delta_m = ParticleType::find(pdg::Delta_m);
585 
586  /* Store the isospin ratio of the given reaction relative to all other
587  * possible isospin-symmetric reactions. */
588  auto add_to_ratios = [&](const ParticleType& a, const ParticleType& b,
589  const ParticleType& c, const ParticleType& d,
590  double weight_numerator, double weight_other) {
591  assert(weight_numerator + weight_other != 0);
592  const auto key =
593  std::make_pair(pack(a.pdgcode().code(), b.pdgcode().code()),
594  pack(c.pdgcode().code(), d.pdgcode().code()));
595  const double ratio = weight_numerator / (weight_numerator + weight_other);
596  ratios[key] = ratio;
597  };
598 
599  /* All inelastic channels are K N -> K Delta -> K pi N or charge exchange,
600  * with identical cross section, weighted by the isospin factor.
601  *
602  * For charge exchange, the isospin factors are 1,
603  * so they are excluded here. */
604  {
605  const auto weight1 = isospin_clebsch_gordan_sqr_2to2(
606  type_p, type_K_p, type_K_z, type_Delta_pp);
607  const auto weight2 = isospin_clebsch_gordan_sqr_2to2(
608  type_p, type_K_p, type_K_p, type_Delta_p);
609 
610  add_to_ratios(type_p, type_K_p, type_K_z, type_Delta_pp, weight1, weight2);
611  add_to_ratios(type_p, type_K_p, type_K_p, type_Delta_p, weight2, weight1);
612  }
613  {
614  const auto weight1 = isospin_clebsch_gordan_sqr_2to2(
615  type_n, type_K_p, type_K_z, type_Delta_p);
616  const auto weight2 = isospin_clebsch_gordan_sqr_2to2(
617  type_n, type_K_p, type_K_p, type_Delta_z);
618 
619  add_to_ratios(type_n, type_K_p, type_K_z, type_Delta_p, weight1, weight2);
620  add_to_ratios(type_n, type_K_p, type_K_p, type_Delta_z, weight2, weight1);
621  }
622  /* K+ and K0 have the same mass and spin, their cross sections are assumed to
623  * only differ in isospin factors. */
624  {
625  const auto weight1 = isospin_clebsch_gordan_sqr_2to2(
626  type_p, type_K_z, type_K_z, type_Delta_p);
627  const auto weight2 = isospin_clebsch_gordan_sqr_2to2(
628  type_p, type_K_z, type_K_p, type_Delta_z);
629 
630  add_to_ratios(type_p, type_K_z, type_K_z, type_Delta_p, weight1, weight2);
631  add_to_ratios(type_p, type_K_z, type_K_p, type_Delta_z, weight2, weight1);
632  }
633  {
634  const auto weight1 = isospin_clebsch_gordan_sqr_2to2(
635  type_n, type_K_z, type_K_z, type_Delta_z);
636  const auto weight2 = isospin_clebsch_gordan_sqr_2to2(
637  type_n, type_K_z, type_K_p, type_Delta_m);
638 
639  add_to_ratios(type_n, type_K_z, type_K_z, type_Delta_z, weight1, weight2);
640  add_to_ratios(type_n, type_K_z, type_K_p, type_Delta_m, weight2, weight1);
641  }
642 }
643 
645  const ParticleType& b,
646  const ParticleType& c,
647  const ParticleType& d) const {
648  /* If this method is called with anti-nucleons, flip all particles to
649  * anti-particles;
650  * the ratio is equal */
651  int flip = 0;
652  for (const auto& p : {&a, &b, &c, &d}) {
653  if (p->is_nucleon()) {
654  if (flip == 0) {
655  flip = p->antiparticle_sign();
656  } else {
657  assert(p->antiparticle_sign() == flip);
658  }
659  }
660  }
661  const auto key = std::make_pair(
662  pack(a.pdgcode().code() * flip, b.pdgcode().code() * flip),
663  pack(c.pdgcode().code() * flip, d.pdgcode().code() * flip));
664  if (ratios_.empty()) {
666  }
667  return ratios_.at(key);
668 }
669 
670 /*thread_local (see #3075)*/ KaonNucleonRatios kaon_nucleon_ratios;
671 
672 double kminusp_kbar0n(double mandelstam_s) {
673  constexpr double a0 = 100; // mb GeV^2
674  constexpr double a1 = 0.15; // GeV
675  constexpr unsigned a2 = 2;
676 
677  const double p_lab = plab_from_s(mandelstam_s, kaon_mass, nucleon_mass);
678  const double p_i = p_lab;
679  const double p_f = p_lab;
680 
681  return a0 * p_f / (p_i * mandelstam_s) *
682  pow_int(a1 * a1 / (a1 * a1 + p_f * p_f), a2);
683 }
684 
685 double kminusp_piminussigmaplus(double sqrts) {
686  return 0.0788265 / smash::square(sqrts - 1.38841);
687 }
688 
689 double kminusp_piplussigmaminus(double sqrts) {
690  return 0.0196741 / smash::square(sqrts - 1.42318);
691 }
692 
693 double kminusp_pi0sigma0(double sqrts) {
694  return 0.0403364 / smash::square(sqrts - 1.39830305);
695 }
696 
697 double kminusp_pi0lambda(double sqrts) {
698  return 0.05932562 / smash::square(sqrts - 1.38786692);
699 }
700 
701 double kminusn_piminussigma0(double sqrts) {
702  return kminusp_piminussigmaplus(sqrts) + kminusp_piplussigmaminus(sqrts) -
703  2. * kminusp_pi0sigma0(sqrts);
704 }
705 
706 double kminusn_piminuslambda(double sqrts) {
707  return 2. * kminusp_pi0lambda(sqrts);
708 }
709 
710 // All K+ p and K+ n channels are forbidden by isospin.
711 
712 double lambdalambda_ximinusp(double sqrts_sqrts0, double p_N, double p_lambda) {
713  assert(p_lambda != 0);
714  assert(sqrts_sqrts0 >= 0);
715  return 37.15 / 2 * p_N / p_lambda * std::pow(sqrts_sqrts0, -0.16);
716 }
717 
718 double lambdalambda_xi0n(double sqrts_sqrts0, double p_N, double p_lambda) {
719  return lambdalambda_ximinusp(sqrts_sqrts0, p_N, p_lambda);
720 }
721 
722 double lambdasigmaplus_xi0p(double sqrts_sqrts0) {
723  assert(sqrts_sqrts0 >= 0);
724  return 24.3781 * std::pow(sqrts_sqrts0, -0.479);
725 }
726 
727 double lambdasigmaminus_ximinusn(double sqrts_sqrts0) {
728  return lambdasigmaplus_xi0p(sqrts_sqrts0);
729 }
730 
731 double lambdasigma0_ximinusp(double sqrts_sqrts0) {
732  assert(sqrts_sqrts0 >= 0);
733  if (sqrts_sqrts0 < 0.03336) {
734  return 6.475 * std::pow(sqrts_sqrts0, -0.4167);
735  } else {
736  return 14.5054 * std::pow(sqrts_sqrts0, -0.1795);
737  }
738 }
739 
740 double lambdasigma0_xi0n(double sqrts_sqrts0) {
741  return lambdasigma0_ximinusp(sqrts_sqrts0);
742 }
743 
744 double sigma0sigma0_ximinusp(double sqrts_sqrts0) {
745  assert(sqrts_sqrts0 >= 0);
746  if (sqrts_sqrts0 < 0.09047) {
747  return 5.625 * std::pow(sqrts_sqrts0, -0.318);
748  } else {
749  return 4.174 * std::pow(sqrts_sqrts0, -0.4421);
750  }
751 }
752 
753 double sigma0sigma0_xi0n(double sqrts_sqrts0) {
754  return sigma0sigma0_ximinusp(sqrts_sqrts0);
755 }
756 
757 double sigmaplussigmaminus_xi0p(double sqrts_sqrts0) {
758  return 4 * sigma0sigma0_ximinusp(sqrts_sqrts0);
759 }
760 
761 double sigma0sigmaminus_ximinusn(double sqrts_sqrts0) {
762  return 4 * sigma0sigma0_ximinusp(sqrts_sqrts0);
763 }
764 
765 double sigmaplussigmaminus_ximinusp(double sqrts_sqrts0) {
766  return 14.194 * std::pow(sqrts_sqrts0, -0.442);
767 }
768 
769 double sigmaplussigmaminus_xi0n(double sqrts_sqrts0) {
770  return sigmaplussigmaminus_ximinusp(sqrts_sqrts0);
771 }
772 
773 } // namespace smash
const std::initializer_list< double > KPLUSN_TOT_SIG
PDG data on K+ n total cross section: cross section.
double piminusp_elastic(double mandelstam_s)
pi-p elastic cross section parametrization Source: GiBUU:parametrizationBarMes_HighEnergy.f90
double xs_high_energy(double mandelstam_s, bool is_opposite_charge, double ma, double mb, double P, double R1, double R2)
total hadronic cross sections at high energies parametrized in the 2016 PDG book(http://pdg.lbl.gov/2016/reviews/rpp2016-rev-cross-section-plots.pdf)
const std::initializer_list< double > PIMINUSP_ELASTIC_P_LAB
PDG data on pi- p elastic cross section: momentum in lab frame.
const std::initializer_list< double > KPLUSP_TOT_PLAB
PDG data on K+ p total cross section: momentum in lab frame.
static std::unique_ptr< InterpolateDataLinear< double > > piminusp_sigma0k0_interpolation
An interpolation that gets lazily filled using the PIMINUSP_SIGMA0K0_RES data.
double lambdasigma0_ximinusp(double sqrts_sqrts0)
Lambda Sigma0 <-> Xi- p cross section parametrization Two hyperon exchange, based on effective model ...
double kminusp_kbar0n(double mandelstam_s)
K- p <-> Kbar0 n cross section parametrization.
double kminusp_elastic_background(double mandelstam_s)
K- p elastic background cross section parametrization Source: Buss:2011mx, B.3.9. ...
double sigma0sigma0_ximinusp(double sqrts_sqrts0)
Sigma0 Sigma0 <-> Xi- p cross section parametrization Two hyperon exchange, based on effective model ...
double get_ratio(const ParticleType &a, const ParticleType &b, const ParticleType &c, const ParticleType &d) const
Return the isospin ratio of the given K N -> K Delta cross section.
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:34
double deuteron_nucleon_elastic(double mandelstam_s)
Deuteron nucleon elastic cross-section [mb] parametrized by Oh:2009gx.
double piplusp_high_energy(double mandelstam_s)
pi+p total cross section at high energies
double piminusp_lambdak0_pdg(double mandelstam_s)
pi- p -> Lambda K0 cross section parametrization, PDG data.
double ppbar_elastic(double mandelstam_s)
ppbar elastic cross section parametrization Source: Bass:1998ca
double piminusp_sigma0k0_res(double mandelstam_s)
pi- p -> Sigma0 K0 cross section parametrization, resonance contribution.
double np_total(double mandelstam_s)
np total cross section parametrization Sources: low-p: Cugnon:1996kh highest-p: Buss:2011mx ...
double piminusp_sigmaminuskplus_pdg(double mandelstam_s)
pi- p -> Sigma- K+ cross section parametrization, PDG data.
double kminusp_piminussigmaplus(double sqrts)
K- p <-> pi- Sigma+ cross section parametrization Taken from UrQMD (Graef:2014mra).
constexpr int K_p
K⁺.
constexpr uint64_t pack(int32_t x, int32_t y)
Pack two int32_t into an uint64_t.
double kplusp_elastic_background(double mandelstam_s)
K+ p elastic background cross section parametrization.
Collection of useful constants that are known at compile time.
double pp_high_energy(double mandelstam_s)
pp total cross section at high energies
constexpr int Delta_m
Δ⁻.
double xs_string_hard(double mandelstam_s, double xs_0, double e_0, double lambda_pow)
Utility function called by specific other parametrizations Parametrized hard scattering cross section...
double kminusn_piminussigma0(double sqrts)
K- n <-> pi- Sigma0 cross section parametrization Follow from the parametrization with the same stran...
double k0n_elastic_background(double mandelstam_s)
K0 n elastic background cross section parametrization Source: Buss:2011mx, B.3.9. ...
double sigmaplussigmaminus_ximinusp(double sqrts_sqrts0)
Sigma+ Sigma- <-> Xi- p cross section parametrization Two hyperon exchange, based on effective model ...
double sigmaplussigmaminus_xi0p(double sqrts_sqrts0)
Sigma+ Sigma- <-> Xi0 p cross section parametrization Two hyperon exchange, based on effective model ...
const std::initializer_list< double > PIPLUSP_ELASTIC_P_LAB
PDG data on pi+ p elastic cross section: momentum in lab frame.
const std::initializer_list< double > PIPLUSP_RES_SIG
Elastic π⁺ N⁺ cross section contributions from decays.
constexpr T pow_int(const T base, unsigned const exponent)
Efficient template for calculating integer powers using squaring.
Definition: pow.h:23
double pp_elastic_high_energy(double mandelstam_s, double m1, double m2)
pp elastic cross section parametrization, with only the high energy part generalized to all energy re...
const std::initializer_list< double > PIMINUSP_RES_SIG
Elastic π⁻ N⁺ cross section contributions from decays.
double NN_string_hard(double mandelstam_s)
nucleon-nucleon hard scattering cross section (with partonic scattering)
double ppbar_high_energy(double mandelstam_s)
ppbar total cross section at high energies
constexpr double nucleon_mass
Nucleon mass in GeV.
Definition: constants.h:52
double k0p_elastic_background(double mandelstam_s)
K0 p elastic background cross section parametrization Source: Buss:2011mx, B.3.9. ...
const std::initializer_list< double > PIMINUSP_RES_SQRTS
Center-of-mass energy list for π⁻ N⁺
static std::unique_ptr< InterpolateDataLinear< double > > kplusn_total_interpolation
An interpolation that gets lazily filled using the KPLUSN_TOT data.
std::unordered_map< std::pair< uint64_t, uint64_t >, double, pair_hash > ratios_
Internal representation of isospin weights once calculated.
const std::initializer_list< double > PIMINUSP_LAMBDAK0_SIG
PDG data on pi- p to Lambda K0 cross section: cross section.
static std::unique_ptr< InterpolateDataLinear< double > > piplusp_elastic_interpolation
An interpolation that gets lazily filled using the PIPLUSP_ELASTIC_SIG data.
const std::initializer_list< double > PIPLUSP_SIGMAPLUSKPLUS_SIG
PDG data on pi+ p to Sigma+ K+ section: cross section.
double kplusn_k0p(double mandelstam_s)
K+ n charge exchange cross section parametrization.
KaonNucleonRatios kaon_nucleon_ratios
double plab_from_s(double mandelstam_s, double mass)
Convert Mandelstam-s to p_lab in a fixed-target collision.
Definition: kinematics.h:157
constexpr int Delta_z
Δ⁰.
const std::initializer_list< double > PIMINUSP_SIGMAMINUSKPLUS_P_LAB
PDG data on pi- p to Sigma- K+ cross section: momentum in lab frame.
constexpr double pion_mass
Pion mass in GeV.
Definition: constants.h:59
double kbar0n_elastic_background(double mandelstam_s)
Kbar0 n elastic background cross section parametrization Source: Buss:2011mx, B.3.9.
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
const std::initializer_list< double > PIMINUSP_SIGMAMINUSKPLUS_SIG
PDG data on pi- p to Sigma- K+ cross section: cross section.
constexpr int Delta_p
Δ⁺.
double kminusn_piminuslambda(double sqrts)
K- n <-> pi- Lambda cross section parametrization Follow from the parametrization with the same stran...
const std::initializer_list< double > PIMINUSP_ELASTIC_SIG
PDG data on pi- p elastic cross section: cross section.
double pipi_string_hard(double mandelstam_s)
pion-pion hard scattering cross section (with partonic scattering)
double kminusp_pi0sigma0(double sqrts)
K- p <-> pi0 Sigma0 cross section parametrization Fit to Landolt-Börnstein instead of UrQMD values...
constexpr int K_z
K⁰.
static std::unique_ptr< InterpolateDataSpline > kminusp_elastic_res_interpolation
An interpolation that gets lazily filled using the KMINUSP_RES data.
static double piplusp_elastic_pdg(double mandelstam_s)
constexpr int Delta_pp
Δ⁺⁺.
const std::initializer_list< double > KMINUSP_RES_SIG
Elastic K̅⁻ N⁺ cross section contributions from decays.
Particle type contains the static properties of a particle species.
Definition: particletype.h:87
double lambdasigma0_xi0n(double sqrts_sqrts0)
Lambda Sigma0 <-> Xi0 n cross section parametrization Two hyperon exchange, based on effective model ...
double lambdasigmaminus_ximinusn(double sqrts_sqrts0)
Lambda Sigma- <-> Xi- n cross section parametrization Two hyperon exchange, based on effective model ...
double piminusp_high_energy(double mandelstam_s)
pi-p total cross section at high energies
double np_elastic(double mandelstam_s)
np elastic cross section parametrization Source: Weil:2013mya, eq.
Calculate and store isospin ratios for K N -> K Delta reactions.
std::pair< std::vector< T >, std::vector< T > > dedup_avg(const std::vector< T > &x, const std::vector< T > &y)
Remove duplicates from data (x, y) by averaging y.
Definition: average.h:65
double lambdalambda_ximinusp(double sqrts_sqrts0, double p_N, double p_lambda)
Lambda Lambda <-> Xi- p cross section parametrization Two hyperon exchange, based on effective model ...
double kplusn_elastic_background(double mandelstam_s)
K+ n elastic background cross section parametrization sigma(K+n->K+n) = sigma(K+n->K0p) = 0...
static std::unique_ptr< InterpolateDataLinear< double > > piminusp_sigmaminuskplus_interpolation
An interpolation that gets lazily filled using the PIMINUSP_SIGMAMINUSKPLUS data. ...
std::vector< T > smooth(const std::vector< T > &x, const std::vector< T > &y, T span=2./3, size_t iter=3, T delta=0)
Apply the LOWESS smoother (see the reference below) to the given data (x, y).
Definition: lowess.h:289
constexpr T square(const T base)
Efficient template for calculating the square.
Definition: pow.h:38
double piplusp_elastic(double mandelstam_s)
pi+p elastic cross section parametrization, PDG data.
constexpr double kaon_mass
Kaon mass in GeV.
Definition: constants.h:66
const std::initializer_list< double > PIMINUSP_SIGMA0K0_RES_SIG
pi- p to Sigma0 K0 cross section: cross section
const std::initializer_list< double > PIPLUSP_RES_SQRTS
Center-of-mass energy list for π⁺ N⁺
static void initialize(std::unordered_map< std::pair< uint64_t, uint64_t >, double, pair_hash > &ratios)
Calculate and store isospin ratios for K N -> K Delta reactions.
static double kminusp_elastic_pdg(double mandelstam_s)
static std::unique_ptr< InterpolateDataLinear< double > > kminusp_elastic_interpolation
An interpolation that gets lazily filled using the KMINUSP_ELASTIC data.
double pp_elastic(double mandelstam_s)
pp elastic cross section parametrization Source: Weil:2013mya, eq.
const std::initializer_list< double > KMINUSP_ELASTIC_P_LAB
PDG data on K- p elastic cross section: momentum in lab frame.
Hash a pair of integers.
double kminusp_piplussigmaminus(double sqrts)
K- p <-> pi+ Sigma- cross section parametrization Taken from UrQMD (Graef:2014mra).
const std::initializer_list< double > KPLUSP_TOT_SIG
PDG data on K+ p total cross section: cross section.
const std::initializer_list< double > PIMINUSP_SIGMA0K0_RES_SQRTS
pi- p to Sigma0 K0 cross section: square root s
double sigma0sigmaminus_ximinusn(double sqrts_sqrts0)
Sigma0 Sigma- <-> Xi- n cross section parametrization Two hyperon exchange, based on effective model ...
double kminusn_elastic_background(double mandelstam_s)
K- n elastic background cross section parametrization Source: Buss:2011mx, B.3.9. ...
static double piminusp_elastic_pdg(double mandelstam_s)
const std::initializer_list< double > PIPLUSP_SIGMAPLUSKPLUS_P_LAB
PDG data on pi+ p to Sigma+ K+ cross section: momentum in lab frame.
double lambdasigmaplus_xi0p(double sqrts_sqrts0)
Lambda Sigma+ <-> Xi0 p cross section parametrization Two hyperon exchange, based on effective model ...
const std::initializer_list< double > KMINUSP_RES_SQRTS
Center-of-mass energy list for K̅⁻ N⁺
constexpr int p
Proton.
static std::unique_ptr< InterpolateDataLinear< double > > kplusp_total_interpolation
An interpolation that gets lazily filled using the KPLUSP_TOT data.
const std::initializer_list< double > PIPLUSP_ELASTIC_SIG
PDG data on pi+ p elastic cross section: cross section.
double isospin_clebsch_gordan_sqr_2to2(const ParticleType &p_a, const ParticleType &p_b, const ParticleType &p_c, const ParticleType &p_d, const int I=-1)
Calculate the squared isospin Clebsch-Gordan coefficient for a 2-to-2 reaction A + B -> C + D...
double kplusn_inelastic_background(double mandelstam_s)
K+ n inelastic background cross section parametrization Source: Buss:2011mx, B.3.8.
const std::initializer_list< double > PIMINUSP_LAMBDAK0_P_LAB
PDG data on pi- p to Lambda K0 cross section: momentum in lab frame.
double lambdalambda_xi0n(double sqrts_sqrts0, double p_N, double p_lambda)
Lambda Lambda <-> Xi0 n cross section parametrization Two hyperon exchange, based on effective model ...
double ppbar_total(double mandelstam_s)
ppbar total cross section parametrization Source: Bass:1998ca
static std::unique_ptr< InterpolateDataSpline > piminusp_elastic_res_interpolation
An interpolation that gets lazily filled using the PIMINUSP_RES data.
double kminusp_pi0lambda(double sqrts)
K- p <-> pi0 Lambda cross section parametrization Fit to Landolt-Börnstein instead of UrQMD values...
double piplusp_sigmapluskplus_pdg(double mandelstam_s)
pi+ p to Sigma+ K+ cross section parametrization, PDG data.
double sigma0sigma0_xi0n(double sqrts_sqrts0)
Sigma0 Sigma0 <-> Xi0 n cross section parametrization Two hyperon exchange, based on effective model ...
double npbar_high_energy(double mandelstam_s)
npbar total cross section at high energies
double piplusp_elastic_high_energy(double mandelstam_s, double m1, double m2)
pi+p elactic cross section parametrization.
const std::initializer_list< double > KMINUSP_ELASTIC_SIG
PDG data on K- p elastic cross section: cross section.
double deuteron_pion_elastic(double mandelstam_s)
Deuteron pion elastic cross-section [mb] parametrized to fit pi-d elastic scattering data (the data c...
double sigmaplussigmaminus_xi0n(double sqrts_sqrts0)
Sigma+ Sigma- <-> Xi0 n cross section parametrization Two hyperon exchange, based on effective model ...
static std::unique_ptr< InterpolateDataSpline > piplusp_elastic_res_interpolation
A null interpolation that gets filled using the PIPLUSP_RES data.
constexpr int n
Neutron.
double kbar0p_elastic_background(double mandelstam_s)
Kbar0 p elastic background cross section parametrization Source: Buss:2011mx, B.3.9.
const std::initializer_list< double > KPLUSN_TOT_PLAB
PDG data on K+ n total cross section: momentum in lab frame.
static std::unique_ptr< InterpolateDataLinear< double > > piminusp_elastic_interpolation
An interpolation that gets lazily filled using the PIMINUSP_ELASTIC data.
static std::unique_ptr< InterpolateDataLinear< double > > piplusp_sigmapluskplus_interpolation
An interpolation that gets lazily filled using the PIPLUSP_SIGMAPLUSKPLUS_SIG data.
double np_high_energy(double mandelstam_s)
np total cross section at high energies
PdgCode pdgcode() const
Definition: particletype.h:146
double pp_total(double mandelstam_s)
pp total cross section parametrization Sources: low-p: Cugnon:1996kh highest-p: Buss:2011mx ...
std::int32_t code() const
Definition: pdgcode.h:249
Definition: action.h:24
static std::unique_ptr< InterpolateDataLinear< double > > piminusp_lambdak0_interpolation
An interpolation that gets lazily filled using the PIMINUSP_LAMBDAK0 data.
double piplusp_elastic_AQM(double mandelstam_s, double m1, double m2)
An overload of piplusp_elastic_high_energy in which the very low part is replaced by a flat 5 mb cros...
double kplusp_inelastic_background(double mandelstam_s)
K+ p inelastic background cross section parametrization Source: Buss:2011mx, B.3.8.
double Npi_string_hard(double mandelstam_s)
nucleon-pion hard scattering cross section (with partonic scattering)
double xs_ppbar_annihilation(double mandelstam_s)
parametrized cross-section for proton-antiproton annihilation used in the UrQMD model ...