Version: SMASH-3.3
bremsstrahlungaction.cc
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2019-2020,2022,2024-2025
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
11 
13 #include "smash/outputinterface.h"
14 #include "smash/random.h"
15 
16 namespace smash {
17 static constexpr int LScatterAction = LogArea::ScatterAction::id;
18 
20  const ParticleList &in, const double time, const int n_frac_photons,
21  const double hadronic_cross_section_input,
22  const SpinInteractionType spin_interaction_type)
23  : ScatterAction(in[0], in[1], time),
24  reac_(bremsstrahlung_reaction_type(in)),
25  number_of_fractional_photons_(n_frac_photons),
26  hadronic_cross_section_(hadronic_cross_section_input),
27  spin_interaction_type_(spin_interaction_type) {}
28 
31  if (in.size() != 2) {
33  }
34 
35  PdgCode a = in[0].pdgcode();
36  PdgCode b = in[1].pdgcode();
37 
38  switch (pack(a.code(), b.code())) {
39  case (pack(pdg::pi_z, pdg::pi_m)):
40  case (pack(pdg::pi_m, pdg::pi_z)):
42 
43  case (pack(pdg::pi_z, pdg::pi_p)):
44  case (pack(pdg::pi_p, pdg::pi_z)):
46 
47  case (pack(pdg::pi_m, pdg::pi_p)):
48  case (pack(pdg::pi_p, pdg::pi_m)):
50 
51  case (pack(pdg::pi_m, pdg::pi_m)):
53 
54  case (pack(pdg::pi_p, pdg::pi_p)):
56 
57  case (pack(pdg::pi_z, pdg::pi_z)):
59 
60  default:
62  }
63 }
64 
65 void BremsstrahlungAction::perform_bremsstrahlung(const OutputsList &outputs) {
66  for (int i = 0; i < number_of_fractional_photons_; i++) {
68  for (const auto &output : outputs) {
69  if (output->is_photon_output()) {
70  // we do not care about the local density
71  output->at_interaction(*this, 0.0);
72  }
73  }
74  }
75 }
76 
78  // we have only one reaction per incoming particle pair
79  if (collision_processes_bremsstrahlung_.size() != 1) {
80  logg[LScatterAction].fatal()
81  << "Problem in BremsstrahlungAction::generate_final_state().\nThe "
82  "brocess branch has "
84  << " entries. It should however have 1.";
85  throw std::runtime_error("");
86  }
87 
88  auto *proc = collision_processes_bremsstrahlung_[0].get();
89 
90  outgoing_particles_ = proc->particle_list();
91  process_type_ = proc->get_type();
92  FourVector interaction_point = get_interaction_point();
93 
94  // Sample k and theta:
95  // minimum cutoff for k to be in accordance with cross section calculations
96  double delta_k; // k-range
97  double k_min = 0.001;
98  double k_max =
99  (sqrt_s() * sqrt_s() - 2 * outgoing_particles_[0].type().mass() * 2 *
100  outgoing_particles_[1].type().mass()) /
101  (2 * sqrt_s());
102 
103  if ((k_max - k_min) < 0.0) {
104  // Make sure it is kinematically even possible to create a photon that is
105  // in accordance with the cross section cutoff
106  k_ = 0.0;
107  delta_k = 0.0;
108  } else {
109  k_ = random::uniform(k_min, k_max);
110  delta_k = (k_max - k_min);
111  }
112  theta_ = random::uniform(0.0, M_PI);
113 
114  // Sample the phase space anisotropically in the local rest frame
116 
117  // Get differential cross sections
118  std::pair<double, double> diff_xs_pair = brems_diff_cross_sections();
119  double diff_xs_k = diff_xs_pair.first;
120  double diff_xs_theta = diff_xs_pair.second;
121 
122  // Assign weighting factor
123  const double W_theta = diff_xs_theta * (M_PI - 0.0);
124  const double W_k = diff_xs_k * delta_k;
125  weight_ = std::sqrt(W_theta * W_k) /
127 
128  // Scale weight by cross section scaling factor of incoming particles
129  weight_ *= incoming_particles_[0].xsec_scaling_factor() *
130  incoming_particles_[1].xsec_scaling_factor();
131 
132  // Set position and formation time and boost back to computational frame
133  for (auto &new_particle : outgoing_particles_) {
134  // assuming decaying particles are always fully formed
135  new_particle.set_formation_time(time_of_execution_);
136  new_particle.set_4position(interaction_point);
137  new_particle.boost_momentum(
139  }
140 
141  // Set unpolarized spin vectors
144  }
145 
146  // Photons are not really part of the normal processes, so we have to set a
147  // constant arbitrary number.
148  const auto id_process = ID_PROCESS_PHOTON;
149  Action::check_conservation(id_process);
150 }
151 
153  assert(outgoing_particles_.size() == 3);
154  const double m_a = outgoing_particles_[0].type().mass(),
155  m_b = outgoing_particles_[1].type().mass(),
156  m_c = outgoing_particles_[2].type().mass();
157  const double sqrts = sqrt_s();
158  const double E_ab = sqrts - m_c - k_; // Ekin of the pion pair in cm frame
159  const double pcm = pCM(sqrts, E_ab, m_c); // cm momentum of (π pair - photon)
160  const double pcm_pions = pCM(E_ab, m_a, m_b); // cm momentum within pion pair
161 
162  // Photon angle: Phi random, theta from theta_ sampled above
163  const Angles phitheta_photon(random::uniform(0.0, twopi), std::cos(theta_));
164  outgoing_particles_[2].set_4momentum(m_c, pcm * phitheta_photon.threevec());
165  // Boost velocity to cm frame of the two pions
166  const ThreeVector beta_cm_pion_pair_photon =
167  pcm * phitheta_photon.threevec() / std::sqrt(pcm * pcm + E_ab * E_ab);
168 
169  // Sample pion pair isotropically
170  Angles phitheta;
171  phitheta.distribute_isotropically();
172  outgoing_particles_[0].set_4momentum(m_a, pcm_pions * phitheta.threevec());
173  outgoing_particles_[1].set_4momentum(m_b, -pcm_pions * phitheta.threevec());
174  outgoing_particles_[0].boost_momentum(beta_cm_pion_pair_photon);
175  outgoing_particles_[1].boost_momentum(beta_cm_pion_pair_photon);
176 }
177 
179  double reaction_cross_section) {
180  CollisionBranchPtr dummy_process = std::make_unique<CollisionBranch>(
181  incoming_particles_[0].type(), incoming_particles_[1].type(),
182  reaction_cross_section, ProcessType::Bremsstrahlung);
183  add_collision(std::move(dummy_process));
184 }
185 
187  CollisionBranchList process_list;
188  // ParticleList final_state_particles;
189  static const ParticleTypePtr photon_particle =
191  static const ParticleTypePtr pi_z_particle = &ParticleType::find(pdg::pi_z);
192  static const ParticleTypePtr pi_p_particle = &ParticleType::find(pdg::pi_p);
193  static const ParticleTypePtr pi_m_particle = &ParticleType::find(pdg::pi_m);
194 
195  // Create interpolation objects, if not yet existent; only trigger for one
196  // of them as either all or none is created
199  }
200 
201  // Find cross section corresponding to given sqrt(s)
202  double sqrts = sqrt_s();
203  double xsection;
204 
206  // Here the final state is determined by the the final state provided by the
207  // sampled process using Monte Carlo techniqus
208 
209  // In the case of two oppositely charged pions as incoming particles,
210  // there are two potential final states: pi+ + pi- and pi0 + pi0
211  double xsection_pipi = (*pipi_pipi_opp_interpolation)(sqrts);
212  double xsection_pi0pi0 = (*pipi_pi0pi0_interpolation)(sqrts);
213 
214  // Prevent negative cross sections due to numerics in interpolation
215  xsection_pipi = (xsection_pipi <= 0.0) ? really_small : xsection_pipi;
216  xsection_pi0pi0 = (xsection_pi0pi0 <= 0.0) ? really_small : xsection_pi0pi0;
217 
218  // Necessary only to decide for a final state with pi+ and pi- as incoming
219  // particles.
220  CollisionBranchList process_list_pipi;
221 
222  // Add both processes to the process_list
223  process_list_pipi.push_back(std::make_unique<CollisionBranch>(
224  incoming_particles_[0].type(), incoming_particles_[1].type(),
225  *photon_particle, xsection_pipi, ProcessType::Bremsstrahlung));
226  process_list_pipi.push_back(std::make_unique<CollisionBranch>(
227  *pi_z_particle, *pi_z_particle, *photon_particle, xsection_pi0pi0,
229 
230  // Decide for one of the possible final states
231  double total_cross_section = xsection_pipi + xsection_pi0pi0;
232  const CollisionBranch *proc =
233  choose_channel<CollisionBranch>(process_list_pipi, total_cross_section);
234 
235  xsection = proc->weight();
236 
237  process_list.push_back(std::make_unique<CollisionBranch>(
238  proc->particle_list()[0].type(), proc->particle_list()[1].type(),
239  *photon_particle, xsection, ProcessType::Bremsstrahlung));
240 
241  } else if (reac_ == ReactionType::pi_m_pi_m ||
245  // Here the final state hadrons are identical to the initial state hadrons
247  xsection = (*pipi_pipi_same_interpolation)(sqrts);
248  } else {
249  // One pi0 in initial and final state
250  xsection = (*pipi0_pipi0_interpolation)(sqrts);
251  }
252 
253  // Prevent negative cross sections due to numerics in interpolation
254  xsection = (xsection <= 0.0) ? really_small : xsection;
255 
256  process_list.push_back(std::make_unique<CollisionBranch>(
257  incoming_particles_[0].type(), incoming_particles_[1].type(),
258  *photon_particle, xsection, ProcessType::Bremsstrahlung));
259 
260  } else if (reac_ == ReactionType::pi_z_pi_z) {
261  // Here we have a hard-coded final state that differs from the initial
262  // state, namely: pi0 + pi0 -> pi+- + pi-+ + gamma
263  xsection = (*pi0pi0_pipi_interpolation)(sqrts);
264 
265  // Prevent negative cross sections due to numerics in interpolation
266  xsection = (xsection <= 0.0) ? really_small : xsection;
267 
268  process_list.push_back(std::make_unique<CollisionBranch>(
269  *pi_p_particle, *pi_m_particle, *photon_particle, xsection,
271  } else {
272  throw std::runtime_error("Unknown ReactionType in BremsstrahlungAction.");
273  }
274 
275  return process_list;
276 }
277 
279  static const ParticleTypePtr pi_z_particle = &ParticleType::find(pdg::pi_z);
280  const double collision_energy = sqrt_s();
281  double dsigma_dk;
282  double dsigma_dtheta;
283 
285  if (outgoing_particles_[0].type() != *pi_z_particle) {
286  // pi+- + pi+-- -> pi+- + pi+- + gamma
287  dsigma_dk =
288  (*pipi_pipi_opp_dsigma_dk_interpolation)(k_, collision_energy);
289  dsigma_dtheta = (*pipi_pipi_opp_dsigma_dtheta_interpolation)(
290  theta_, collision_energy);
291  } else {
292  // pi+- + pi+-- -> pi0 + pi0 + gamma
293  dsigma_dk = (*pipi_pi0pi0_dsigma_dk_interpolation)(k_, collision_energy);
294  dsigma_dtheta =
295  (*pipi_pi0pi0_dsigma_dtheta_interpolation)(theta_, collision_energy);
296  }
297  } else if (reac_ == ReactionType::pi_p_pi_p ||
299  dsigma_dk = (*pipi_pipi_same_dsigma_dk_interpolation)(k_, collision_energy);
300  dsigma_dtheta =
301  (*pipi_pipi_same_dsigma_dtheta_interpolation)(theta_, collision_energy);
302  } else if (reac_ == ReactionType::pi_z_pi_p ||
304  dsigma_dk = (*pipi0_pipi0_dsigma_dk_interpolation)(k_, collision_energy);
305  dsigma_dtheta =
306  (*pipi0_pipi0_dsigma_dtheta_interpolation)(theta_, collision_energy);
307  } else if (reac_ == ReactionType::pi_z_pi_z) {
308  dsigma_dk = (*pi0pi0_pipi_dsigma_dk_interpolation)(k_, collision_energy);
309  dsigma_dtheta =
310  (*pi0pi0_pipi_dsigma_dtheta_interpolation)(theta_, collision_energy);
311  } else {
312  throw std::runtime_error(
313  "Unkown channel when computing differential cross sections for "
314  "bremsstrahlung processes.");
315  }
316 
317  // Prevent negative cross sections due to numerics in interpolation
318  dsigma_dk = (dsigma_dk < 0.0) ? really_small : dsigma_dk;
319  dsigma_dtheta = (dsigma_dtheta < 0.0) ? really_small : dsigma_dtheta;
320 
321  // Combine differential cross sections to a pair
322  std::pair<double, double> diff_x_sections = {dsigma_dk, dsigma_dtheta};
323 
324  return diff_x_sections;
325 }
326 
328  // Read in tabularized values for sqrt(s), k and theta
329  std::vector<double> sqrts = BREMS_SQRTS;
330  std::vector<double> photon_momentum = BREMS_K;
331  std::vector<double> photon_angle = BREMS_THETA;
332 
333  // Read in tabularized total cross sections
334  std::vector<double> sigma_pipi_pipi_opp = BREMS_PIPI_PIPI_OPP_SIG;
335  std::vector<double> sigma_pipi_pipi_same = BREMS_PIPI_PIPI_SAME_SIG;
336  std::vector<double> sigma_pipi0_pipi0 = BREMS_PIPI0_PIPI0_SIG;
337  std::vector<double> sigma_pipi_pi0pi0 = BREMS_PIPI_PI0PI0_SIG;
338  std::vector<double> sigma_pi0pi0_pipi = BREMS_PI0PI0_PIPI_SIG;
339 
340  // Read in tabularized differential cross sections dSigma/dk
341  std::vector<double> dsigma_dk_pipi_pipi_opp = BREMS_PIPI_PIPI_OPP_DIFF_SIG_K;
342  std::vector<double> dsigma_dk_pipi_pipi_same =
344  std::vector<double> dsigma_dk_pipi0_pipi0 = BREMS_PIPI0_PIPI0_DIFF_SIG_K;
345  std::vector<double> dsigma_dk_pipi_pi0pi0 = BREMS_PIPI_PI0PI0_DIFF_SIG_K;
346  std::vector<double> dsigma_dk_pi0pi0_pipi = BREMS_PI0PI0_PIPI_DIFF_SIG_K;
347 
348  // Read in tabularized differential cross sections dSigma/dtheta
349  std::vector<double> dsigma_dtheta_pipi_pipi_opp =
351  std::vector<double> dsigma_dtheta_pipi_pipi_same =
353  std::vector<double> dsigma_dtheta_pipi0_pipi0 =
355  std::vector<double> dsigma_dtheta_pipi_pi0pi0 =
357  std::vector<double> dsigma_dtheta_pi0pi0_pipi =
359 
360  // Create interpolation objects containing linear interpolations for
361  // total cross sections
362  pipi_pipi_opp_interpolation = std::make_unique<InterpolateDataLinear<double>>(
363  sqrts, sigma_pipi_pipi_opp);
365  std::make_unique<InterpolateDataLinear<double>>(sqrts,
366  sigma_pipi_pipi_same);
368  std::make_unique<InterpolateDataLinear<double>>(sqrts, sigma_pipi0_pipi0);
370  std::make_unique<InterpolateDataLinear<double>>(sqrts, sigma_pipi_pi0pi0);
372  std::make_unique<InterpolateDataLinear<double>>(sqrts, sigma_pi0pi0_pipi);
373 
374  // Create interpolation objects containing bicubic interpolations for
375  // differential dSigma/dk
377  std::make_unique<InterpolateData2DSpline>(photon_momentum, sqrts,
378  dsigma_dk_pipi_pipi_opp);
380  std::make_unique<InterpolateData2DSpline>(photon_momentum, sqrts,
381  dsigma_dk_pipi_pipi_same);
383  std::make_unique<InterpolateData2DSpline>(photon_momentum, sqrts,
384  dsigma_dk_pipi0_pipi0);
386  std::make_unique<InterpolateData2DSpline>(photon_momentum, sqrts,
387  dsigma_dk_pipi_pi0pi0);
389  std::make_unique<InterpolateData2DSpline>(photon_momentum, sqrts,
390  dsigma_dk_pi0pi0_pipi);
391 
392  // Create interpolation objects containing bicubic interpolations for
393  // differential dSigma/dtheta
395  std::make_unique<InterpolateData2DSpline>(photon_angle, sqrts,
396  dsigma_dtheta_pipi_pipi_opp);
398  std::make_unique<InterpolateData2DSpline>(photon_angle, sqrts,
399  dsigma_dtheta_pipi_pipi_same);
401  std::make_unique<InterpolateData2DSpline>(photon_angle, sqrts,
402  dsigma_dtheta_pipi0_pipi0);
404  std::make_unique<InterpolateData2DSpline>(photon_angle, sqrts,
405  dsigma_dtheta_pipi_pi0pi0);
407  std::make_unique<InterpolateData2DSpline>(photon_angle, sqrts,
408  dsigma_dtheta_pi0pi0_pipi);
409 }
410 } // namespace smash
FourVector total_momentum_of_outgoing_particles() const
Calculate the total kinetic momentum of the outgoing particles.
Definition: action.cc:160
ParticleList outgoing_particles_
Initially this stores only the PDG codes of final-state particles.
Definition: action.h:372
const double time_of_execution_
Time at which the action is supposed to be performed (absolute time in the lab frame in fm).
Definition: action.h:378
void assign_unpolarized_spin_vector_to_outgoing_particles()
Assign an unpolarized spin vector to all outgoing particles.
Definition: action.cc:478
virtual double check_conservation(const uint32_t id_process) const
Check various conservation laws.
Definition: action.cc:484
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
Definition: action.h:271
ParticleList incoming_particles_
List with data of incoming particles.
Definition: action.h:364
FourVector get_interaction_point() const
Get the interaction point.
Definition: action.cc:68
ProcessType process_type_
type of process
Definition: action.h:381
Angles provides a common interface for generating directions: i.e., two angles that should be interpr...
Definition: angles.h:59
ThreeVector threevec() const
Definition: angles.h:288
void distribute_isotropically()
Populate the object with a new direction.
Definition: angles.h:199
void sample_3body_phasespace()
Sample the final state anisotropically, considering the differential cross sections with respect to t...
ReactionType
Enum for encoding the photon process.
const ReactionType reac_
Reaction process as determined from incoming particles.
double k_
Sampled value of k (photon momentum)
const SpinInteractionType spin_interaction_type_
Type of spin interaction to use.
const int number_of_fractional_photons_
Number of photons created for each hadronic scattering, needed for correct weighting.
void add_dummy_hadronic_process(double reaction_cross_section)
Adds one hadronic process with a given cross-section.
void generate_final_state() override
Generate the final-state for the Bremsstrahlung process.
static ReactionType bremsstrahlung_reaction_type(const ParticleList &in)
Determine photon process from incoming particles.
void perform_bremsstrahlung(const OutputsList &outputs)
Create the final state and write to output.
CollisionBranchList brems_cross_sections()
Computes the total cross section of the bremsstrahlung process.
BremsstrahlungAction(const ParticleList &in, const double time, const int n_frac_photons, const double hadronic_cross_section_input, const SpinInteractionType spin_interaction_type=SpinInteractionType::Off)
Construct a ScatterActionBrems object.
double hadronic_cross_section() const
Return the total cross section of the underlying hadronic scattering It is necessary for the weightin...
CollisionBranchList collision_processes_bremsstrahlung_
Holds the bremsstrahlung branch.
double weight_
Weight of the produced photon.
void create_interpolations()
Create interpolation objects for tabularized cross sections: total cross section, differential dSigma...
double theta_
Sampled value of theta (angle of the photon)
std::pair< double, double > brems_diff_cross_sections()
Computes the differential cross sections dSigma/dk and dSigma/dtheta of the bremsstrahlung process.
CollisionBranch is a derivative of ProcessBranch, which is used to represent particular final-state c...
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
Definition: fourvector.h:33
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
Definition: fourvector.h:333
A pointer-like interface to global references to ParticleType objects.
Definition: particletype.h:682
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
Definition: particletype.cc:99
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
Definition: pdgcode.h:108
std::int32_t code() const
Definition: pdgcode.h:319
ParticleList particle_list() const
double weight() const
ScatterAction is a special action which takes two incoming particles and performs a scattering,...
Definition: scatteraction.h:30
void add_collision(CollisionBranchPtr p)
Add a new collision channel.
The ThreeVector class represents a physical three-vector with the components .
Definition: threevector.h:31
SpinInteractionType
Possible spin interaction types.
@ Off
No spin interactions.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:40
constexpr int pi_p
π⁺.
constexpr int pi_z
π⁰.
constexpr int photon
Photon.
constexpr int pi_m
π⁻.
T uniform(T min, T max)
Definition: random.h:88
Definition: action.h:24
static std::unique_ptr< InterpolateData2DSpline > pipi_pi0pi0_dsigma_dk_interpolation
const std::initializer_list< double > BREMS_PIPI_PIPI_SAME_DIFF_SIG_K
dSigma/dk for π+ + π+ -> π+ + π+ + γ or π- + π- -> π- + π- + γ
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
Definition: kinematics.h:79
const std::initializer_list< double > BREMS_PIPI_PI0PI0_SIG
Total π+- + π-+ -> π0 + π0 + γ cross section.
static std::unique_ptr< InterpolateDataLinear< double > > pipi_pipi_opp_interpolation
const std::initializer_list< double > BREMS_PIPI_PIPI_OPP_SIG
Total π+- + π-+ -> π+- + π-+ + γ cross section.
const std::initializer_list< double > BREMS_PIPI_PI0PI0_DIFF_SIG_THETA
dSigma/dtheta for π+- + π-+ -> π0 + π0 + γ
const std::initializer_list< double > BREMS_PIPI0_PIPI0_DIFF_SIG_K
dSigma/dk for π0 + π -> π0 + π + γ
const std::initializer_list< double > BREMS_PIPI_PIPI_SAME_SIG
Total π+ + π+ -> π+ + π+ + γ or π- + π- -> π- + π- + γ cross section.
static std::unique_ptr< InterpolateDataLinear< double > > pipi_pipi_same_interpolation
static std::unique_ptr< InterpolateData2DSpline > pipi0_pipi0_dsigma_dtheta_interpolation
constexpr std::uint32_t ID_PROCESS_PHOTON
Process ID for any photon process.
Definition: constants.h:122
static std::unique_ptr< InterpolateData2DSpline > pipi_pipi_opp_dsigma_dk_interpolation
static std::unique_ptr< InterpolateData2DSpline > pi0pi0_pipi_dsigma_dk_interpolation
static std::unique_ptr< InterpolateDataLinear< double > > pipi0_pipi0_interpolation
@ Bremsstrahlung
See here for a short description.
const std::initializer_list< double > BREMS_PIPI_PIPI_OPP_DIFF_SIG_THETA
dSigma/dtheta for π+- + π-+ -> π+- + π-+ + γ
static std::unique_ptr< InterpolateData2DSpline > pipi0_pipi0_dsigma_dk_interpolation
const std::initializer_list< double > BREMS_PI0PI0_PIPI_SIG
Total π0 + π0 -> π+- + π-+ + γ cross section.
const std::initializer_list< double > BREMS_K
photon momentum
constexpr double twopi
.
Definition: constants.h:49
const std::initializer_list< double > BREMS_SQRTS
Center-of-mass energy.
const std::initializer_list< double > BREMS_THETA
theta angle with respect to collision axis of incoming pions
const std::initializer_list< double > BREMS_PI0PI0_PIPI_DIFF_SIG_THETA
dSigma/dtheta for π0 + π0 -> π+- + π-+ + γ
static std::unique_ptr< InterpolateDataLinear< double > > pipi_pi0pi0_interpolation
static std::unique_ptr< InterpolateData2DSpline > pipi_pipi_same_dsigma_dtheta_interpolation
constexpr uint64_t pack(int32_t x, int32_t y)
Pack two int32_t into an uint64_t.
const std::initializer_list< double > BREMS_PIPI0_PIPI0_SIG
Total π0 + π -> π0 + π + γ cross section.
const std::initializer_list< double > BREMS_PI0PI0_PIPI_DIFF_SIG_K
dSigma/dk for π0 + π0 -> π+- + π-+ + γ
const std::initializer_list< double > BREMS_PIPI_PI0PI0_DIFF_SIG_K
dSigma/dk for π+- + π-+ -> π0 + π0 + γ
static std::unique_ptr< InterpolateDataLinear< double > > pi0pi0_pipi_interpolation
static std::unique_ptr< InterpolateData2DSpline > pipi_pipi_same_dsigma_dk_interpolation
static std::unique_ptr< InterpolateData2DSpline > pipi_pi0pi0_dsigma_dtheta_interpolation
static std::unique_ptr< InterpolateData2DSpline > pi0pi0_pipi_dsigma_dtheta_interpolation
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:41
static constexpr int LScatterAction
const std::initializer_list< double > BREMS_PIPI_PIPI_OPP_DIFF_SIG_K
dSigma/dk for π+- + π-+ -> π+- + π-+ + γ
static std::unique_ptr< InterpolateData2DSpline > pipi_pipi_opp_dsigma_dtheta_interpolation
const std::initializer_list< double > BREMS_PIPI0_PIPI0_DIFF_SIG_THETA
dSigma/dtheta for π0 + π -> π0 + π + γ
const std::initializer_list< double > BREMS_PIPI_PIPI_SAME_DIFF_SIG_THETA
dSigma/dtheta for π+ + π+ -> π+ + π+ + γ or π- + π- -> π- + π- + γ