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