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