10 #ifndef SRC_INCLUDE_PROCESSSTRING_H_ 11 #define SRC_INCLUDE_PROCESSSTRING_H_ 18 #include "Pythia8/Pythia.h" 66 std::array<FourVector, 2>
plab_;
68 std::array<FourVector, 2>
pcom_;
234 double gluon_beta,
double gluon_pmin,
235 double quark_alpha,
double quark_beta,
236 double strange_supp,
double diquark_supp,
237 double sigma_perp,
double stringz_a_leading,
238 double stringz_b_leading,
double stringz_a,
239 double stringz_b,
double string_sigma_T,
240 double factor_t_form,
241 bool mass_dependent_formation_times,
242 double prob_proton_to_d_uu,
243 bool separate_fragment_baryon,
double popcorn_rate);
266 double diquark_supp,
double popcorn_rate,
267 double stringz_a,
double stringz_b,
268 double string_sigma_T);
277 const auto &log = logger<LogArea::Pythia>();
281 pythia_hadron_->rndm.init(seed_new);
282 log.debug(
"pythia_hadron_ : rndm is initialized with seed ", seed_new);
305 double sqrts_threshold = 2. * (1. + 1.0e-6);
308 const int pdg_a_mod =
309 (std::abs(pdg_a) > 1000) ? pdg_a : 10 * (std::abs(pdg_a) / 10) + 3;
310 const int pdg_b_mod =
311 (std::abs(pdg_b) > 1000) ? pdg_b : 10 * (std::abs(pdg_b) / 10) + 3;
312 sqrts_threshold += pythia_hadron_->particleData.m0(pdg_a_mod) +
313 pythia_hadron_->particleData.m0(pdg_b_mod);
316 if (sqrt_s < sqrts_threshold) {
317 sqrt_s = sqrts_threshold;
319 pythia_sigmatot_.calc(pdg_a, pdg_b, sqrt_s);
320 return {pythia_sigmatot_.sigmaAX(), pythia_sigmatot_.sigmaXB(),
321 pythia_sigmatot_.sigmaXX()};
339 pmin_gluon_lightcone_ = p_light_cone_min;
358 pow_fquark_alpha_ = alphapow;
359 pow_fquark_beta_ = betapow;
371 kappa_tension_string_ = kappa_string;
383 void init(
const ParticleList &incoming,
double tcoll);
393 std::array<ThreeVector, 3> &evec_basis);
409 const std::array<std::array<int, 2>, 2> &quarks,
410 const std::array<FourVector, 2> &pstr_com,
411 std::array<double, 2> &m_str,
412 std::array<ThreeVector, 2> &evec_str);
425 const std::array<std::array<int, 2>, 2> &quarks,
426 const std::array<FourVector, 2> &pstr_com,
427 const std::array<double, 2> &m_str,
428 const std::array<ThreeVector, 2> &evec_str,
429 bool flip_string_ends,
bool separate_fragment_baryon);
499 std::array<int, 5> &excess_quark,
500 std::array<int, 5> &excess_antiq);
523 std::array<int, 5> &excess_constituent);
539 std::array<int, 5> &nquark_total,
540 std::array<int, 5> &nantiq_total);
568 std::array<int, 5> &nquark_total, std::array<int, 5> &nantiq_total,
569 bool sign_constituent,
570 std::array<std::array<int, 5>, 2> &excess_constituent);
600 std::array<std::array<int, 5>, 2> &excess_quark,
601 std::array<std::array<int, 5>, 2> &excess_antiq);
641 std::array<std::array<int, 5>, 2> &excess_quark,
642 std::array<std::array<int, 5>, 2> &excess_antiq);
669 Pythia8::Event &event_intermediate,
670 Pythia8::Event &event_hadronize);
703 Pythia8::Event &event_intermediate,
704 Pythia8::Event &event_hadronize);
728 Pythia8::Event &event_intermediate,
729 Pythia8::Event &event_hadronize);
748 Pythia8::Event &event) {
750 for (
int ip = 2; ip <
event.size() - np_end; ip++) {
751 const double y_quark_current =
event[ip].y();
752 const double y_quark_forward =
event[iforward].y();
753 if ((find_forward && y_quark_current > y_quark_forward) ||
754 (!find_forward && y_quark_current < y_quark_forward)) {
793 ParticleList &intermediate_particles) {
794 const std::string s = std::to_string(pdgid);
798 intermediate_particles.push_back(new_particle);
807 if (pythia_id == 310 || pythia_id == 130) {
845 return Pythia8::Vec4(mom.
x1(), mom.
x2(), mom.
x3(), energy);
860 std::array<ThreeVector, 3> &evec_basis) {
861 ThreeVector three_momentum = evec_basis[0] * particle.pz() +
862 evec_basis[1] * particle.px() +
863 evec_basis[2] * particle.py();
864 return FourVector(particle.e(), three_momentum);
886 bool separate_fragment_baryon,
887 ParticleList &intermediate_particles);
926 bool separate_fragment_baryon,
927 std::array<ThreeVector, 3> &evec_basis,
928 double &ppos_string,
double &pneg_string,
929 double &QTrx_string_pos,
double &QTrx_string_neg,
930 double &QTry_string_pos,
double &QTry_string_neg,
931 Pythia8::FlavContainer &flav_string_pos,
932 Pythia8::FlavContainer &flav_string_neg,
933 std::vector<int> &pdgid_frag,
934 std::vector<FourVector> &momentum_frag);
976 bool separate_fragment_hadron,
977 double ppos_string,
double pneg_string,
978 double mTrn_had_forward,
double mTrn_had_backward,
979 double &ppos_had_forward,
double &ppos_had_backward,
980 double &pneg_had_forward,
double &pneg_had_backward);
1012 Pythia8::Event &event_fragments, std::array<ThreeVector, 3> &evec_basis,
1013 double ppos_string,
double pneg_string,
1014 double QTrx_string,
double QTry_string,
1015 double QTrx_add_pos,
double QTry_add_pos,
1016 double QTrx_add_neg,
double QTry_add_neg);
1029 Pythia8::Event &event, std::array<ThreeVector, 3> &evec_basis,
1030 double factor_yrapid,
double diff_yrapid) {
1031 Pythia8::Vec4 pvec_string_now = Pythia8::Vec4(0., 0., 0., 0.);
1033 for (
int ipyth = 1; ipyth <
event.size(); ipyth++) {
1034 if (!event[ipyth].isFinal()) {
1039 event[ipyth].e(), event[ipyth].px(),
1040 event[ipyth].py(), event[ipyth].pz());
1041 const double E_frag = p_frag.
x0();
1042 const double pl_frag = p_frag.
threevec() * evec_basis[0];
1043 const double ppos_frag = (E_frag + pl_frag) * M_SQRT1_2;
1044 const double pneg_frag = (E_frag - pl_frag) * M_SQRT1_2;
1045 const double mTrn_frag = std::sqrt(2. * ppos_frag * pneg_frag);
1047 const double y_frag = 0.5 * std::log(ppos_frag / pneg_frag);
1050 const double y_new_frag = factor_yrapid * y_frag + diff_yrapid;
1052 const double E_new_frag = mTrn_frag * std::cosh(y_new_frag);
1053 const double pl_new_frag = mTrn_frag * std::sinh(y_new_frag);
1055 p_frag.
threevec() + (pl_new_frag - pl_frag) * evec_basis[0];
1056 Pythia8::Vec4 pvec_new_frag =
set_Vec4(E_new_frag, mom_new_frag);
1057 event[ipyth].p(pvec_new_frag);
1058 pvec_string_now += pvec_new_frag;
1060 event[0].p(pvec_string_now);
1061 event[0].m(pvec_string_now.mCalc());
1080 ParticleList& outgoing_particles,
1082 double suppression_factor);
1098 static std::pair<int, int>
find_leading(
int nq1,
int nq2, ParticleList& list);
1113 double suppression_factor);
1139 #endif // SRC_INCLUDE_PROCESSSTRING_H_ double sample_zLund(double a, double b, double mTrn)
Sample lightcone momentum fraction according to the LUND fragmentation function.
double additional_xsec_supp_
additional cross-section suppression factor to take coherence effect into account.
bool remake_kinematics_fragments(Pythia8::Event &event_fragments, std::array< ThreeVector, 3 > &evec_basis, double ppos_string, double pneg_string, double QTrx_string, double QTry_string, double QTrx_add_pos, double QTry_add_pos, double QTrx_add_neg, double QTry_add_neg)
std::array< int, 2 > NpartString_
number of particles fragmented from strings
bool next_NDiffSoft()
Soft Non-diffractive process is modelled in accordance with dual-topological approach Capella:1978ig...
The ThreeVector class represents a physical three-vector with the components .
double PPosA_
forward lightcone momentum p^{+} of incoming particle A in CM-frame [GeV]
static FourVector reorient(Pythia8::Particle &particle, std::array< ThreeVector, 3 > &evec_basis)
compute the four-momentum properly oriented in the lab frame.
static Pythia8::Vec4 set_Vec4(double energy, const ThreeVector &mom)
Easy setter of Pythia Vec4 from SMASH.
double pmin_gluon_lightcone_
the minimum lightcone momentum scale carried by a gluon [GeV]
String excitation processes used in SMASH.
bool make_final_state_2strings(const std::array< std::array< int, 2 >, 2 > &quarks, const std::array< FourVector, 2 > &pstr_com, const std::array< double, 2 > &m_str, const std::array< ThreeVector, 2 > &evec_str, bool flip_string_ends, bool separate_fragment_baryon)
Prepare kinematics of two strings, fragment them and append to final_state.
void set_pmin_gluon_lightcone(double p_light_cone_min)
set the minimum lightcone momentum scale carried by gluon.
void set_sigma_qperp_(double sigma_qperp)
set the average amount of transverse momentum transfer sigma_qperp_.
int NpartFinal_
total number of final state particles
static void make_string_ends(const PdgCode &pdgcode_in, int &idq1, int &idq2, double xi)
make a random selection to determine partonic contents at the string ends.
int fragment_string(int idq1, int idq2, double mString, ThreeVector &evecLong, bool flip_string_ends, bool separate_fragment_baryon, ParticleList &intermediate_particles)
perform string fragmentation to determine species and momenta of hadrons by implementing PYTHIA 8...
Collection of useful constants that are known at compile time.
void replace_constituent(Pythia8::Particle &particle, std::array< int, 5 > &excess_constituent)
Convert a partonic PYTHIA partice into the desired species and update the excess of constituents...
std::unique_ptr< Pythia8::Pythia > pythia_hadron_
PYTHIA object used in fragmentation.
std::array< FourVector, 2 > plab_
momenta of incoming particles in the lab frame [GeV]
int append_final_state(ParticleList &intermediate_particles, const FourVector &uString, const ThreeVector &evecLong)
compute the formation time and fill the arrays with final-state particles as described in Andersson:1...
static int pdg_map_for_pythia(PdgCode &pdg)
Take pdg code and map onto particle specie which can be handled by PYTHIA.
std::array< PdgCode, 2 > PDGcodes_
PdgCodes of incoming particles.
double stringz_b_leading_
parameter (StringZ:bLund) for the fragmentation function of leading baryon in soft non-diffractive st...
void find_total_number_constituent(Pythia8::Event &event_intermediate, std::array< int, 5 > &nquark_total, std::array< int, 5 > &nantiq_total)
Compute how many quarks and antiquarks we have in the system, and update the correspoing arrays with ...
ThreeVector threevec() const
double prob_proton_to_d_uu_
Probability of splitting a nucleon into the quark flavour it has only once and a diquark it has twice...
Pythia8::SigmaTotal pythia_sigmatot_
An object to compute cross-sections.
double PNegB_
backward lightcone momentum p^{-} of incoming particle B in CM-frame [GeV]
std::array< FourVector, 2 > pcom_
momenta of incoming particles in the center of mass frame [GeV]
std::unique_ptr< Pythia8::Pythia > pythia_parton_
PYTHIA object used in hard string routine.
double sqrtsAB_
sqrt of Mandelstam variable s of collision [GeV]
static const ParticleType & find(PdgCode pdgcode)
Returns the ParticleType object for the given pdgcode.
void rearrange_excess(std::array< int, 5 > &nquark_total, std::array< std::array< int, 5 >, 2 > &excess_quark, std::array< std::array< int, 5 >, 2 > &excess_antiq)
Take total number of quarks and check if the system has enough constitents that need to be converted ...
int get_resonance_from_quark(int idq1, int idq2, double mass)
std::array< double, 3 > cross_sections_diffractive(int pdg_a, int pdg_b, double sqrt_s)
Interface to pythia_sigmatot_ to compute cross-sections of A+B-> different final states Schuler:1993w...
static void quarks_from_diquark(int diquark, int &q1, int &q2, int °_spin)
find two quarks from a diquark.
static bool append_intermediate_list(int pdgid, FourVector momentum, ParticleList &intermediate_particles)
append new particle from PYTHIA to a specific particle list
double stringz_b_produce_
parameter (StringZ:bLund) for the fragmentation function of other (produced) hadrons in soft non-diff...
int get_hadrontype_from_quark(int idq1, int idq2)
Determines hadron type from valence quark constituents.
void init(const ParticleList &incoming, double tcoll)
initialization feed intial particles, time of collision and gamma factor of the center of mass...
static int diquark_from_quarks(int q1, int q2)
Construct diquark from two quarks.
static void make_orthonormal_basis(ThreeVector &evec_polar, std::array< ThreeVector, 3 > &evec_basis)
compute three orthonormal basis vectors from unit vector in the longitudinal direction ...
T uniform_int(T min, T max)
Pythia8::Event event_intermediate_
event record for intermediate partonic state in the hard string routine
double kappa_tension_string_
string tension [GeV/fm]
double pow_fquark_alpha_
parameter for the quark distribution function
static void assign_all_scaling_factors(int baryon_string, ParticleList &outgoing_particles, const ThreeVector &evecLong, double suppression_factor)
Assign a cross section scaling factor to all outgoing particles.
void compute_incoming_lightcone_momenta()
compute the lightcone momenta of incoming particles where the longitudinal direction is set to be sam...
bool next_SDiff(bool is_AB_to_AX)
Single-diffractive process is based on single pomeron exchange described in Ingelman:1984ns.
double stringz_a_produce_
parameter (StringZ:aLund) for the fragmentation function of other (produced) hadrons in soft non-diff...
static void assign_scaling_factor(int nquark, ParticleData &data, double suppression_factor)
Assign a cross section scaling factor to the given particle.
void shift_rapidity_event(Pythia8::Event &event, std::array< ThreeVector, 3 > &evec_basis, double factor_yrapid, double diff_yrapid)
Shift the momentum rapidity of all particles in a given event record.
void compose_string_parton(bool find_forward_string, Pythia8::Event &event_intermediate, Pythia8::Event &event_hadronize)
Identify a set of partons, which are connected to form a color-neutral string, from a given PYTHIA ev...
bool next_DDiff()
Double-diffractive process ( A + B -> X + X ) is similar to the single-diffractive process...
bool separate_fragment_baryon_
Whether to use a separate fragmentation function for leading baryons.
void set_4momentum(const FourVector &momentum_vector)
Set the particle's 4-momentum directly.
ThreeVector vcomAB_
velocity three vector of the center of mass in the lab frame
ParticleList final_state_
final state array which must be accessed after the collision
bool pythia_parton_initialized_
Remembers if Pythia is initialized or not.
static std::pair< int, int > find_leading(int nq1, int nq2, ParticleList &list)
Find the leading string fragments.
PdgCode stores a Particle Data Group Particle Numbering Scheme particle type number.
void common_setup_pythia(Pythia8::Pythia *pythia_in, double strange_supp, double diquark_supp, double popcorn_rate, double stringz_a, double stringz_b, double string_sigma_T)
Common setup of PYTHIA objects for soft and hard string routines.
double PNegA_
backward lightcone momentum p^{-} of incoming particle A in CM-frame [GeV]
double PPosB_
forward lightcone momentum p^{+} of incoming particle B in CM-frame [GeV]
Pythia8::StringFlav pythia_stringflav_
An object for the flavor selection in string fragmentation in the case of separate fragmentation func...
static void convert_KaonLS(int &pythia_id)
convert Kaon-L or Kaon-S into K0 or Anti-K0
void find_junction_leg(bool sign_color, std::vector< int > &col, Pythia8::Event &event_intermediate, Pythia8::Event &event_hadronize)
Identify partons, which are associated with junction legs, from a given PYTHIA event record...
Pythia8::Pythia * get_ptr_pythia_parton()
Function to get the PYTHIA object for hard string routine.
void compose_string_junction(bool &find_forward_string, Pythia8::Event &event_intermediate, Pythia8::Event &event_hadronize)
Identify a set of partons and junction(s), which are connected to form a color-neutral string...
FourVector ucomAB_
velocity four vector of the center of mass in the lab frame
double time_collision_
time of collision in the computational frame [fm]
void set_tension_string(double kappa_string)
set the string tension, which is used in append_final_state.
ParticleList get_final_state()
a function to get the final state particle list which is called after the collision ...
bool set_mass_and_direction_2strings(const std::array< std::array< int, 2 >, 2 > &quarks, const std::array< FourVector, 2 > &pstr_com, std::array< double, 2 > &m_str, std::array< ThreeVector, 2 > &evec_str)
Determine string masses and directions in which strings are stretched.
double soft_t_form_
factor to be multiplied to formation times in soft strings
int fragment_off_hadron(bool from_forward, bool separate_fragment_baryon, std::array< ThreeVector, 3 > &evec_basis, double &ppos_string, double &pneg_string, double &QTrx_string_pos, double &QTrx_string_neg, double &QTry_string_pos, double &QTry_string_neg, Pythia8::FlavContainer &flav_string_pos, Pythia8::FlavContainer &flav_string_neg, std::vector< int > &pdgid_frag, std::vector< FourVector > &momentum_frag)
Fragment one hadron from a given string configuration if the string mass is above threshold (given by...
void set_pow_fquark(double alphapow, double betapow)
lightcone momentum fraction of quark is sampled according to probability distribution in non-diffrac...
bool restore_constituent(Pythia8::Event &event_intermediate, std::array< std::array< int, 5 >, 2 > &excess_quark, std::array< std::array< int, 5 >, 2 > &excess_antiq)
Take the intermediate partonic state from PYTHIA event with mapped hadrons and convert constituents i...
double pow_fgluon_beta_
parameter for the gluon distribution function
constexpr int maximum_rndm_seed_in_pythia
The maximum value of the random seed used in PYTHIA.
bool mass_dependent_formation_times_
Whether the formation time should depend on the mass of the fragment according to Andersson:1983ia eq...
StringProcess(double string_tension, double time_formation, double gluon_beta, double gluon_pmin, double quark_alpha, double quark_beta, double strange_supp, double diquark_supp, double sigma_perp, double stringz_a_leading, double stringz_b_leading, double stringz_a, double stringz_b, double string_sigma_T, double factor_t_form, bool mass_dependent_formation_times, double prob_proton_to_d_uu, bool separate_fragment_baryon, double popcorn_rate)
Constructor, initializes pythia.
std::array< ThreeVector, 3 > evecBasisAB_
Orthonormal basis vectors in the center of mass frame, where the 0th one is parallel to momentum of i...
void init_pythia_hadron_rndm()
Set PYTHIA random seeds to be desired values.
int get_index_forward(bool find_forward, int np_end, Pythia8::Event &event)
Obtain index of the most forward or backward particle in a given PYTHIA event record.
bool next_BBbarAnn()
Baryon-antibaryon annihilation process Based on what UrQMD Bass:1998ca, Bleicher:1999xi does...
void set_pow_fgluon(double betapow)
lightcone momentum fraction of gluon is sampled according to probability distribution P(x) = 1/x * (1...
double sigma_qperp_
Transverse momentum spread of the excited strings.
The FourVector class holds relevant values in Minkowski spacetime with (+, −, −, −) metric signature.
static void find_excess_constituent(PdgCode &pdg_actual, PdgCode &pdg_mapped, std::array< int, 5 > &excess_quark, std::array< int, 5 > &excess_antiq)
Compare the valence quark contents of the actual and mapped hadrons and evaluate how many more consti...
bool next_NDiffHard()
Hard Non-diffractive process is based on PYTHIA 8 with partonic showers and interactions.
double stringz_a_leading_
parameter (StringZ:aLund) for the fragmentation function of leading baryon in soft non-diffractive st...
double massB_
mass of incoming particle B [GeV]
ParticleData contains the dynamic information of a certain particle.
double time_formation_const_
constant proper time in the case of constant formation time [fm]
double massA_
mass of incoming particle A [GeV]
bool splitting_gluon_qqbar(Pythia8::Event &event_intermediate, std::array< int, 5 > &nquark_total, std::array< int, 5 > &nantiq_total, bool sign_constituent, std::array< std::array< int, 5 >, 2 > &excess_constituent)
Take total number of quarks and check if the system has enough constituents that need to be converted...
bool make_lightcone_final_two(bool separate_fragment_hadron, double ppos_string, double pneg_string, double mTrn_had_forward, double mTrn_had_backward, double &ppos_had_forward, double &ppos_had_backward, double &pneg_had_forward, double &pneg_had_backward)
Determines lightcone momenta of two final hadrons fragmented from a string in the same way as StringF...
double pow_fquark_beta_
parameter for the quark distribution function