PYTHIA  8.311
ProcessLevel.h
1 // ProcessLevel.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2024 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // This file contains the main class for process-level event generation.
7 // ProcessLevel: administrates the selection of "hard" process.
8 
9 #ifndef Pythia8_ProcessLevel_H
10 #define Pythia8_ProcessLevel_H
11 
12 #include "Pythia8/Basics.h"
13 #include "Pythia8/BeamParticle.h"
14 #include "Pythia8/Event.h"
15 #include "Pythia8/Info.h"
16 #include "Pythia8/ParticleData.h"
17 #include "Pythia8/PartonDistributions.h"
18 #include "Pythia8/PhysicsBase.h"
19 #include "Pythia8/ProcessContainer.h"
20 #include "Pythia8/PythiaStdlib.h"
21 #include "Pythia8/ResonanceDecays.h"
22 #include "Pythia8/Settings.h"
23 #include "Pythia8/SigmaTotal.h"
24 #include "Pythia8/SusyCouplings.h"
25 #include "Pythia8/SLHAinterface.h"
26 #include "Pythia8/StandardModel.h"
27 #include "Pythia8/UserHooks.h"
28 
29 namespace Pythia8 {
30 
31 //==========================================================================
32 
33 // The ProcessLevel class contains the top-level routines to generate
34 // the characteristic "hard" process of an event.
35 
36 class ProcessLevel : public PhysicsBase {
37 
38 public:
39 
40  // Constructor.
41  ProcessLevel() = default;
42 
43  // Destructor to delete processes in containers.
44  ~ProcessLevel();
45 
46  // Initialization.
47  bool init(bool doLHAin, SLHAinterface* slhaInterfacePtrIn,
48  vector<SigmaProcessPtr>& sigmaPtrs, vector<PhaseSpacePtr>& phaseSpacePtrs);
49 
50  // Store or replace Les Houches pointer.
51  void setLHAPtr( LHAupPtr lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn;
52  if (iLHACont >= 0) containerPtrs[iLHACont]->setLHAPtr(lhaUpPtr);}
53 
54  // Switch to new beam particle identities; for similar hadrons only.
55  void updateBeamIDs();
56 
57  // Generate the next "hard" process.
58  bool next( Event& process, int procTypeIn = 0);
59 
60  // Special case: LHA input of resonance decay only.
61  bool nextLHAdec( Event& process);
62 
63  // Accumulate and update statistics (after possible user veto).
64  void accumulate( bool doAccumulate = true);
65 
66  // Print statistics on cross sections and number of events.
67  void statistics(bool reset = false);
68 
69  // Reset statistics.
70  void resetStatistics();
71 
72  // Add any junctions to the process event record list.
73  void findJunctions( Event& junEvent);
74 
75  // Initialize and call resonance decays separately.
76  void initDecays( LHAupPtr lhaUpPtrIn) {
77  containerLHAdec.setLHAPtr(lhaUpPtrIn, particleDataPtr, settingsPtr,
78  rndmPtr); }
79 
80  bool nextDecays( Event& process) { return resonanceDecays.next( process);}
81 
82 protected:
83 
84  virtual void onInitInfoPtr() override {
85  registerSubObject(resonanceDecays);
86  registerSubObject(gammaKin);
87  }
88 
89 private:
90 
91  // Constants: could only be changed in the code itself.
92  static const int MAXLOOP;
93 
94  // Generic info for process generation.
95  bool doVarEcm, allowIDAswitch, doSecondHard, doSameCuts, allHardSame,
96  noneHardSame, someHardSame, cutsAgree, cutsOverlap, doResDecays,
97  doISR, doMPI, doWt2, switchedID, switchedEcm, useHVcols;
98  int startColTag, procType;
99  double maxPDFreweight, mHatMin1, mHatMax1, pTHatMin1, pTHatMax1, mHatMin2,
100  mHatMax2, pTHatMin2, pTHatMax2, sigmaND, eCMold;
101 
102  // Info for process generation with photon beams.
103  bool beamHasGamma;
104  int gammaMode;
105 
106  // Vector of containers of internally-generated processes.
107  vector<ProcessContainer*> containerPtrs;
108  int iContainer, iLHACont = -1;
109  double sigmaMaxSum;
110 
111  // Ditto for optional choice of a second hard process.
112  vector<ProcessContainer*> container2Ptrs;
113  int i2Container;
114  double sigma2MaxSum;
115 
116  // Single half-dummy container for LHA input of resonance decay only.
117  ProcessContainer containerLHAdec;
118 
119  // Pointer to SusyLesHouches object for interface to SUSY spectra.
120  SLHAinterface* slhaInterfacePtr;
121 
122  // Pointer to LHAup for generating external events.
123  LHAupPtr lhaUpPtr;
124 
125  // ResonanceDecay object does sequential resonance decays.
126  ResonanceDecays resonanceDecays;
127 
128  // Samples photon kinematics from leptons.
129  GammaKinematics gammaKin;
130 
131  // Generate the next event with one interaction.
132  bool nextOne( Event& process);
133 
134  // Generate the next event with two hard interactions.
135  bool nextTwo( Event& process);
136 
137  // Check that enough room for beam remnants in photon beam.
138  bool roomForRemnants();
139 
140  // Append the second to the first process list.
141  void combineProcessRecords( Event& process, Event& process2);
142 
143  // Check that colours match up.
144  bool checkColours( Event& process);
145 
146  // Print statistics when two hard processes allowed.
147  void statistics2(bool reset);
148 
149 };
150 
151 //==========================================================================
152 
153 } // end namespace Pythia8
154 
155 #endif // Pythia8_ProcessLevel_H
bool init(bool doLHAin, SLHAinterface *slhaInterfacePtrIn, vector< SigmaProcessPtr > &sigmaPtrs, vector< PhaseSpacePtr > &phaseSpacePtrs)
Initialization.
Definition: ProcessLevel.cc:44
Settings * settingsPtr
Pointer to the settings database.
Definition: PhysicsBase.h:81
Definition: PhysicsBase.h:27
void registerSubObject(PhysicsBase &pb)
Register a sub object that should have its information in sync with this.
Definition: PhysicsBase.cc:56
The Event class holds all info on the generated event.
Definition: Event.h:453
bool next(Event &process, int procTypeIn=0)
Generate the next "hard" process.
Definition: ProcessLevel.cc:383
Definition: ProcessContainer.h:39
Definition: SLHAinterface.h:27
virtual void onInitInfoPtr() override
Definition: ProcessLevel.h:84
Rndm * rndmPtr
Pointer to the random number generator.
Definition: PhysicsBase.h:93
Definition: ProcessLevel.h:36
void resetStatistics()
Reset statistics.
Definition: ProcessLevel.cc:626
void initDecays(LHAupPtr lhaUpPtrIn)
Initialize and call resonance decays separately.
Definition: ProcessLevel.h:76
~ProcessLevel()
Destructor to delete processes in containers.
Definition: ProcessLevel.cc:28
ParticleData * particleDataPtr
Pointer to the particle data table.
Definition: PhysicsBase.h:84
ProcessLevel()=default
Constructor.
void setLHAPtr(LHAupPtr lhaUpPtrIn, ParticleData *particleDataPtrIn=0, Settings *settingsPtrIn=0, Rndm *rndmPtrIn=0)
Store or replace Les Houches pointer.
Definition: ProcessContainer.h:66
void findJunctions(Event &junEvent)
Add any junctions to the process event record list.
Definition: ProcessLevel.cc:1228
void statistics(bool reset=false)
Print statistics on cross sections and number of events.
Definition: ProcessLevel.cc:510
Class to sample the virtuality and transverse momentum of emitted photons.
Definition: GammaKinematics.h:23
void setLHAPtr(LHAupPtr lhaUpPtrIn)
Store or replace Les Houches pointer.
Definition: ProcessLevel.h:51
void accumulate(bool doAccumulate=true)
Accumulate and update statistics (after possible user veto).
Definition: ProcessLevel.cc:422
void updateBeamIDs()
Switch to new beam particle identities; for similar hadrons only.
Definition: ProcessLevel.cc:366
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
Definition: ResonanceDecays.h:28
bool next(Event &process, int iDecNow=0)
Generate the next decay sequence.
Definition: ResonanceDecays.cc:48
bool nextLHAdec(Event &process)
Special case: LHA input of resonance decay only.
Definition: ProcessLevel.cc:402