Back to home page

EIC code displayed by LXR

 
 

    


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

0001 // LowEnergyProcess.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 // Header file for low-energy hadronic collisions, as used for rescattering.
0007 
0008 #ifndef Pythia8_LowEnergyProcess_H
0009 #define Pythia8_LowEnergyProcess_H
0010 
0011 #include "Pythia8/Basics.h"
0012 #include "Pythia8/Event.h"
0013 #include "Pythia8/FragmentationSystems.h"
0014 #include "Pythia8/SigmaLowEnergy.h"
0015 #include "Pythia8/MiniStringFragmentation.h"
0016 #include "Pythia8/HadronWidths.h"
0017 #include "Pythia8/NucleonExcitations.h"
0018 #include "Pythia8/PythiaStdlib.h"
0019 #include "Pythia8/StringFragmentation.h"
0020 
0021 namespace Pythia8 {
0022 
0023 //==========================================================================
0024 
0025 // LowEnergyProcess class.
0026 // Is used to describe the low-energy collision between two hadrons.
0027 
0028 class LowEnergyProcess : public PhysicsBase {
0029 
0030 public:
0031 
0032   // Constructor.
0033   LowEnergyProcess() = default;
0034 
0035   // Initialize the class.
0036   void init( StringFlav* flavSelPtrIn, StringFragmentation* stringFragPtrIn,
0037     MiniStringFragmentation* ministringFragPtrIn,
0038     SigmaLowEnergy* sigmaLowEnergyPtrIn,
0039     NucleonExcitations* nucleonExcitationsPtrIn);
0040 
0041   // Produce outgoing primary hadrons from collision of incoming pair.
0042   bool collide( int i1, int i2, int typeIn, Event& event, Vec4 vtx = Vec4(),
0043     Vec4 vtx1 = Vec4(), Vec4 vtx2 = Vec4());
0044 
0045   // Event record to handle hadronization.
0046   Event         leEvent;
0047 
0048   // Give access to b slope in elastic and diffractive interactions.
0049   double bSlope( int id1In, int id2In, double eCMIn, double mAIn, double mBIn,
0050     int typeIn = 2) { id1 = id1In; id2 = id2In; eCM = eCMIn, sCM = eCM * eCM;
0051     mA = mAIn; mB = mBIn; type = typeIn; return bSlope();}
0052 
0053 private:
0054 
0055   // Initialization flag.
0056   bool isInit = false;
0057 
0058   // Parameters of the generation process.
0059   double probStoUD, fracEtass, fracEtaPss, xPowMes, xPowBar, xDiqEnhance,
0060          sigmaQ, mStringMin, sProton, probDoubleAnn;
0061 
0062   // Properties of the current collision. 1 or 2 is two incoming hadrons.
0063   // "c" or "ac" is colour or anticolour component of hadron.
0064   bool   isBaryon1, isBaryon2;
0065   int    type, sizeOld, id1, id2, idc1, idac1, idc2, idac2, nHadron,
0066          id1sv = {}, id2sv = {};
0067   double m1, m2, eCM, sCM, mThr1, mThr2, z1, z2, mT1, mT2, mA, mB,
0068          mc1, mac1, px1, py1, pTs1, mTsc1, mTsac1, mTc1, mTac1,
0069          mc2, mac2, px2, py2, pTs2, mTsc2, mTsac2, mTc2, mTac2, bA, bB;
0070 
0071   // Pointer to class for flavour generation.
0072   StringFlav* flavSelPtr;
0073 
0074   // Pointer to the generator for normal string fragmentation.
0075   StringFragmentation* stringFragPtr;
0076 
0077   // Pointer to the generator for special low-mass string fragmentation.
0078   MiniStringFragmentation* ministringFragPtr;
0079 
0080   // Separate configuration for simple collisions.
0081   ColConfig simpleColConfig;
0082 
0083   // Cross sections for low-energy processes.
0084   SigmaLowEnergy* sigmaLowEnergyPtr;
0085 
0086   // Pointer to class for handling nucleon excitations
0087   NucleonExcitations* nucleonExcitationsPtr;
0088 
0089   // Handle inelastic nondiffractive collision.
0090   bool nondiff();
0091 
0092   // Handle elastic and diffractive collisions.
0093   bool eldiff();
0094 
0095   // Handle excitation collisions.
0096   bool excitation();
0097 
0098   // Handle annihilation collisions.
0099   bool annihilation();
0100 
0101   // Handle resonant collisions.
0102   bool resonance();
0103 
0104   // Simple version of hadronization for low-energy hadronic collisions.
0105   bool simpleHadronization();
0106 
0107   // Special case with isotropic two-body final state.
0108   bool twoBody();
0109 
0110   // Special case with isotropic three-body final state.
0111   bool threeBody();
0112 
0113   // Split up hadron A or B into a colour pair, with masses and pT values.
0114   bool splitA( double mMax, double redMpT, bool splitFlavour = true);
0115   bool splitB( double mMax, double redMpT, bool splitFlavour = true);
0116 
0117   // Split a hadron inte a colour and an anticolour part.
0118   pair< int, int> splitFlav( int id);
0119 
0120   // Choose relative momentum of colour and anticolour constituents in hadron.
0121   double splitZ( int iq1, int iq2, double mRat1, double mRat2);
0122 
0123   // Estimate lowest possible mass state for flavour combination.
0124   double mThreshold( int iq1, int iq2);
0125 
0126   // Estimate lowest possible mass state for diffractive excitation.
0127   double mDiffThr( int idNow, double mNow);
0128 
0129   // Pick slope b of exp(b * t) for elastic and diffractive events.
0130   double bSlope();
0131 
0132 };
0133 
0134 //==========================================================================
0135 
0136 } // end namespace Pythia8
0137 
0138 #endif // Pythia8_LowEnergyProcess_H