Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:06:24

0001 // HadronLevel.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2024 Torbjorn Sjostrand.
0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
0004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
0005 
0006 // This file contains the main class for hadron-level generation.
0007 // HadronLevel: handles administration of fragmentation and decay.
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 // The HadronLevel class contains the top-level routines to generate
0043 // the transition from the partonic to the hadronic stage of an event.
0044 
0045 class HadronLevel : public PhysicsBase {
0046 
0047 public:
0048 
0049   // Constructor.
0050   HadronLevel() = default;
0051 
0052   // Initialize HadronLevel classes as required.
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   // Get pointer to StringFlav instance (needed by BeamParticle).
0060   StringFlav* getStringFlavPtr() {return &flavSel;}
0061 
0062   // Generate the next event.
0063   bool next(Event& event);
0064 
0065   // Try to decay the specified particle. Returns false if decay failed.
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   // Special routine to allow more decays if on/off switches changed.
0071   bool moreDecays(Event& event);
0072 
0073   // Perform rescattering. Return true if new strings must be hadronized.
0074   bool rescatter(Event& event);
0075 
0076   // Prepare and pick process for a low-energy hadron-hadron scattering.
0077   bool initLowEnergyProcesses();
0078   int pickLowEnergyProcess(int idA, int idB, double eCM, double mA, double mB);
0079 
0080   // Special routine to do a low-energy hadron-hadron scattering.
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   // Tell if we did an early user-defined veto of the event.
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   // Constants: could only be changed in the code itself.
0111   static const double MTINY;
0112 
0113   // Initialization data, read from Settings.
0114   bool doHadronize{}, doDecay{}, doPartonVertex{}, doBoseEinstein{},
0115     doDeuteronProd{}, allowRH{}, closePacking{}, doNonPertAll{};
0116   double mStringMin{}, pNormJunction{}, widthSepBE{}, widthSepRescatter{};
0117   vector<int> nonPertProc{};
0118 
0119   // Configuration of colour-singlet systems.
0120   ColConfig      colConfig;
0121 
0122   // Colour and mass information.
0123   vector<int>    iParton{}, iJunLegA{}, iJunLegB{}, iJunLegC{},
0124                  iAntiLegA{}, iAntiLegB{}, iAntiLegC{}, iGluLeg{};
0125   vector<double> m2Pair{};
0126 
0127   // The generator class for normal string fragmentation.
0128   StringFragmentation stringFrag;
0129 
0130   // The generator class for special low-mass string fragmentation.
0131   MiniStringFragmentation ministringFrag;
0132 
0133   // Try ministring fragmentation also if normal fails.
0134   bool tryMiniAfterFailedFrag{};
0135 
0136   // The generator class for normal decays.
0137   ParticleDecays decays;
0138 
0139   // The generator class for Bose-Einstein effects.
0140   BoseEinstein boseEinstein;
0141 
0142   // The generator class for deuteron production.
0143   DeuteronProduction deuteronProd;
0144 
0145   // Classes for flavour, pT and z generation.
0146   StringFlav flavSel;
0147   StringPT   pTSel;
0148   StringZ    zSel;
0149 
0150   // Class for colour tracing.
0151   ColourTracing colTrace;
0152 
0153   // Junction splitting class.
0154   JunctionSplitting junctionSplitting;
0155 
0156   // The RHadrons class is used to fragment off and decay R-hadrons.
0157   RHadrons*  rHadronsPtr;
0158 
0159   // Special class for Hidden-Valley hadronization. Not always used.
0160   HiddenValleyFragmentation hiddenvalleyFrag;
0161   bool useHiddenValley{};
0162 
0163   // Special case: colour-octet onium decays, to be done initially.
0164   bool decayOctetOnia(Event& event);
0165 
0166   // Trace colour flow in the event to form colour singlet subsystems.
0167   // Option to keep junctions, needed for rope hadronization.
0168   bool findSinglets(Event& event, bool keepJunctions = false);
0169 
0170   // Class to displace hadron vertices from parton impact-parameter picture.
0171   PartonVertexPtr partonVertexPtr;
0172 
0173   // Hadronic rescattering.
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   // User veto performed right after string hadronization.
0186   bool doHadronizeVeto;
0187 
0188   // The generator class for low-energy hadron-hadron collisions.
0189   LowEnergyProcess lowEnergyProcess;
0190   int    impactModel{};
0191   double impactOpacity{};
0192 
0193   // Cross sections for low-energy processes.
0194   SigmaLowEnergy* sigmaLowEnergyPtr = {};
0195 
0196   // Nucleon excitations data.
0197   NucleonExcitations* nucleonExcitationsPtr = {};
0198 
0199   // Class for event geometry for Rope Hadronization. Production vertices.
0200   StringRepPtr stringRepulsionPtr;
0201   FragModPtr fragmentationModifierPtr;
0202 
0203   // Extract rapidity pairs.
0204   vector< vector< pair<double,double> > > rapidityPairs(Event& event);
0205 
0206   // Calculate the rapidity for string ends, protected against too large y.
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 } // end namespace Pythia8
0221 
0222 #endif // Pythia8_HadronLevel_H