Back to home page

EIC code displayed by LXR

 
 

    


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

0001 // ProcessLevel.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 process-level event generation.
0007 // ProcessLevel: administrates the selection of "hard" process.
0008 
0009 #ifndef Pythia8_ProcessLevel_H
0010 #define Pythia8_ProcessLevel_H
0011 
0012 #include "Pythia8/Basics.h"
0013 #include "Pythia8/BeamParticle.h"
0014 #include "Pythia8/Event.h"
0015 #include "Pythia8/Info.h"
0016 #include "Pythia8/ParticleData.h"
0017 #include "Pythia8/PartonDistributions.h"
0018 #include "Pythia8/PhysicsBase.h"
0019 #include "Pythia8/ProcessContainer.h"
0020 #include "Pythia8/PythiaStdlib.h"
0021 #include "Pythia8/ResonanceDecays.h"
0022 #include "Pythia8/Settings.h"
0023 #include "Pythia8/SigmaTotal.h"
0024 #include "Pythia8/SusyCouplings.h"
0025 #include "Pythia8/SLHAinterface.h"
0026 #include "Pythia8/StandardModel.h"
0027 #include "Pythia8/UserHooks.h"
0028 
0029 namespace Pythia8 {
0030 
0031 //==========================================================================
0032 
0033 // The ProcessLevel class contains the top-level routines to generate
0034 // the characteristic "hard" process of an event.
0035 
0036 class ProcessLevel : public PhysicsBase {
0037 
0038 public:
0039 
0040   // Constructor.
0041   ProcessLevel() = default;
0042 
0043   // Destructor to delete processes in containers.
0044   ~ProcessLevel();
0045 
0046   // Initialization.
0047   bool init(bool doLHAin, SLHAinterface* slhaInterfacePtrIn,
0048     vector<SigmaProcessPtr>& sigmaPtrs, vector<PhaseSpacePtr>& phaseSpacePtrs);
0049 
0050   // Store or replace Les Houches pointer.
0051   void setLHAPtr( LHAupPtr lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn;
0052     if (iLHACont >= 0) containerPtrs[iLHACont]->setLHAPtr(lhaUpPtr);}
0053 
0054   // Switch to new beam particle identities; for similar hadrons only.
0055   void updateBeamIDs();
0056 
0057   // Generate the next "hard" process.
0058   bool next( Event& process, int procTypeIn = 0);
0059 
0060   // Special case: LHA input of resonance decay only.
0061   bool nextLHAdec( Event& process);
0062 
0063   // Accumulate and update statistics (after possible user veto).
0064   void accumulate( bool doAccumulate = true);
0065 
0066   // Print statistics on cross sections and number of events.
0067   void statistics(bool reset = false);
0068 
0069   // Reset statistics.
0070   void resetStatistics();
0071 
0072   // Add any junctions to the process event record list.
0073   void findJunctions( Event& junEvent);
0074 
0075   // Initialize and call resonance decays separately.
0076   void initDecays( LHAupPtr lhaUpPtrIn) {
0077     containerLHAdec.setLHAPtr(lhaUpPtrIn, particleDataPtr, settingsPtr,
0078       rndmPtr); }
0079 
0080   bool nextDecays( Event& process) { return resonanceDecays.next( process);}
0081 
0082 protected:
0083 
0084   virtual void onInitInfoPtr() override {
0085     registerSubObject(resonanceDecays);
0086     registerSubObject(gammaKin);
0087   }
0088 
0089 private:
0090 
0091   // Constants: could only be changed in the code itself.
0092   static const int MAXLOOP;
0093 
0094   // Generic info for process generation.
0095   bool   doVarEcm, allowIDAswitch, doSecondHard, doSameCuts, allHardSame,
0096          noneHardSame, someHardSame, cutsAgree, cutsOverlap, doResDecays,
0097          doISR, doMPI, doWt2, switchedID, switchedEcm, useHVcols;
0098   int    startColTag, procType;
0099   double maxPDFreweight, mHatMin1, mHatMax1, pTHatMin1, pTHatMax1, mHatMin2,
0100          mHatMax2, pTHatMin2, pTHatMax2, sigmaND, eCMold;
0101 
0102   // Info for process generation with photon beams.
0103   bool   beamHasGamma;
0104   int    gammaMode;
0105 
0106   // Vector of containers of internally-generated processes.
0107   vector<ProcessContainer*> containerPtrs;
0108   int    iContainer, iLHACont = -1;
0109   double sigmaMaxSum;
0110 
0111   // Ditto for optional choice of a second hard process.
0112   vector<ProcessContainer*> container2Ptrs;
0113   int    i2Container;
0114   double sigma2MaxSum;
0115 
0116   // Single half-dummy container for LHA input of resonance decay only.
0117   ProcessContainer containerLHAdec;
0118 
0119   // Pointer to SusyLesHouches object for interface to SUSY spectra.
0120   SLHAinterface*  slhaInterfacePtr;
0121 
0122   // Pointer to LHAup for generating external events.
0123   LHAupPtr          lhaUpPtr;
0124 
0125   // ResonanceDecay object does sequential resonance decays.
0126   ResonanceDecays resonanceDecays;
0127 
0128   // Samples photon kinematics from leptons.
0129   GammaKinematics gammaKin;
0130 
0131   // Generate the next event with one interaction.
0132   bool nextOne( Event& process);
0133 
0134   // Generate the next event with two hard interactions.
0135   bool nextTwo( Event& process);
0136 
0137   // Check that enough room for beam remnants in photon beam.
0138   bool roomForRemnants();
0139 
0140   // Append the second to the first process list.
0141   void combineProcessRecords( Event& process, Event& process2);
0142 
0143   // Check that colours match up.
0144   bool checkColours( Event& process);
0145 
0146   // Print statistics when two hard processes allowed.
0147   void statistics2(bool reset);
0148 
0149 };
0150 
0151 //==========================================================================
0152 
0153 } // end namespace Pythia8
0154 
0155 #endif // Pythia8_ProcessLevel_H