PYTHIA  8.312
LowEnergyProcess.h
1 // LowEnergyProcess.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2024 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Header file for low-energy hadronic collisions, as used for rescattering.
7 
8 #ifndef Pythia8_LowEnergyProcess_H
9 #define Pythia8_LowEnergyProcess_H
10 
11 #include "Pythia8/Basics.h"
12 #include "Pythia8/Event.h"
13 #include "Pythia8/FragmentationSystems.h"
14 #include "Pythia8/SigmaLowEnergy.h"
15 #include "Pythia8/MiniStringFragmentation.h"
16 #include "Pythia8/HadronWidths.h"
17 #include "Pythia8/NucleonExcitations.h"
18 #include "Pythia8/PythiaStdlib.h"
19 #include "Pythia8/StringFragmentation.h"
20 
21 namespace Pythia8 {
22 
23 //==========================================================================
24 
25 // LowEnergyProcess class.
26 // Is used to describe the low-energy collision between two hadrons.
27 
28 class LowEnergyProcess : public PhysicsBase {
29 
30 public:
31 
32  // Constructor.
33  LowEnergyProcess() = default;
34 
35  // Initialize the class.
36  void init( StringFlav* flavSelPtrIn, StringFragmentation* stringFragPtrIn,
37  MiniStringFragmentation* ministringFragPtrIn,
38  SigmaLowEnergy* sigmaLowEnergyPtrIn,
39  NucleonExcitations* nucleonExcitationsPtrIn);
40 
41  // Produce outgoing primary hadrons from collision of incoming pair.
42  bool collide( int i1, int i2, int typeIn, Event& event, Vec4 vtx = Vec4(),
43  Vec4 vtx1 = Vec4(), Vec4 vtx2 = Vec4());
44 
45  // Event record to handle hadronization.
47 
48  // Give access to b slope in elastic and diffractive interactions.
49  double bSlope( int id1In, int id2In, double eCMIn, double mAIn, double mBIn,
50  int typeIn = 2) { id1 = id1In; id2 = id2In; eCM = eCMIn, sCM = eCM * eCM;
51  mA = mAIn; mB = mBIn; type = typeIn; return bSlope();}
52 
53 private:
54 
55  // Initialization flag.
56  bool isInit = false;
57 
58  // Parameters of the generation process.
59  double probStoUD, fracEtass, fracEtaPss, xPowMes, xPowBar, xDiqEnhance,
60  sigmaQ, mStringMin, sProton, probDoubleAnn;
61 
62  // Properties of the current collision. 1 or 2 is two incoming hadrons.
63  // "c" or "ac" is colour or anticolour component of hadron.
64  bool isBaryon1, isBaryon2;
65  int type, sizeOld, id1, id2, idc1, idac1, idc2, idac2, nHadron,
66  id1sv = {}, id2sv = {};
67  double m1, m2, eCM, sCM, mThr1, mThr2, z1, z2, mT1, mT2, mA, mB,
68  mc1, mac1, px1, py1, pTs1, mTsc1, mTsac1, mTc1, mTac1,
69  mc2, mac2, px2, py2, pTs2, mTsc2, mTsac2, mTc2, mTac2, bA, bB;
70 
71  // Pointer to class for flavour generation.
72  StringFlav* flavSelPtr;
73 
74  // Pointer to the generator for normal string fragmentation.
75  StringFragmentation* stringFragPtr;
76 
77  // Pointer to the generator for special low-mass string fragmentation.
78  MiniStringFragmentation* ministringFragPtr;
79 
80  // Separate configuration for simple collisions.
81  ColConfig simpleColConfig;
82 
83  // Cross sections for low-energy processes.
84  SigmaLowEnergy* sigmaLowEnergyPtr;
85 
86  // Pointer to class for handling nucleon excitations
87  NucleonExcitations* nucleonExcitationsPtr;
88 
89  // Handle inelastic nondiffractive collision.
90  bool nondiff();
91 
92  // Handle elastic and diffractive collisions.
93  bool eldiff();
94 
95  // Handle excitation collisions.
96  bool excitation();
97 
98  // Handle annihilation collisions.
99  bool annihilation();
100 
101  // Handle resonant collisions.
102  bool resonance();
103 
104  // Simple version of hadronization for low-energy hadronic collisions.
105  bool simpleHadronization();
106 
107  // Special case with isotropic two-body final state.
108  bool twoBody();
109 
110  // Special case with isotropic three-body final state.
111  bool threeBody();
112 
113  // Split up hadron A or B into a colour pair, with masses and pT values.
114  bool splitA( double mMax, double redMpT, bool splitFlavour = true);
115  bool splitB( double mMax, double redMpT, bool splitFlavour = true);
116 
117  // Split a hadron inte a colour and an anticolour part.
118  pair< int, int> splitFlav( int id);
119 
120  // Choose relative momentum of colour and anticolour constituents in hadron.
121  double splitZ( int iq1, int iq2, double mRat1, double mRat2);
122 
123  // Estimate lowest possible mass state for flavour combination.
124  double mThreshold( int iq1, int iq2);
125 
126  // Estimate lowest possible mass state for diffractive excitation.
127  double mDiffThr( int idNow, double mNow);
128 
129  // Pick slope b of exp(b * t) for elastic and diffractive events.
130  double bSlope();
131 
132 };
133 
134 //==========================================================================
135 
136 } // end namespace Pythia8
137 
138 #endif // Pythia8_LowEnergyProcess_H
Definition: LowEnergyProcess.h:28
bool collide(int i1, int i2, int typeIn, Event &event, Vec4 vtx=Vec4(), Vec4 vtx1=Vec4(), Vec4 vtx2=Vec4())
Produce outgoing primary hadrons from collision of incoming pair.
Definition: LowEnergyProcess.cc:105
double bSlope(int id1In, int id2In, double eCMIn, double mAIn, double mBIn, int typeIn=2)
Give access to b slope in elastic and diffractive interactions.
Definition: LowEnergyProcess.h:49
Definition: NucleonExcitations.h:23
Definition: PhysicsBase.h:27
The Event class holds all info on the generated event.
Definition: Event.h:453
Gets cross sections for hadron-hadron collisions at low energies.
Definition: SigmaLowEnergy.h:22
Definition: StringFragmentation.h:119
LowEnergyProcess()=default
Constructor.
The ColConfig class describes the colour configuration of the whole event.
Definition: FragmentationSystems.h:60
Event leEvent
Event record to handle hadronization.
Definition: LowEnergyProcess.h:46
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:84
Definition: MiniStringFragmentation.h:30
Definition: Basics.h:32
void init(StringFlav *flavSelPtrIn, StringFragmentation *stringFragPtrIn, MiniStringFragmentation *ministringFragPtrIn, SigmaLowEnergy *sigmaLowEnergyPtrIn, NucleonExcitations *nucleonExcitationsPtrIn)
Initialize the class.
Definition: LowEnergyProcess.cc:50