9 #ifndef Pythia8_HadronLevel_H
10 #define Pythia8_HadronLevel_H
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/LowEnergySigma.h"
25 #include "Pythia8/MiniStringFragmentation.h"
26 #include "Pythia8/NucleonExcitations.h"
27 #include "Pythia8/ParticleData.h"
28 #include "Pythia8/ParticleDecays.h"
29 #include "Pythia8/PartonVertex.h"
30 #include "Pythia8/PhysicsBase.h"
31 #include "Pythia8/PythiaStdlib.h"
32 #include "Pythia8/RHadrons.h"
33 #include "Pythia8/Settings.h"
34 #include "Pythia8/StringFragmentation.h"
35 #include "Pythia8/TimeShower.h"
36 #include "Pythia8/UserHooks.h"
45 class HadronLevel :
public PhysicsBase {
50 HadronLevel() =
default;
53 bool init( TimeShowerPtr timesDecPtr, RHadrons* rHadronsPtrIn,
54 DecayHandlerPtr decayHandlePtr, vector<int> handledParticles,
55 StringIntPtr stringInteractionsPtrIn, PartonVertexPtr partonVertexPtrIn);
58 StringFlav* getStringFlavPtr() {
return &flavSel;}
61 bool next(
Event& event);
64 bool moreDecays(
Event& event);
67 bool initLowEnergyProcesses();
68 int pickLowEnergyProcess(
int idA,
int idB,
double eCM,
double mA,
double mB);
71 bool doLowEnergyProcess(
int i1,
int i2,
int type,
Event& event) {
72 if (!lowEnergyProcess.collide( i1, i2, type, event)) {
73 infoPtr->errorMsg(
"Error in HadronLevel::doLowEnergyProcess: "
74 "Low energy collision failed");
81 double getLowEnergySigma(
int idA,
int idB,
double eCM,
double mA,
82 double mB,
int type = 0) {
83 return lowEnergySigma.sigmaPartial( idA, idB, eCM, mA, mB, type); }
86 double getLowEnergySlope(
int idA,
int idB,
double eCM,
double mA,
87 double mB,
int type = 2) {
88 return lowEnergyProcess.bSlope( idA, idB, eCM, mA, mB, type); }
92 virtual void onInitInfoPtr()
override{
93 registerSubObject(flavSel);
94 registerSubObject(pTSel);
95 registerSubObject(zSel);
96 registerSubObject(stringFrag);
97 registerSubObject(ministringFrag);
98 registerSubObject(decays);
99 registerSubObject(lowEnergyProcess);
100 registerSubObject(lowEnergySigma);
101 registerSubObject(nucleonExcitations);
102 registerSubObject(boseEinstein);
103 registerSubObject(hiddenvalleyFrag);
104 registerSubObject(junctionSplitting);
105 registerSubObject(deuteronProd);
111 static const double MTINY;
114 bool doHadronize{}, doDecay{}, doPartonVertex{}, doBoseEinstein{},
115 doDeuteronProd{}, allowRH{}, closePacking{}, doNonPertAll{};
116 double mStringMin{}, eNormJunction{}, widthSepBE{}, widthSepRescatter{};
117 vector<int> nonPertProc{};
123 vector<int> iParton{}, iJunLegA{}, iJunLegB{}, iJunLegC{},
124 iAntiLegA{}, iAntiLegB{}, iAntiLegC{}, iGluLeg{};
125 vector<double> m2Pair{};
128 StringFragmentation stringFrag;
131 MiniStringFragmentation ministringFrag;
134 ParticleDecays decays;
137 BoseEinstein boseEinstein;
140 DeuteronProduction deuteronProd;
148 ColourTracing colTrace;
151 JunctionSplitting junctionSplitting;
154 RHadrons* rHadronsPtr;
157 HiddenValleyFragmentation hiddenvalleyFrag;
158 bool useHiddenValley{};
161 bool decayOctetOnia(
Event& event);
165 bool findSinglets(
Event& event,
bool keepJunctions =
false);
168 PartonVertexPtr partonVertexPtr;
172 bool doRescatter{}, scatterManyTimes{}, scatterQuickCheck{},
173 scatterNeighbours{}, delayRegeneration{};
174 double b2Max, tauRegeneration{};
175 void queueDecResc(
Event& event,
int iStart,
176 priority_queue<HadronLevel::PriorityNode>& queue);
180 bool useVelocityFrame;
183 LowEnergyProcess lowEnergyProcess;
185 double impactOpacity{};
188 LowEnergySigma lowEnergySigma;
191 NucleonExcitations nucleonExcitations = {};
194 StringRepPtr stringRepulsionPtr;
195 FragModPtr fragmentationModifierPtr;
198 vector< vector< pair<double,double> > > rapidityPairs(
Event& event);
201 double yMax(Particle pIn,
double mTiny) {
202 double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, pIn.mT()) );
203 return (pIn.pz() > 0) ? temp : -temp; }
211 #endif // Pythia8_HadronLevel_H