PYTHIA  8.317
HadronLevel.h
1 // HadronLevel.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2026 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 // This file contains the main class for hadron-level generation.
7 // HadronLevel: handles administration of fragmentation and decay.
8 
9 #ifndef Pythia8_HadronLevel_H
10 #define Pythia8_HadronLevel_H
11 
12 #include "Pythia8/Basics.h"
13 #include "Pythia8/BoseEinstein.h"
14 #include "Pythia8/ColourTracing.h"
15 #include "Pythia8/DeuteronProduction.h"
16 #include "Pythia8/Event.h"
17 #include "Pythia8/FragmentationFlavZpT.h"
18 #include "Pythia8/FragmentationSystems.h"
19 #include "Pythia8/HadronWidths.h"
20 #include "Pythia8/HiddenValleyFragmentation.h"
21 #include "Pythia8/Info.h"
22 #include "Pythia8/JunctionSplitting.h"
23 #include "Pythia8/LowEnergyProcess.h"
24 #include "Pythia8/SigmaLowEnergy.h"
25 #include "Pythia8/NucleonExcitations.h"
26 #include "Pythia8/ParticleData.h"
27 #include "Pythia8/ParticleDecays.h"
28 #include "Pythia8/PartonVertex.h"
29 #include "Pythia8/PhysicsBase.h"
30 #include "Pythia8/PythiaStdlib.h"
31 #include "Pythia8/RHadrons.h"
32 #include "Pythia8/Settings.h"
33 #include "Pythia8/StringFragmentation.h"
34 #include "Pythia8/TimeShower.h"
35 #include "Pythia8/UserHooks.h"
36 
37 namespace Pythia8 {
38 
39 //==========================================================================
40 
41 // The HadronLevel class contains the top-level routines to generate
42 // the transition from the partonic to the hadronic stage of an event.
43 
44 class HadronLevel : public PhysicsBase {
45 
46 public:
47 
48  // Constructor.
49  HadronLevel() = default;
50 
51  // Initialize HadronLevel classes as required.
52  bool init( TimeShowerPtr timesDecPtrIn, RHadronsPtr rHadronsPtrIn,
53  LundFragmentationPtr fragPtrIn, vector<FragmentationModelPtr>* fragPtrsIn,
54  DecayHandlerPtr decayHandlePtr, vector<int> handledParticles,
55  StringIntPtr stringInteractionsPtrIn, PartonVertexPtr partonVertexPtrIn,
56  SigmaLowEnergy& sigmaLowEnergyIn,
57  NucleonExcitations& nucleonExcitationsIn);
58 
59  // Get pointer to StringFlav instance (needed by BeamParticle).
60  StringFlav* getStringFlavPtr() {return &flavSel;}
61 
62  // Generate the next event.
63  bool next(Event& event);
64 
65  // Try to decay the specified particle. Returns false if decay
66  // failed. If allowPartons is true the decay products may consist
67  // of partons, which need to be hadronised.
68  bool decay( int iDec, Event& event, bool allowPartons = true) { return
69  (event[iDec].isFinal() && event[iDec].canDecay() && event[iDec].mayDecay())
70  ? decays.decay( iDec, event, allowPartons) : true;}
71 
72  // Special routine to allow more decays if on/off switches changed.
73  bool moreDecays(Event& event, bool allowPartons = false);
74 
75  // Perform rescattering. Return true if new strings must be hadronized.
76  bool rescatter(Event& event);
77 
78  // Prepare and pick process for a low-energy hadron-hadron scattering.
80  int pickLowEnergyProcess(int idA, int idB, double eCM, double mA, double mB);
81 
82  // Special routine to do a low-energy hadron-hadron scattering.
83  bool doLowEnergyProcess(int i1, int i2, int procTypeIn, Event& event) {
84  if (!lowEnergyProcess.collide( i1, i2, procTypeIn, event)) {
85  loggerPtr->ERROR_MSG("low energy collision failed");
86  return false;
87  }
88  return true;
89  }
90 
91  // Tell if we did an early user-defined veto of the event.
92  bool hasVetoedHadronize() const {return doHadronizeVeto; }
93 
94 protected:
95 
96  virtual void onInitInfoPtr() override{
97  registerSubObject(flavSel);
98  registerSubObject(pTSel);
99  registerSubObject(zSel);
100  registerSubObject(decays);
101  registerSubObject(lowEnergyProcess);
102  registerSubObject(boseEinstein);
103  registerSubObject(junctionSplitting);
104  registerSubObject(deuteronProd);
105  }
106 
107 private:
108 
109  // Constants: could only be changed in the code itself.
110  static const double MTINY;
111 
112  // Initialization data, read from Settings.
113  bool doHadronize{}, doDecay{}, doPartonVertex{}, doBoseEinstein{},
114  doDeuteronProd{}, allowRH{}, closePacking{}, doNonPertAll{},
115  doQED;
116  double pNormJunction{}, widthSepBE{}, widthSepRescatter{};
117  vector<int> nonPertProc{};
118 
119  // Configuration of colour-singlet systems.
120  ColConfig colConfig;
121 
122  // Colour and mass information.
123  vector<int> iParton{}, iJunLegA{}, iJunLegB{}, iJunLegC{},
124  iAntiLegA{}, iAntiLegB{}, iAntiLegC{}, iGluLeg{};
125  vector<double> m2Pair{};
126 
127  // The generator class for normal decays.
128  ParticleDecays decays;
129 
130  // Pointer to TimeShower for interleaved QED in hadron decays.
131  TimeShowerPtr timesDecPtr;
132 
133  // The generator class for Bose-Einstein effects.
134  BoseEinstein boseEinstein;
135 
136  // The generator class for deuteron production.
137  DeuteronProduction deuteronProd;
138 
139  // Classes for flavour, pT and z generation.
140  StringFlav flavSel;
141  StringPT pTSel;
142  StringZ zSel;
143 
144  // Class for colour tracing.
145  ColourTracing colTrace;
146 
147  // Junction splitting class.
148  JunctionSplitting junctionSplitting;
149 
150  // The RHadrons class is used to fragment off and decay R-hadrons.
151  RHadronsPtr rHadronsPtr{};
152 
153  // The LundFragmentation class is used for fragmentation and low
154  // energy processes.
155  LundFragmentationPtr fragPtr{};
156 
157  // The fragmentation model pointer vector allows multiple
158  // fragmentation models to be used sequentially.
159  vector<FragmentationModelPtr>* fragPtrs{};
160 
161  // Special case: colour-octet onium decays, to be done initially.
162  bool decayOctetOnia(Event& event);
163 
164  // Trace colour flow in the event to form colour singlet subsystems.
165  // Option to keep junctions, needed for rope hadronization.
166  bool findSinglets(Event& event, bool keepJunctions = false);
167 
168  // Class to displace hadron vertices from parton impact-parameter picture.
169  PartonVertexPtr partonVertexPtr;
170 
171  // Hadronic rescattering.
172  class PriorityNode;
173  bool doRescatter{}, scatterManyTimes{}, scatterQuickCheck{},
174  scatterNeighbours{}, delayRegeneration{};
175  double b2Max, tauRegeneration{};
176  void queueDecResc(Event& event, int iStart,
177  priority_queue<HadronLevel::PriorityNode>& queue);
178  int boostDir;
179  double boost;
180  bool doBoost;
181  bool useVelocityFrame;
182 
183  // User veto performed right after string hadronization.
184  bool doHadronizeVeto;
185 
186  // The generator class for low-energy hadron-hadron collisions.
187  LowEnergyProcess lowEnergyProcess;
188  int impactModel{};
189  double impactOpacity{};
190 
191  // Cross sections for low-energy processes.
192  SigmaLowEnergy* sigmaLowEnergyPtr = {};
193 
194  // Nucleon excitations data.
195  NucleonExcitations* nucleonExcitationsPtr = {};
196 
197  // Class for event geometry for Rope Hadronization. Production vertices.
198  StringRepPtr stringRepulsionPtr;
199  FragModPtr fragmentationModifierPtr;
200 
201  // Extract rapidity pairs.
202  vector< vector< pair<double,double> > > rapidityPairs(Event& event);
203 
204  // Calculate the rapidity for string ends, protected against too large y.
205  double yMax(Particle pIn, double mTiny) {
206  double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, pIn.mT()) );
207  return (pIn.pz() > 0) ? temp : -temp; }
208  double yMax(Vec4 pIn, double mTiny) {
209  double mTemp = pIn.m2Calc() + pIn.pT2();
210  mTemp = (mTemp >= 0.) ? sqrt(mTemp) : -sqrt(-mTemp);
211  double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, mTemp) );
212  return (pIn.pz() > 0) ? temp : -temp; }
213 
214  // Fragmentation weights container.
215  WeightsFragmentation* wgtsPtr{};
216 
217 };
218 
219 //==========================================================================
220 
221 } // end namespace Pythia8
222 
223 #endif // Pythia8_HadronLevel_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
Definition: NucleonExcitations.h:23
Definition: PhysicsBase.h:26
void registerSubObject(PhysicsBase &pb)
Register a sub object that should have its information in sync with this.
Definition: PhysicsBase.cc:57
The ParticleDecays class contains the routines to decay a particle.
Definition: ParticleDecays.h:57
The Event class holds all info on the generated event.
Definition: Event.h:408
bool moreDecays(Event &event, bool allowPartons=false)
Special routine to allow more decays if on/off switches changed.
Definition: HadronLevel.cc:341
Definition: HadronLevel.cc:20
Definition: JunctionSplitting.h:30
The StringPT class is used to select select transverse momenta.
Definition: FragmentationFlavZpT.h:300
Gets cross sections for hadron-hadron collisions at low energies.
Definition: SigmaLowEnergy.h:22
StringFlav * getStringFlavPtr()
Get pointer to StringFlav instance (needed by BeamParticle).
Definition: HadronLevel.h:60
int pickLowEnergyProcess(int idA, int idB, double eCM, double mA, double mB)
Pick process for a low-energy hadron-hadron scattering.
Definition: HadronLevel.cc:391
bool decay(int iDec, Event &event, bool allowPartons=true)
Definition: HadronLevel.h:68
HadronLevel()=default
Constructor.
The StringZ class is used to sample the fragmentation function f(z).
Definition: FragmentationFlavZpT.h:219
Definition: HadronLevel.h:44
bool next(Event &event)
Generate the next event.
Definition: HadronLevel.cc:175
bool rescatter(Event &event)
Perform rescattering. Return true if new strings must be hadronized.
Definition: HadronLevel.cc:624
bool doLowEnergyProcess(int i1, int i2, int procTypeIn, Event &event)
Special routine to do a low-energy hadron-hadron scattering.
Definition: HadronLevel.h:83
bool hasVetoedHadronize() const
Tell if we did an early user-defined veto of the event.
Definition: HadronLevel.h:92
ColourTracing class. It is used to trace colours within the event record.
Definition: ColourTracing.h:22
bool initLowEnergyProcesses()
Prepare and pick process for a low-energy hadron-hadron scattering.
Definition: HadronLevel.cc:361
Definition: Event.h:32
The ColConfig class describes the colour configuration of the whole event.
Definition: FragmentationSystems.h:60
Definition: BoseEinstein.h:48
bool init(TimeShowerPtr timesDecPtrIn, RHadronsPtr rHadronsPtrIn, LundFragmentationPtr fragPtrIn, vector< FragmentationModelPtr > *fragPtrsIn, DecayHandlerPtr decayHandlePtr, vector< int > handledParticles, StringIntPtr stringInteractionsPtrIn, PartonVertexPtr partonVertexPtrIn, SigmaLowEnergy &sigmaLowEnergyIn, NucleonExcitations &nucleonExcitationsIn)
Initialize HadronLevel classes as required.
Definition: HadronLevel.cc:57
Logger * loggerPtr
Pointer to logger.
Definition: PhysicsBase.h:91
virtual void onInitInfoPtr() override
Definition: HadronLevel.h:96
Definition: Weights.h:355
bool decay(int iDec, Event &event, bool allowPartons=true)
Definition: ParticleDecays.cc:120
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:76
The DeuteronProduction class.
Definition: DeuteronProduction.h:22
Definition: Basics.h:32