File indexing completed on 2025-01-18 10:06:24
0001
0002
0003
0004
0005
0006
0007
0008
0009 #ifndef Pythia8_HadronLevel_H
0010 #define Pythia8_HadronLevel_H
0011
0012 #include "Pythia8/Basics.h"
0013 #include "Pythia8/BoseEinstein.h"
0014 #include "Pythia8/ColourTracing.h"
0015 #include "Pythia8/DeuteronProduction.h"
0016 #include "Pythia8/Event.h"
0017 #include "Pythia8/FragmentationFlavZpT.h"
0018 #include "Pythia8/FragmentationSystems.h"
0019 #include "Pythia8/HadronWidths.h"
0020 #include "Pythia8/HiddenValleyFragmentation.h"
0021 #include "Pythia8/Info.h"
0022 #include "Pythia8/JunctionSplitting.h"
0023 #include "Pythia8/LowEnergyProcess.h"
0024 #include "Pythia8/SigmaLowEnergy.h"
0025 #include "Pythia8/MiniStringFragmentation.h"
0026 #include "Pythia8/NucleonExcitations.h"
0027 #include "Pythia8/ParticleData.h"
0028 #include "Pythia8/ParticleDecays.h"
0029 #include "Pythia8/PartonVertex.h"
0030 #include "Pythia8/PhysicsBase.h"
0031 #include "Pythia8/PythiaStdlib.h"
0032 #include "Pythia8/RHadrons.h"
0033 #include "Pythia8/Settings.h"
0034 #include "Pythia8/StringFragmentation.h"
0035 #include "Pythia8/TimeShower.h"
0036 #include "Pythia8/UserHooks.h"
0037
0038 namespace Pythia8 {
0039
0040
0041
0042
0043
0044
0045 class HadronLevel : public PhysicsBase {
0046
0047 public:
0048
0049
0050 HadronLevel() = default;
0051
0052
0053 bool init( TimeShowerPtr timesDecPtr, RHadrons* rHadronsPtrIn,
0054 DecayHandlerPtr decayHandlePtr, vector<int> handledParticles,
0055 StringIntPtr stringInteractionsPtrIn, PartonVertexPtr partonVertexPtrIn,
0056 SigmaLowEnergy& sigmaLowEnergyIn,
0057 NucleonExcitations& nucleonExcitationsIn);
0058
0059
0060 StringFlav* getStringFlavPtr() {return &flavSel;}
0061
0062
0063 bool next(Event& event);
0064
0065
0066 bool decay( int iDec, Event& event) { return
0067 (event[iDec].isFinal() && event[iDec].canDecay() && event[iDec].mayDecay())
0068 ? decays.decay( iDec, event) : true;}
0069
0070
0071 bool moreDecays(Event& event);
0072
0073
0074 bool rescatter(Event& event);
0075
0076
0077 bool initLowEnergyProcesses();
0078 int pickLowEnergyProcess(int idA, int idB, double eCM, double mA, double mB);
0079
0080
0081 bool doLowEnergyProcess(int i1, int i2, int procTypeIn, Event& event) {
0082 if (!lowEnergyProcess.collide( i1, i2, procTypeIn, event)) {
0083 loggerPtr->ERROR_MSG("low energy collision failed");
0084 return false;
0085 }
0086 return true;
0087 }
0088
0089
0090 bool hasVetoedHadronize() const {return doHadronizeVeto; }
0091
0092 protected:
0093
0094 virtual void onInitInfoPtr() override{
0095 registerSubObject(flavSel);
0096 registerSubObject(pTSel);
0097 registerSubObject(zSel);
0098 registerSubObject(stringFrag);
0099 registerSubObject(ministringFrag);
0100 registerSubObject(decays);
0101 registerSubObject(lowEnergyProcess);
0102 registerSubObject(boseEinstein);
0103 registerSubObject(hiddenvalleyFrag);
0104 registerSubObject(junctionSplitting);
0105 registerSubObject(deuteronProd);
0106 }
0107
0108 private:
0109
0110
0111 static const double MTINY;
0112
0113
0114 bool doHadronize{}, doDecay{}, doPartonVertex{}, doBoseEinstein{},
0115 doDeuteronProd{}, allowRH{}, closePacking{}, doNonPertAll{};
0116 double mStringMin{}, pNormJunction{}, widthSepBE{}, widthSepRescatter{};
0117 vector<int> nonPertProc{};
0118
0119
0120 ColConfig colConfig;
0121
0122
0123 vector<int> iParton{}, iJunLegA{}, iJunLegB{}, iJunLegC{},
0124 iAntiLegA{}, iAntiLegB{}, iAntiLegC{}, iGluLeg{};
0125 vector<double> m2Pair{};
0126
0127
0128 StringFragmentation stringFrag;
0129
0130
0131 MiniStringFragmentation ministringFrag;
0132
0133
0134 bool tryMiniAfterFailedFrag{};
0135
0136
0137 ParticleDecays decays;
0138
0139
0140 BoseEinstein boseEinstein;
0141
0142
0143 DeuteronProduction deuteronProd;
0144
0145
0146 StringFlav flavSel;
0147 StringPT pTSel;
0148 StringZ zSel;
0149
0150
0151 ColourTracing colTrace;
0152
0153
0154 JunctionSplitting junctionSplitting;
0155
0156
0157 RHadrons* rHadronsPtr;
0158
0159
0160 HiddenValleyFragmentation hiddenvalleyFrag;
0161 bool useHiddenValley{};
0162
0163
0164 bool decayOctetOnia(Event& event);
0165
0166
0167
0168 bool findSinglets(Event& event, bool keepJunctions = false);
0169
0170
0171 PartonVertexPtr partonVertexPtr;
0172
0173
0174 class PriorityNode;
0175 bool doRescatter{}, scatterManyTimes{}, scatterQuickCheck{},
0176 scatterNeighbours{}, delayRegeneration{};
0177 double b2Max, tauRegeneration{};
0178 void queueDecResc(Event& event, int iStart,
0179 priority_queue<HadronLevel::PriorityNode>& queue);
0180 int boostDir;
0181 double boost;
0182 bool doBoost;
0183 bool useVelocityFrame;
0184
0185
0186 bool doHadronizeVeto;
0187
0188
0189 LowEnergyProcess lowEnergyProcess;
0190 int impactModel{};
0191 double impactOpacity{};
0192
0193
0194 SigmaLowEnergy* sigmaLowEnergyPtr = {};
0195
0196
0197 NucleonExcitations* nucleonExcitationsPtr = {};
0198
0199
0200 StringRepPtr stringRepulsionPtr;
0201 FragModPtr fragmentationModifierPtr;
0202
0203
0204 vector< vector< pair<double,double> > > rapidityPairs(Event& event);
0205
0206
0207 double yMax(Particle pIn, double mTiny) {
0208 double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, pIn.mT()) );
0209 return (pIn.pz() > 0) ? temp : -temp; }
0210 double yMax(Vec4 pIn, double mTiny) {
0211 double mTemp = pIn.m2Calc() + pIn.pT2();
0212 mTemp = (mTemp >= 0.) ? sqrt(mTemp) : -sqrt(-mTemp);
0213 double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, mTemp) );
0214 return (pIn.pz() > 0) ? temp : -temp; }
0215
0216 };
0217
0218
0219
0220 }
0221
0222 #endif