PYTHIA  8.317
ProcessContainer.h
1 // ProcessContainer.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2026 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 collected machinery of a process.
7 // ProcessContainer: contains information on a particular process.
8 // SetupContainers: administrates the selection/creation of processes.
9 
10 #ifndef Pythia8_ProcessContainer_H
11 #define Pythia8_ProcessContainer_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/BeamParticle.h"
15 #include "Pythia8/Event.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/ParticleData.h"
18 #include "Pythia8/PartonDistributions.h"
19 #include "Pythia8/PhaseSpace.h"
20 #include "Pythia8/PhysicsBase.h"
21 #include "Pythia8/PythiaStdlib.h"
22 #include "Pythia8/ResonanceDecays.h"
23 #include "Pythia8/Settings.h"
24 #include "Pythia8/SigmaProcess.h"
25 #include "Pythia8/SigmaTotal.h"
26 #include "Pythia8/SigmaOnia.h"
27 #include "Pythia8/StandardModel.h"
28 #include "Pythia8/SusyCouplings.h"
29 #include "Pythia8/SLHAinterface.h"
30 #include "Pythia8/UserHooks.h"
31 
32 namespace Pythia8 {
33 
34 //==========================================================================
35 
36 // The ProcessContainer class combines pointers to matrix element and
37 // phase space generator with general generation info.
38 
39 class ProcessContainer : public PhysicsBase {
40 
41 public:
42 
43  // Constructor.
44  ProcessContainer(SigmaProcessPtr sigmaProcessPtrIn = 0,
45  PhaseSpacePtr phaseSpacePtrIn = 0) :
46  sigmaProcessPtr(sigmaProcessPtrIn),
47  phaseSpacePtr(phaseSpacePtrIn), resDecaysPtr(), gammaKinPtr(),
48  matchInOut(), idRenameBeams(), setLifetime(), setQuarkMass(),
49  setLeptonMass(), idNewM(), smearHadronMass(), idSmearHadIn(),
50  idSmearHadrons(), mRecalculate(), mNewM(), isLHA(), isNonDiff(),
51  isResolved(), isDiffA(), isDiffB(), isDiffC(), isQCD3body(),
52  allowNegSig(), isSameSave(), increaseMaximum(), canVetoResDecay(),
53  lhaStrat(), lhaStratAbs(), processCode(), useStrictLHEFscales(),
54  isAsymLHA(), betazLHA(), newSigmaMx(), nTry(), nSel(), nAcc(),
55  nTryStat(), sigmaMx(), sigmaSgn(), sigmaSum(), sigma2Sum(), sigmaNeg(),
56  sigmaAvg(), sigmaFin(), deltaFin(), weightNow(), wtAccSum(),
57  beamAhasResGamma(), beamBhasResGamma(), beamHasResGamma(),
58  beamHasGamma(), beamAgammaMode(), beamBgammaMode(), gammaModeEvent(),
59  approximatedGammaFlux(), reweightNLO(), doMerging(), nTryRequested(),
60  nSelRequested(), nAccRequested(), sigmaTemp(), sigma2Temp(),
61  normVar3() {}
62 
63  // Initialize phase space and counters.
64  bool init(bool isFirst, ResonanceDecays* resDecaysPtrIn,
65  SLHAinterface* slhaInterfacePtr, GammaKinematics* gammaKinPtrIn);
66 
67  // Store or replace Les Houches pointer.
68  void setLHAPtr( LHAupPtr lhaUpPtrIn, ParticleData* particleDataPtrIn = 0,
69  Settings* settingsPtrIn = 0, Rndm* rndmPtrIn = 0)
70  {lhaUpPtr = lhaUpPtrIn; setLifetime = 0;
71  if (settingsPtrIn && rndmPtrIn) {
72  rndmPtr = rndmPtrIn;
73  setLifetime = settingsPtrIn->mode("LesHouches:setLifetime");
74  }
75  if (particleDataPtrIn != 0) particleDataPtr = particleDataPtrIn;
76  if (sigmaProcessPtr != 0) sigmaProcessPtr->setLHAPtr(lhaUpPtr);
77  if (phaseSpacePtr != 0) phaseSpacePtr->setLHAPtr(lhaUpPtr);}
78 
79  // Switch to new beam particle identities; for similar hadrons only.
80  void updateBeamIDs() {phaseSpacePtr->updateBeamIDs();}
81 
82  // Update the CM energy of the event.
83  void newECM(double eCM) {phaseSpacePtr->newECM(eCM);}
84 
85  // Generate a trial event; accepted or not.
86  bool trialProcess();
87 
88  // Pick flavours and colour flow of process.
89  bool constructState();
90 
91  // Give the hard subprocess (with option for a second hard subprocess).
92  bool constructProcess( Event& process, bool isHardest = true);
93 
94  // Give the Les Houches decay chain for external resonance.
95  bool constructDecays( Event& process);
96 
97  // Do resonance decays.
98  bool decayResonances( Event& process);
99 
100  // Accumulate statistics after user veto.
101  void accumulate();
102 
103  // Reset statistics on events generated so far.
104  void reset();
105 
106  // Set whether (photon) beam is resolved or unresolved.
107  // Method propagates the choice of photon process type to beam pointers.
108  void setBeamModes(bool setVMD = false, bool isSampled = true);
109 
110  // Process name and code, and the number of final-state particles.
111  string name() const {return sigmaProcessPtr->name();}
112  int code() const {return sigmaProcessPtr->code();}
113  int nFinal() const {return sigmaProcessPtr->nFinal();}
114  bool isSUSY() const {return sigmaProcessPtr->isSUSY();}
115  bool isNonDiffractive() const {return isNonDiff;}
116  bool isSoftQCD() const {return (code() > 100 && code() < 107);}
117  bool isElastic() const {return (code() == 102);}
118 
119  // Member functions for info on generation process.
120  bool newSigmaMax() const {return newSigmaMx;}
121  double sigmaMax() const {return sigmaMx;}
122  long nTried() const {return nTry;}
123  long nSelected() const {return nSel;}
124  long nAccepted() const {return nAcc;}
125  double weightSum() const {return wtAccSum;}
126  double sigmaSelMC( bool doAccumulate = true)
127  { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return sigmaAvg;}
128  double sigmaMC( bool doAccumulate = true)
129  { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return sigmaFin;}
130  double deltaMC( bool doAccumulate = true)
131  { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return deltaFin;}
132 
133  // New value for switched beam identity or energy (for SoftQCD processes).
134  double sigmaMaxSwitch() {sigmaMx = phaseSpacePtr->sigmaMaxSwitch();
135  return sigmaMx;}
136 
137  // Some kinematics quantities.
138  int id1() const {return sigmaProcessPtr->id(1);}
139  int id2() const {return sigmaProcessPtr->id(2);}
140  double x1() const {return phaseSpacePtr->x1();}
141  double x2() const {return phaseSpacePtr->x2();}
142  double Q2Fac() const {return sigmaProcessPtr->Q2Fac();}
143  double mHat() const {return sqrtpos(phaseSpacePtr->sHat());}
144  double pTHat() const {return phaseSpacePtr->pTHat();}
145 
146  // Tell whether container is for Les Houches events.
147  bool isLHAContainer() const {return isLHA;}
148  int lhaStrategy() const {return lhaStrat;}
149 
150  // Info on Les Houches events.
151  int codeLHASize() const {return codeLHA.size();}
152  int subCodeLHA(int i) const {return codeLHA[i];}
153  long nTriedLHA(int i) const {return nTryLHA[i];}
154  long nSelectedLHA(int i) const {return nSelLHA[i];}
155  long nAcceptedLHA(int i) const {return nAccLHA[i];}
156 
157  // When two hard processes set or get info whether process is matched.
158  void isSame( bool isSameIn) { isSameSave = isSameIn;}
159  bool isSame() const {return isSameSave;}
160 
161 private:
162 
163  // Constants: could only be changed in the code itself.
164  static const int N12SAMPLE, N3SAMPLE;
165 
166  // Pointer to the subprocess matrix element. Mark if external.
167  SigmaProcessPtr sigmaProcessPtr;
168 
169  // Pointer to the phase space generator.
170  PhaseSpacePtr phaseSpacePtr;
171 
172  // Pointer to ResonanceDecays object for sequential resonance decays.
173  ResonanceDecays* resDecaysPtr;
174 
175  // Pointer to LHAup for generating external events.
176  LHAupPtr lhaUpPtr;
177 
178  // Pointer to the phase space generator of photons from leptons.
179  GammaKinematics* gammaKinPtr;
180 
181  // Possibility to modify Les Houches input.
182  bool matchInOut;
183  int idRenameBeams, setLifetime, setQuarkMass, setLeptonMass, idNewM[9],
184  smearHadronMass;
185  vector<int> idSmearHadIn, idSmearHadrons;
186  double mRecalculate, mNewM[9];
187 
188  // Info on process.
189  bool isLHA, isNonDiff, isResolved, isDiffA, isDiffB, isDiffC, isQCD3body,
190  allowNegSig, isSameSave, increaseMaximum, canVetoResDecay;
191  int lhaStrat, lhaStratAbs, processCode;
192  bool useStrictLHEFscales;
193 
194  // Boost Les Houches events to CM frame (when originally asymmetric).
195  bool isAsymLHA;
196  double betazLHA;
197 
198  // Statistics on generation process. (Long integers just in case.)
199  bool newSigmaMx;
200  long nTry, nSel, nAcc, nTryStat;
201  double sigmaMx, sigmaSgn, sigmaSum, sigma2Sum, sigmaNeg, sigmaAvg,
202  sigmaFin, deltaFin, weightNow, wtAccSum;
203 
204  // Flags to store whether beam has a (un)resolved photon.
205  bool beamAhasResGamma, beamBhasResGamma, beamHasResGamma, beamHasGamma;
206  int beamAgammaMode, beamBgammaMode, gammaModeEvent;
207 
208  // Use approximated photon flux for process sampling.
209  bool approximatedGammaFlux;
210 
211  // Whether the cross section should be reweighted to the NLO cross section.
212  bool reweightNLO;
213 
214  // Check if merging is enabled.
215  bool doMerging;
216 
217  // Statistics for Les Houches event classification.
218  vector<int> codeLHA;
219  vector<long> nTryLHA, nSelLHA, nAccLHA;
220 
221  // More fine-grained event counting.
222  long nTryRequested, nSelRequested, nAccRequested;
223 
224  // Temporary summand for handling (weighted) events when vetoes are applied.
225  double sigmaTemp, sigma2Temp, normVar3;
226 
227  // Estimate integrated cross section and its uncertainty.
228  void sigmaDelta();
229 
230 };
231 
232 //==========================================================================
233 
234 // The SetupContainers class turns the list of user-requested processes
235 // into a vector of ProcessContainer objects, each with a process.
236 
238 
239 public:
240 
241  // Constructor.
242  SetupContainers() : nVecA(), nVecB() {}
243 
244  // Initialization assuming all necessary data already read.
245  bool init(vector<ProcessContainer*>& containerPtrs, Info* infoPtr);
246 
247  // Initialization of a second hard process.
248  bool init2(vector<ProcessContainer*>& container2Ptrs, Info* infoPtr);
249 
250 private:
251 
252  // Methods to check that outgoing SUSY particles are allowed ones.
253  void setupIdVecs( Settings& settings);
254  bool allowIdVals( int idCheck1, int idCheck2);
255 
256  // Arrays of allowed outgoing SUSY particles and their lengths.
257  vector<int> idVecA, idVecB;
258  int nVecA, nVecB;
259 
260  // Helper class to setup onia production.
261  SigmaOniaSetup charmonium, bottomonium;
262 
263 };
264 
265 //==========================================================================
266 
267 } // end namespace Pythia8
268 
269 #endif // Pythia8_ProcessContainer_H
ProcessContainer(SigmaProcessPtr sigmaProcessPtrIn=0, PhaseSpacePtr phaseSpacePtrIn=0)
Constructor.
Definition: ProcessContainer.h:44
void accumulate()
Accumulate statistics after user veto.
Definition: ProcessContainer.cc:461
double sigmaMaxSwitch()
New value for switched beam identity or energy (for SoftQCD processes).
Definition: ProcessContainer.h:134
Definition: PhysicsBase.h:26
Definition: Info.h:45
void newECM(double eCM)
Update the CM energy of the event.
Definition: ProcessContainer.h:83
The Event class holds all info on the generated event.
Definition: Event.h:408
bool trialProcess()
Generate a trial event; accepted or not.
Definition: ProcessContainer.cc:276
Definition: ProcessContainer.h:39
Definition: SLHAinterface.h:27
bool constructProcess(Event &process, bool isHardest=true)
Give the hard subprocess (with option for a second hard subprocess).
Definition: ProcessContainer.cc:544
void isSame(bool isSameIn)
When two hard processes set or get info whether process is matched.
Definition: ProcessContainer.h:158
bool init(bool isFirst, ResonanceDecays *resDecaysPtrIn, SLHAinterface *slhaInterfacePtr, GammaKinematics *gammaKinPtrIn)
Initialize phase space and counters.
Definition: ProcessContainer.cc:48
Rndm * rndmPtr
Pointer to the random number generator.
Definition: PhysicsBase.h:97
SetupContainers()
Constructor.
Definition: ProcessContainer.h:242
A helper class used to setup the onia processes.
Definition: SigmaOnia.h:63
int id1() const
Some kinematics quantities.
Definition: ProcessContainer.h:138
ParticleData * particleDataPtr
Pointer to the particle data table.
Definition: PhysicsBase.h:88
bool newSigmaMax() const
Member functions for info on generation process.
Definition: ProcessContainer.h:120
Definition: Basics.h:393
int codeLHASize() const
Info on Les Houches events.
Definition: ProcessContainer.h:151
void setLHAPtr(LHAupPtr lhaUpPtrIn, ParticleData *particleDataPtrIn=0, Settings *settingsPtrIn=0, Rndm *rndmPtrIn=0)
Store or replace Les Houches pointer.
Definition: ProcessContainer.h:68
void reset()
Reset statistics on events generated so far.
Definition: ProcessContainer.cc:1394
Class to sample the virtuality and transverse momentum of emitted photons.
Definition: GammaKinematics.h:23
string name() const
Process name and code, and the number of final-state particles.
Definition: ProcessContainer.h:111
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
void setBeamModes(bool setVMD=false, bool isSampled=true)
Set the photon modes according to the process.
Definition: ProcessContainer.cc:507
Info * infoPtr
Definition: PhysicsBase.h:82
bool decayResonances(Event &process)
Do resonance decays.
Definition: ProcessContainer.cc:1334
Definition: ProcessContainer.h:237
Definition: ResonanceDecays.h:28
bool constructDecays(Event &process)
Give the Les Houches decay chain for external resonance.
Definition: ProcessContainer.cc:1210
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:423
void updateBeamIDs()
Switch to new beam particle identities; for similar hadrons only.
Definition: ProcessContainer.h:80
bool constructState()
Pick flavours and colour flow of process.
Definition: ProcessContainer.cc:487
Definition: Settings.h:196
bool isLHAContainer() const
Tell whether container is for Les Houches events.
Definition: ProcessContainer.h:147
double sqrtpos(const double &x)
Avoid problem with negative square root argument (from roundoff).
Definition: PythiaStdlib.h:184