14 #include "Pythia8/Pythia.h"
30 bool isotropic,
double string_formation_time,
31 double box_length,
bool is_total_parametrized)
32 :
Action({in_part_a, in_part_b}, time),
33 sum_of_partial_cross_sections_(0.),
34 isotropic_(isotropic),
35 string_formation_time_(string_formation_time),
36 is_total_parametrized_(is_total_parametrized) {
37 box_length_ = box_length;
38 if (is_total_parametrized_) {
39 parametrized_total_cross_section_ = NAN;
100 "ScatterAction::generate_final_state: Invalid process type " +
101 std::to_string(
static_cast<int>(
process_type_)) +
" was requested. " +
108 new_particle.boost_momentum(
112 new_particle.set_4position(middle_point);
121 throw std::logic_error(
122 "add_all_scatterings should be called only once per ScatterAction "
129 CollisionBranchList processes =
143 const double xs_diff =
154 if (pseudoresonance && finder_parameters.
two_to_one) {
162 auto pseudoresonance_branch = std::make_unique<CollisionBranch>(
168 << pseudoresonance->
name() <<
" with cross section " << xs_gap
182 <<
"Trying to rescale branches before adding processes.";
183 throw std::logic_error(
184 "This function can only be called after having added processes.");
191 const std::set<ParticleTypePtr> pair{type_a, type_b};
194 <<
"Total cross section between " << type_a->
name() <<
" and "
195 << type_b->
name() <<
"is roughly zero at sqrt(s) = " <<
sqrt_s()
196 <<
" GeV, and no rescaling to match the parametrized value will be "
197 "done.\nAn elastic process will be added, instead, to match the "
198 "total cross section.\nFor this pair of particles, this warning "
199 "will be subsequently suppressed.";
202 auto elastic_branch = std::make_unique<CollisionBranch>(
207 const double reweight =
212 proc->set_weight(proc->weight() * reweight);
222 const double desired_mass =
sqrt_s();
231 if (nucleon_and_pion) {
233 }
else if (two_nucleons) {
235 }
else if (two_pions) {
237 }
else if (nucleon_and_kaon) {
261 if (std::empty(list)) {
267 auto largest = *std::max_element(list.begin(), list.end(),
269 return a->mass() < b->mass();
275 return std::abs(a->
mass() - desired_mass) <
276 std::abs(b->mass() - desired_mass);
278 auto closest = *std::min_element(list.begin(), list.end(), comparison);
281 throw std::logic_error(
"Unknown method for selecting pseudoresonance.");
295 <<
"Trying to parametrize total cross section when it shouldn't be.";
296 throw std::logic_error(
297 "This function can only be called on ScatterAction objects with "
298 "parametrized cross section.");
318 return (1. / std::sqrt(1.0 -
beta_cm().sqr()));
339 const double lamb =
lambda_tilde(m_s, m1 * m1, m2 * m2);
358 " position difference [fm]: ", pos_diff,
359 ", momentum difference [GeV]: ", mom_diff);
361 const double dp2 = mom_diff.
sqr();
362 const double dr2 = pos_diff.
sqr();
367 const double dpdr = pos_diff * mom_diff;
376 const double result = dr2 - dpdr * dpdr / dp2;
377 return result > 0.0 ? result : 0.0;
386 const double mom_diff_sqr =
388 const double x_sqr = delta_x.
sqr();
396 const double p_a_dot_x = p_a.
momentum().
Dot(delta_x);
397 const double p_b_dot_x = p_b.
momentum().
Dot(delta_x);
402 (p_a_sqr * std::pow(p_b_dot_x, 2) + p_b_sqr * std::pow(p_a_dot_x, 2) -
403 2 * p_a_dot_p_b * p_a_dot_x * p_b_dot_x) /
404 (std::pow(p_a_dot_p_b, 2) - p_a_sqr * p_b_sqr);
405 return b_sqr > 0.0 ? b_sqr : 0.0;
418 return 7.6 + 0.66 * std::log(mandelstam_s);
436 return 5.5 * p8 / (7.7 + p8);
455 }
else if (plab < 0.6) {
456 return 16.53 * (plab - 0.225);
457 }
else if (plab < 1.6) {
458 return -1.63 * plab + 7.16;
465 double kinetic_energy_cm) {
488 const double mass_a = masses.first;
489 const double mass_b = masses.second;
491 const std::array<double, 2> t_range = get_t_range<double>(
492 kinetic_energy_cm, mass_in_a, mass_in_b, mass_a, mass_b);
498 double mandelstam_s_new = 0.;
514 double bb, a, plab =
plab_from_s(mandelstam_s_new);
521 a = (plab < 0.8) ? 1. : 0.64 / (plab * plab);
531 t = t_range[0] + t_range[1] - t;
535 1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
548 t = t_range[0] + t_range[1] - t;
551 1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
555 const std::array<double, 4>
p{1.46434, 5.80311, -6.89358, 1.94302};
556 const double a =
p[0] + mass_a * (
p[1] + mass_a * (
p[2] + mass_a *
p[3]));
560 double t = t_range[0];
565 t = t_range[0] + t_range[1] - t;
568 1. - 2. * (t - t_range[0]) / (t_range[1] - t_range[0]));
581 const double p_f =
pCM(kinetic_energy_cm, mass_a, mass_b);
584 " radial momentum: ", p_f);
624 "resonance_formation: "
625 "Incorrect number of particles in final state: ";
649 out_mom += data.momentum();
652 logg[
LPythia].debug(
"Outgoing momenta string:", out_mom);
667 bool success =
false;
669 const int ntry_max = 10000;
670 while (!success && ntry < ntry_max) {
697 logg[
LPythia].error(
"Unknown string process required.");
701 if (ntry == ntry_max) {
706 bool success_newtry =
false;
722 while (!success_newtry && ntry_new < ntry_max) {
731 if (success_newtry) {
735 if (!success_newtry) {
756 out <<
" (not performed)";
Action is the base class for a generic process that takes a number of incoming particles and transfor...
void sample_2body_phasespace()
Sample the full 2-body phase-space (masses, momenta, angles) in the center-of-mass frame for the fina...
FourVector total_momentum_of_outgoing_particles() const
Calculate the total kinetic momentum of the outgoing particles.
void assign_formation_time_to_outgoing_particles()
Assign the formation time to the outgoing particles.
std::pair< FourVector, FourVector > get_potential_at_interaction_point() const
Get the skyrme and asymmetry potential at the interaction point.
ParticleList outgoing_particles_
Initially this stores only the PDG codes of final-state particles.
FourVector total_momentum() const
Sum of 4-momenta of incoming particles.
virtual void sample_manybody_phasespace()
Sample the full n-body phase-space (masses, momenta, angles) in the center-of-mass frame for the fina...
const double time_of_execution_
Time at which the action is supposed to be performed (absolute time in the lab frame in fm).
const ParticleList & incoming_particles() const
Get the list of particles that go into the action.
double sqrt_s() const
Determine the total energy in the center-of-mass frame [GeV].
ParticleList incoming_particles_
List with data of incoming particles.
FourVector get_interaction_point() const
Get the interaction point.
static double lambda_tilde(double a, double b, double c)
Little helper function that calculates the lambda function (sometimes written with a tilde to better ...
ProcessType process_type_
type of process
Angles provides a common interface for generating directions: i.e., two angles that should be interpr...
ThreeVector threevec() const
void distribute_isotropically()
Populate the object with a new direction.
CollisionBranch is a derivative of ProcessBranch, which is used to represent particular final-state c...
ProcessType get_type() const override
The cross section class assembels everything that is needed to calculate the cross section and return...
double high_energy(const ScatterActionsFinderParameters &finder_parameters) const
Determine the parametrized total cross section at high energies for the given collision,...
double string_probability(const ScatterActionsFinderParameters &finder_parameters) const
CollisionBranchList string_excitation(double total_string_xs, StringProcess *string_process, const ScatterActionsFinderParameters &finder_parameters) const
Determine the cross section for string excitations, which is given by the difference between the para...
double parametrized_total(const ScatterActionsFinderParameters &finder_parameters) const
Select the parametrization for the total cross section, given the types of incoming particles.
CollisionBranchList generate_collision_list(const ScatterActionsFinderParameters &finder_parameters, StringProcess *string_process) const
Generate a list of all possible collisions between the incoming particles with the given c....
Guard type that safely disables floating point traps for the scope in which it is placed.
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
double sqr() const
calculate the square of the vector (which is a scalar)
ThreeVector threevec() const
double Dot(const FourVector &a) const
calculate the scalar product with another four-vector
ThreeVector velocity() const
Get the velocity (3-vector divided by zero component).
ParticleData contains the dynamic information of a certain particle.
PdgCode pdgcode() const
Get the pdgcode of the particle.
void set_4momentum(const FourVector &momentum_vector)
Set the particle's 4-momentum directly.
const ParticleType & type() const
Get the type of the particle.
const FourVector & momentum() const
Get the particle's 4-momentum.
void boost(const ThreeVector &v)
Apply a full Lorentz boost of momentum and position.
const FourVector & position() const
Get the particle's position in Minkowski space.
A pointer-like interface to global references to ParticleType objects.
const std::string & name() const
int32_t charge() const
The charge of the particle.
bool is_Deltastar() const
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
int antiparticle_sign() const
int baryon_number() const
ParticleList particle_list() const
Thrown when ScatterAction is called to perform with unknown ProcessType.
bool isotropic_
Do this collision isotropically?
ScatterAction(const ParticleData &in_part1, const ParticleData &in_part2, double time, bool isotropic=false, double string_formation_time=1.0, double box_length=-1.0, bool is_total_parametrized=false)
Construct a ScatterAction object.
std::optional< double > parametrized_total_cross_section_
If cross section is parametrized, store the value.
void rescale_outgoing_branches()
Loop over the possible branches and rescales their weight according to the desired total cross sectio...
void add_collision(CollisionBranchPtr p)
Add a new collision channel.
void resonance_formation()
Perform a 2->1 resonance-formation process.
ThreeVector beta_cm() const
Get the velocity of the center of mass of the scattering/incoming particles in the calculation frame.
double partial_cross_section_
Partial cross-section to the chosen outgoing channel.
double relative_velocity() const
Get the relative velocity of the two incoming particles.
double mandelstam_s() const
Determine the Mandelstam s variable,.
void set_parametrized_total_cross_section(const ScatterActionsFinderParameters &finder_parameters)
Given the incoming particles, assigns the correct parametrization of the total cross section.
double get_partial_weight() const override
Get the partial cross section of the chosen channel.
StringProcess * string_process_
Pointer to interface class for strings.
bool is_total_parametrized_
Whether the total cross section is parametrized.
void create_string_final_state()
Creates the final states for string-processes after they are performed.
double sum_of_partial_cross_sections_
Current sum of partial hadronic cross sections.
double cm_momentum_squared() const
Get the squared momentum of the center of mass of the incoming particles in the calculation frame.
void generate_final_state() override
Generate the final-state of the scattering process.
void add_all_scatterings(const ScatterActionsFinderParameters &finder_parameters)
Add all possible scattering subprocesses for this action object.
void two_to_many_scattering()
Perform an inelastic two-to-many-body scattering (more than 2)
double get_total_weight() const override
Get the total cross section of scattering particles.
void string_excitation()
Todo(ryu): document better - it is not really UrQMD-based, isn't it? Perform the UrQMD-based string e...
void elastic_scattering()
Perform an elastic two-body scattering, i.e. just exchange momentum.
void inelastic_scattering()
Perform an inelastic two-body scattering, i.e. new particles are formed.
void sample_angles(std::pair< double, double > masses, double kinetic_energy_cm) override
Sample final-state angles in a 2->2 collision (possibly anisotropic).
void add_collisions(CollisionBranchList pv)
Add several new collision channels at once.
double gamma_cm() const
Get the gamma factor corresponding to a boost to the center of mass frame of the colliding particles.
double cm_momentum() const
Get the momentum of the center of mass of the incoming particles in the calculation frame.
bool were_processes_added_
Lock for calling add_all_scatterings only once.
static std::set< std::set< ParticleTypePtr > > warned_no_rescaling_available
Warn about zero cross section only once per particle type pair.
double cov_transverse_distance_sqr() const
Calculate the transverse distance of the two incoming particles in their local rest frame written in ...
ParticleTypePtr try_find_pseudoresonance(const PseudoResonance method, const StringTransitionParameters &transition) const
Try to find a pseudo-resonance that can be created from the incoming particles using a given method.
CollisionBranchList collision_channels_
List of possible collisions.
double transverse_distance_sqr() const
Calculate the transverse distance of the two incoming particles in their local rest frame.
Helper class for ScatterActionsFinder.
const bool strings_with_probability
This indicates whether the string fragmentation is swiched on with a probability smoothly increasing ...
const StringTransitionParameters transition_high_energy
Constants related to transition between low collision energies - mediated via resonances - and high c...
const PseudoResonance pseudoresonance_method
Which pseudo-resonance to choose.
const bool two_to_one
Enables resonance production.
bool next_SDiff(bool is_AB_to_AX)
Single-diffractive process is based on single pomeron exchange described in Ingelman:1984ns .
bool next_NDiffSoft()
Soft Non-diffractive process is modelled in accordance with dual-topological approach Capella:1978ig ...
bool next_DDiff()
Double-diffractive process ( A + B -> X + X ) is similar to the single-diffractive process,...
ParticleList get_final_state()
a function to get the final state particle list which is called after the collision
bool next_BBbarAnn()
Baryon-antibaryon annihilation process Based on what UrQMD Bass:1998ca , Bleicher:1999xi does,...
void init(const ParticleList &incoming, double tcoll)
initialization feed intial particles, time of collision and gamma factor of the center of mass.
bool next_NDiffHard()
Hard Non-diffractive process is based on PYTHIA 8 with partonic showers and interactions.
The ThreeVector class represents a physical three-vector with the components .
void rotate_z_axis_to(ThreeVector &r)
Rotate the z-axis onto the vector r.
Collection of useful constants that are known at compile time.
PseudoResonance
Which pseudo-resonance fills the inelastic gap in the transition to string region of cross sections.
@ Closest
Resonance with the pole mass closest from the invariant mass of incoming particles for all processes.
@ ClosestFromUnstable
Closest resonance for a given mass from processes with at least one resonance in the incoming particl...
@ None
No pseudo-resonance is created.
@ LargestFromUnstable
Heaviest possible resonance from processes with at least one resonance in the incoming particles.
@ Largest
Resonance of largest mass for all processes.
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
void format_debug_output(std::ostream &out) const override
Writes information about this scatter action to the out stream.
T power(T n, T xMin, T xMax)
Draws a random number according to a power-law distribution ~ x^n.
T expo(T A, T x1, T x2)
Draws a random number x from an exponential distribution exp(A*x), where A is assumed to be positive,...
double plab_from_s(double mandelstam_s, double mass)
Convert Mandelstam-s to p_lab in a fixed-target collision.
static double Cugnon_bnp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic np...
T pCM(const T sqrts, const T mass_a, const T mass_b) noexcept
ParticleTypePtrList list_possible_resonances(const ParticleTypePtr type_a, const ParticleTypePtr type_b)
Lists the possible resonances that decay into two particles.
T pCM_sqr(const T sqrts, const T mass_a, const T mass_b) noexcept
@ FailedString
See here for a short description.
@ TwoToOne
See here for a short description.
@ StringSoftDoubleDiffractive
See here for a short description.
@ TwoToFive
See here for a short description.
@ StringSoftSingleDiffractiveXB
See here for a short description.
@ TwoToTwo
See here for a short description.
@ Elastic
See here for a short description.
@ TwoToFour
See here for a short description.
@ StringSoftAnnihilation
See here for a short description.
@ StringSoftNonDiffractive
See here for a short description.
@ StringSoftSingleDiffractiveAX
See here for a short description.
@ StringHard
See here for a short description.
@ TwoToThree
See here for a short description.
static double high_energy_bpp(double plab)
Computes the B coefficients from the STAR fit, see fig.
constexpr double nucleon_mass
Nucleon mass in GeV.
constexpr T pow_int(const T base, unsigned const exponent)
Efficient template for calculating integer powers using squaring.
constexpr double pion_mass
Pion mass in GeV.
T pCM_from_s(const T s, const T mass_a, const T mass_b) noexcept
constexpr double really_small
Numerical error tolerance.
static constexpr int LScatterAction
static constexpr int LPythia
bool is_string_soft_process(ProcessType p)
Check if a given process type is a soft string excitation.
double s_from_plab(double plab, double m_P, double m_T)
Convert p_lab to Mandelstam-s for a fixed-target setup, with a projectile of mass m_P and momentum pl...
static double Cugnon_bpp(double plab)
Computes the B coefficients from the Cugnon parametrization of the angular distribution in elastic pp...
Constants related to transition between low and high collision energies.
const std::pair< double, double > sqrts_range_Npi
Transition range in N collisions.
const double pipi_offset
Constant offset as to where to turn on the strings and elastic processes for reactions (this is an e...
const double sqrts_add_lower
Constant for the lower end of transition region in the case of AQM this is added to the sum of masses...
const double KN_offset
Constant offset as to where to shift from 2to2 to string processes (in GeV) in the case of KN reactio...
const std::pair< double, double > sqrts_range_NN
Transition range in NN collisions.