PYTHIA  8.317
Pythia.h
1 // Pythia.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 main class for event generation.
7 // Pythia: provide the main user interface to everything else.
8 
9 #ifndef Pythia8_Pythia_H
10 #define Pythia8_Pythia_H
11 
12 // Version number defined for use in macros and for consistency checks.
13 #define PYTHIA_VERSION 8.317
14 #define PYTHIA_VERSION_INTEGER 8317
15 
16 // Header files for the Pythia class and for what else the user may need.
17 #include "Pythia8/Analysis.h"
18 #include "Pythia8/Basics.h"
19 #include "Pythia8/BeamParticle.h"
20 #include "Pythia8/BeamSetup.h"
21 #include "Pythia8/Event.h"
22 #include "Pythia8/FragmentationFlavZpT.h"
23 #include "Pythia8/HadronLevel.h"
24 #include "Pythia8/HadronWidths.h"
25 #include "Pythia8/Info.h"
26 #include "Pythia8/JunctionSplitting.h"
27 #include "Pythia8/LesHouches.h"
28 #include "Pythia8/Logger.h"
29 #include "Pythia8/SigmaLowEnergy.h"
30 #include "Pythia8/Merging.h"
31 #include "Pythia8/MergingHooks.h"
32 #include "Pythia8/PartonLevel.h"
33 #include "Pythia8/ParticleData.h"
34 #include "Pythia8/PartonDistributions.h"
35 #include "Pythia8/PartonSystems.h"
36 #include "Pythia8/PartonVertex.h"
37 #include "Pythia8/PhysicsBase.h"
38 #include "Pythia8/ProcessLevel.h"
39 #include "Pythia8/PythiaStdlib.h"
40 #include "Pythia8/ResonanceWidths.h"
41 #include "Pythia8/RHadrons.h"
42 #include "Pythia8/Ropewalk.h"
43 #include "Pythia8/Settings.h"
44 #include "Pythia8/ShowerModel.h"
45 #include "Pythia8/SigmaTotal.h"
46 #include "Pythia8/SimpleSpaceShower.h"
47 #include "Pythia8/SimpleTimeShower.h"
48 #include "Pythia8/SpaceShower.h"
49 #include "Pythia8/StandardModel.h"
50 #include "Pythia8/StringInteractions.h"
51 #include "Pythia8/SusyCouplings.h"
52 #include "Pythia8/SLHAinterface.h"
53 #include "Pythia8/ThermalFragmentation.h"
54 #include "Pythia8/TimeShower.h"
55 #include "Pythia8/UserHooks.h"
56 #include "Pythia8/VinciaCommon.h"
57 #include "Pythia8/Weights.h"
58 
59 namespace Pythia8 {
60 
61 //==========================================================================
62 
63 // Forward declaration of the HeavyIons and HIUserHooks classes.
64 class HeavyIons;
65 class HIUserHooks;
66 
67 // Forward declaration of PythiaParallel class, to be friended.
68 class PythiaParallel;
69 
70 // The Pythia class contains the top-level routines to generate an event.
71 
72 class Pythia {
73 
74 public:
75 
76  // Constructor. (See Pythia.cc file.)
77  Pythia(string xmlDir = "../share/Pythia8/xmldoc", bool printBanner = true);
78 
79  // Constructor to copy settings and particle database from another Pythia
80  // object instead of XML files (to speed up multiple initialisations).
81  Pythia(Settings& settingsIn, ParticleData& particleDataIn,
82  bool printBanner = true);
83 
84  // Constructor taking input from streams instead of files.
85  Pythia( istream& settingsStrings, istream& particleDataStrings,
86  bool printBanner = true);
87 
88  // Destructor.
89  ~Pythia() {}
90 
91  // Copy and = constructors cannot be used.
92  Pythia(const Pythia&) = delete;
93  Pythia& operator=(const Pythia&) = delete;
94 
95  // Check consistency of version numbers (called by constructors).
96  bool checkVersion();
97 
98  // Read in one update for a setting or particle data from a single line.
99  inline bool readString(string line, bool warn = true,
100  int subrun = SUBRUNDEFAULT) {
101  return isConstructed ? settings.readString(line, warn, subrun) : false;}
102 
103  // Read in updates for settings or particle data from user-defined file.
104  inline bool readFile(string fileName, bool warn = true,
105  int subrun = SUBRUNDEFAULT) {
106  return isConstructed ? settings.readFile(fileName, warn, subrun) : false;}
107  inline bool readFile(string fileName, int subrun) {
108  return readFile(fileName, true, subrun);}
109  inline bool readFile(istream& is = cin, bool warn = true,
110  int subrun = SUBRUNDEFAULT) {
111  return isConstructed ? settings.readFile(is, warn, subrun) : false;}
112  inline bool readFile(istream& is, int subrun) {
113  return readFile(is, true, subrun);}
114 
115  // Possibility to pass in pointers to PDF's.
116  inline bool setPDFPtr(PDFPtr pdfAPtrIn, PDFPtr pdfBPtrIn,
117  PDFPtr pdfHardAPtrIn = nullptr, PDFPtr pdfHardBPtrIn = nullptr,
118  PDFPtr pdfPomAPtrIn = nullptr, PDFPtr pdfPomBPtrIn = nullptr,
119  PDFPtr pdfGamAPtrIn = nullptr, PDFPtr pdfGamBPtrIn = nullptr,
120  PDFPtr pdfHardGamAPtrIn = nullptr, PDFPtr pdfHardGamBPtrIn = nullptr,
121  PDFPtr pdfUnresAPtrIn = nullptr, PDFPtr pdfUnresBPtrIn = nullptr,
122  PDFPtr pdfUnresGamAPtrIn = nullptr, PDFPtr pdfUnresGamBPtrIn = nullptr,
123  PDFPtr pdfVMDAPtrIn = nullptr, PDFPtr pdfVMDBPtrIn = nullptr) {
124  return beamSetup.setPDFPtr( pdfAPtrIn, pdfBPtrIn, pdfHardAPtrIn,
125  pdfHardBPtrIn, pdfPomAPtrIn, pdfPomBPtrIn, pdfGamAPtrIn, pdfGamBPtrIn,
126  pdfHardGamAPtrIn, pdfHardGamBPtrIn, pdfUnresAPtrIn, pdfUnresBPtrIn,
127  pdfUnresGamAPtrIn, pdfUnresGamBPtrIn, pdfVMDAPtrIn, pdfVMDBPtrIn); }
128  inline bool setPDFAPtr( PDFPtr pdfAPtrIn ) {
129  return beamSetup.setPDFAPtr( pdfAPtrIn); }
130  inline bool setPDFBPtr( PDFPtr pdfBPtrIn ) {
131  return beamSetup.setPDFBPtr( pdfBPtrIn); }
132 
133  // Set photon fluxes externally. Used with option "PDF:lepton2gammaSet = 2".
134  inline bool setPhotonFluxPtr(PDFPtr photonFluxAIn, PDFPtr photonFluxBIn) {
135  return beamSetup.setPhotonFluxPtr( photonFluxAIn, photonFluxBIn); }
136 
137  // Possibility to pass in pointer to external LHA-interfaced generator.
138  inline bool setLHAupPtr(LHAupPtr lhaUpPtrIn) {
139  lhaUpPtr = lhaUpPtrIn;
140  useNewLHA = false;
141  return beamSetup.setLHAupPtr(lhaUpPtrIn);}
142 
143  // Possibility to pass in pointer for external handling of some decays.
144  inline bool setDecayPtr(DecayHandlerPtr decayHandlePtrIn,
145  vector<int> handledParticlesIn = {}) {
146  decayHandlePtr = decayHandlePtrIn;
147  handledParticles = handledParticlesIn.size() == 0 ?
148  decayHandlePtrIn->handledParticles() : handledParticlesIn; return true;}
149 
150  // Possibility to pass in pointer for external random number generation.
151  inline bool setRndmEnginePtr(RndmEnginePtr rndmEnginePtrIn) {
152  return rndm.rndmEnginePtr(rndmEnginePtrIn);}
153 
154  // Possibility to pass in pointer for user hooks.
155  inline bool setUserHooksPtr(UserHooksPtr userHooksPtrIn) {
156  userHooksPtr = userHooksPtrIn; return true;}
157 
158  // Possibility to add further pointers to allow multiple user hooks.
159  inline bool addUserHooksPtr(UserHooksPtr userHooksPtrIn) {
160  if ( !userHooksPtrIn ) return false;
161  if ( !userHooksPtr ) return setUserHooksPtr(userHooksPtrIn);
162  shared_ptr<UserHooksVector> uhv =
163  dynamic_pointer_cast<UserHooksVector>(userHooksPtr);
164  if ( !uhv ) { uhv = make_shared<UserHooksVector>();
165  uhv->hooks.push_back(userHooksPtr); userHooksPtr = uhv; }
166  uhv->hooks.push_back(userHooksPtrIn); return true;}
167 
168  // Possibility to insert a user hook.
169  inline bool insertUserHooksPtr(int idx, UserHooksPtr userHooksPtrIn) {
170  if ( !userHooksPtrIn || !userHooksPtr ) return false;
171  shared_ptr<UserHooksVector> uhv =
172  dynamic_pointer_cast<UserHooksVector>(userHooksPtr);
173  if ( !uhv || idx < 0 || idx > (int)uhv->hooks.size() ) return false;
174  uhv->hooks.insert(uhv->hooks.begin() + idx, userHooksPtrIn); return true;}
175 
176  // Possibility to pass in pointer for full merging class.
177  inline bool setMergingPtr(MergingPtr mergingPtrIn) {
178  mergingPtr = mergingPtrIn; return true;}
179 
180  // Possibility to pass in pointer for merging hooks.
181  inline bool setMergingHooksPtr(MergingHooksPtr mergingHooksPtrIn) {
182  mergingHooksPtr = mergingHooksPtrIn; return true;}
183 
184  // Possibility to pass in pointer for beam shape.
185  inline bool setBeamShapePtr(BeamShapePtr beamShapePtrIn) {
186  return beamSetup.setBeamShapePtr(beamShapePtrIn);}
187 
188  // Possibility to pass in pointer for external cross section,
189  // with option to include external phase-space generator.
190  inline bool setSigmaPtr(SigmaProcessPtr sigmaPtrIn,
191  PhaseSpacePtr phaseSpacePtrIn = nullptr) {
192  sigmaPtrs.resize(0), phaseSpacePtrs.resize(0);
193  sigmaPtrs.push_back(sigmaPtrIn);
194  phaseSpacePtrs.push_back(phaseSpacePtrIn); return true;}
195 
196  // Possibility to add further pointers to allow for multiple cross sections.
197  inline bool addSigmaPtr(SigmaProcessPtr sigmaPtrIn,
198  PhaseSpacePtr phaseSpacePtrIn = nullptr) {
199  sigmaPtrs.push_back(sigmaPtrIn);
200  phaseSpacePtrs.push_back(phaseSpacePtrIn); return true;}
201 
202  // Possibility to insert further pointers to allow for multiple
203  // cross sections.
204  inline bool insertSigmaPtr(int idx, SigmaProcessPtr sigmaPtrIn,
205  PhaseSpacePtr phaseSpacePtrIn = nullptr) {
206  if (idx < 0 || idx > (int)sigmaPtrs.size()) return false;
207  sigmaPtrs.insert(sigmaPtrs.begin() + idx, sigmaPtrIn);
208  phaseSpacePtrs.insert(phaseSpacePtrs.begin() + idx, phaseSpacePtrIn);
209  return true;}
210 
211  // Possibility to pass in pointer for external resonance.
212  inline bool setResonancePtr(ResonanceWidthsPtr resonancePtrIn) {
213  resonancePtrs.resize(0);
214  resonancePtrs.push_back( resonancePtrIn); return true;}
215 
216  // Possibility to add further pointers to allow for multiple resonances.
217  inline bool addResonancePtr(ResonanceWidthsPtr resonancePtrIn) {
218  resonancePtrs.push_back( resonancePtrIn); return true;}
219 
220  // Possibility to insert further pointers to allow for multiple resonances.
221  inline bool insertResonancePtr(int idx, ResonanceWidthsPtr resonancePtrIn) {
222  if (idx < 0 || idx > (int)resonancePtrs.size()) return false;
223  resonancePtrs.insert( resonancePtrs.begin() + idx, resonancePtrIn);
224  return true;}
225 
226  // Possibility to pass in pointer for external showers.
227  inline bool setShowerModelPtr(ShowerModelPtr showerModelPtrIn) {
228  showerModelPtr = showerModelPtrIn; return true;}
229 
230  // Possibility to pass in pointer for external fragmentation model.
231  inline bool setFragmentationPtr(FragmentationModelPtr fragmentationPtrIn) {
232  fragPtrs.resize(0);
233  fragPtrs.push_back(fragmentationPtrIn); return true;}
234 
235  // Possibility to allow for multiple external fragmentation models.
236  inline bool addFragmentationPtr(FragmentationModelPtr fragmentationPtrIn) {
237  fragPtrs.push_back(fragmentationPtrIn); return true;}
238 
239  // Possibility to insert external fragmentation model, in specific position.
240  inline bool insertFragmentationPtr(int idx,
241  FragmentationModelPtr fragmentationPtrIn) {
242  if (idx < 0 || idx > (int)fragPtrs.size()) return false;
243  fragPtrs.insert( fragPtrs.begin() + idx, fragmentationPtrIn);
244  return true;}
245 
246  // Possibility to pass in pointer for modelling of heavy ion collisions.
247  inline bool setHeavyIonsPtr(HeavyIonsPtr heavyIonsPtrIn) {
248  heavyIonsPtr = heavyIonsPtrIn; return true;}
249 
250  // Possibility to pass a HIUserHooks pointer for modifying the
251  // behavior of the heavy ion modelling.
252  inline bool setHIHooks(HIUserHooksPtr hiHooksPtrIn) {
253  hiHooksPtr = hiHooksPtrIn; return true; }
254 
255  // Possibility to get the pointer to a object modelling heavy ion
256  // collisions.
257  inline HeavyIonsPtr getHeavyIonsPtr() { return heavyIonsPtr;}
258 
259  // Possibility to access the pointer to the BeamShape object.
260  inline BeamShapePtr getBeamShapePtr() { return beamSetup.getBeamShapePtr();}
261 
262  // Possibility to get the pointer to the parton-shower model.
263  inline ShowerModelPtr getShowerModelPtr() { return showerModelPtr;}
264 
265  // Possibility to get the pointer to the LHA accessor.
266  inline LHAupPtr getLHAupPtr() { return lhaUpPtr;}
267 
268  // Possibility to pass in pointer for setting of parton space-time vertices.
269  inline bool setPartonVertexPtr( PartonVertexPtr partonVertexPtrIn) {
270  partonVertexPtr = partonVertexPtrIn; return true;}
271 
272  // Initialize.
273  bool init();
274 
275  // Generate the next event.
276  inline bool next() { return next(0); }
277  bool next(int procTypeIn);
278 
279  // Switch to new beam particle identities; for similar hadrons only.
280  bool setBeamIDs( int idAin, int idBin = 0);
281 
282  // Switch beam kinematics.
283  bool setKinematics(double eCMIn);
284  bool setKinematics(double eAIn, double eBIn);
285  bool setKinematics(double pxAIn, double pyAIn, double pzAIn,
286  double pxBIn, double pyBIn, double pzBIn);
287  bool setKinematics(Vec4 pAIn, Vec4 pBIn);
288 
289  // Generate only a single timelike shower as in a decay.
290  inline int forceTimeShower( int iBeg, int iEnd, double pTmax,
291  int nBranchMax = 0) {
292  if (!isInit) {
293  logger.ERROR_MSG("Pythia is not properly initialized");
294  return 0;
295  }
297  infoPrivate.setScalup( 0, pTmax);
298  return timesDecPtr->shower(iBeg, iEnd, event, pTmax, nBranchMax);}
299 
300  // Generate only the hadronization/decay stage.
301  bool forceHadronLevel( bool findJunctions = true);
302 
303  // Special routine to allow more decays if on/off switches changed.
304  bool moreDecays(bool allowPartons = false) {
305  return hadronLevel.moreDecays(event, allowPartons);}
306  bool moreDecays(int index, bool allowPartons = false) {
307  return hadronLevel.decay(index, event, allowPartons);}
308 
309  // Special routine to force R-hadron decay when not done before.
310  bool forceRHadronDecays() {return doRHadronDecays();}
311 
312  // Do a low-energy collision between two hadrons in the event record.
313  inline bool doLowEnergyProcess(int i1, int i2, int procTypeIn) {
314  if (!isInit) {
315  logger.ERROR_MSG("Pythia is not properly initialized"); return false; }
316  return hadronLevel.doLowEnergyProcess( i1, i2, procTypeIn, event); }
317 
318  // Get total cross section for two hadrons in the event record or standalone.
319  inline double getSigmaTotal() {
320  return getSigmaTotal(beamSetup.idA, beamSetup.idB, beamSetup.eCM, 0);}
321  inline double getSigmaTotal(double eCM12, int mixLoHi = 0) {
322  return getSigmaTotal(beamSetup.idA, beamSetup.idB, eCM12, mixLoHi);}
323  inline double getSigmaTotal(int id1, int id2, double eCM12,
324  int mixLoHi = 0) {
325  return getSigmaTotal(id1, id2, eCM12, particleData.m0(id1),
326  particleData.m0(id2), mixLoHi); }
327  inline double getSigmaTotal(int id1, int id2, double eCM12, double m1,
328  double m2, int mixLoHi = 0) {
329  if (!isInit) {
330  logger.ERROR_MSG("Pythia is not properly initialized"); return 0.; }
331  return sigmaCmb.sigmaTotal(id1, id2, eCM12, m1, m2, mixLoHi); }
332 
333  // Get partial (elastic, diffractive, non-diffractive, ...) cross sections
334  // for two hadrons in the event record or standalone.
335  inline double getSigmaPartial(int procTypeIn) {
336  return getSigmaPartial(beamSetup.idA, beamSetup.idB, beamSetup.eCM,
337  procTypeIn, 0); }
338  inline double getSigmaPartial(double eCM12, int procTypeIn,
339  int mixLoHi = 0) {return getSigmaPartial(
340  beamSetup.idA, beamSetup.idB, eCM12, procTypeIn, mixLoHi); }
341  inline double getSigmaPartial(int id1, int id2, double eCM12, int procTypeIn,
342  int mixLoHi = 0) {return getSigmaPartial(id1, id2, eCM12,
343  particleData.m0(id1), particleData.m0(id2), procTypeIn, mixLoHi); }
344  inline double getSigmaPartial(int id1, int id2, double eCM12, double m1,
345  double m2, int procTypeIn, int mixLoHi = 0) {
346  if (!isInit) {
347  logger.ERROR_MSG("Pythia is not properly initialized"); return 0.;}
348  return sigmaCmb.sigmaPartial(
349  id1, id2, eCM12, m1, m2, procTypeIn, mixLoHi);}
350 
351  // Return a parton density set among list of possibilities.
352  inline PDFPtr getPDFPtr(int idIn, int sequence = 1, string beam = "A",
353  bool resolved = true) {
354  return beamSetup.getPDFPtr( idIn, sequence, beam, resolved); }
355 
356  // List the current Les Houches event.
357  inline void LHAeventList() { if (lhaUpPtr != 0) lhaUpPtr->listEvent();}
358 
359  // Skip a number of Les Houches events at input.
360  inline bool LHAeventSkip(int nSkip) {
361  if (lhaUpPtr != 0) return lhaUpPtr->skipEvent(nSkip);
362  return false;}
363 
364  // Main routine to provide final statistics on generation.
365  void stat();
366 
367  // Read in settings values: shorthand, not new functionality.
368  bool flag(string key) {return settings.flag(key);}
369  int mode(string key) {return settings.mode(key);}
370  double parm(string key) {return settings.parm(key);}
371  string word(string key) {return settings.word(key);}
372 
373  // The event record for the parton-level central process.
375 
376  // The event record for the complete event history.
377  Event event = {};
378 
379  // Public information and statistic on the generation.
380  const Info& info = infoPrivate;
381 
382  // Logger: for diagnostic messages, errors, statistics, etc.
384 
385  // Settings: databases of flags/modes/parms/words to control run.
387 
388  // ParticleData: the particle data table/database.
390 
391  // Random number generator.
392  Rndm rndm = {};
393 
394  // Standard Model couplings, including alphaS and alphaEM.
396 
397  // SUSY couplings.
399 
400  // SLHA Interface
402 
403  // The partonic content of each subcollision system (auxiliary to event).
405 
406  // Merging object as wrapper for matrix element merging routines.
407  MergingPtr mergingPtr = {};
408 
409  // Pointer to MergingHooks object for user interaction with the merging.
410  // MergingHooks also more generally steers the matrix element merging.
411  MergingHooksPtr mergingHooksPtr;
412 
413  // Pointer to a HeavyIons object for generating heavy ion collisions
414  HeavyIonsPtr heavyIonsPtr = {};
415 
416  // Pointer to a HIUserHooks object to modify heavy ion modelling.
417  HIUserHooksPtr hiHooksPtr = {};
418 
419  // HadronWidths: the hadron widths data table/database.
421 
422  // The two incoming beams.
423  const BeamParticle& beamA = beamSetup.beamA;
424  const BeamParticle& beamB = beamSetup.beamB;
425 
426 private:
427 
428  // Friend PythiaParallel to give full access to underlying info.
429  friend class PythiaParallel;
430  friend class PhysicsBase;
431 
432  // The collector of all event generation weights that should eventually
433  // be transferred to the final output.
434  WeightContainer weightContainer = {};
435 
436  // The main keeper/collector of information, accessible from all
437  // PhysicsBase objects. The information is available from the
438  // outside through the public info object.
439  Info infoPrivate = {};
440 
441  // Initialise new Pythia object (called by constructors).
442  void initPtrs();
443 
444  // Initialise fragmentation model objects.
445  void initFragPtrs();
446 
447  // Initialise user provided plugins.
448  void initPlugins();
449 
450  // Functions to be called at the beginning and end of a next() call.
451  void beginEvent();
452  void endEvent(PhysicsBase::Status status);
453 
454  // Register a PhysicsBase object and give it a pointer to the info object.
455  inline void registerPhysicsBase(PhysicsBase &pb) {
456  if (find(physicsPtrs.begin(), physicsPtrs.end(), &pb) != physicsPtrs.end())
457  return;
458  pb.initInfoPtr(infoPrivate);
459  physicsPtrs.push_back(&pb);
460  }
461 
462  // If new pointers are set in Info propagate this to all
463  // PhysicsBase objects.
464  void pushInfo();
465 
466  // Constants: could only be changed in the code itself.
467  static const double VERSIONNUMBERHEAD, VERSIONNUMBERCODE;
468  // Maximum number of tries to produce parton level from given input.
469  // Negative integer to denote that no subrun has been set.
470  static const int NTRY = 10;
471 
472  // Initialization data, extracted from database.
473  string xmlPath = {};
474  bool doProcessLevel = {}, doPartonLevel = {}, doHadronLevel = {},
475  doLowEnergy = {}, doSoftQCDall = {}, doSoftQCDinel = {},
476  doCentralDiff = {}, doDiffraction = {}, doSoftQCD = {},
477  doHardDiff = {}, doResDec = {}, doFSRinRes = {}, decayRHadrons = {},
478  doPartonVertex = {}, abortIfVeto = {}, checkEvent = {},
479  checkHistory = {}, doNonPert = {};
480  int nErrList = {};
481  double epTolErr = {}, epTolWarn = {}, mTolErr = {}, mTolWarn = {};
482 
483  // Initialization data, from init(...) call, plus some event-specific.
484  bool isConstructed = {}, isInit = {}, showSaV = {}, showMaD = {},
485  doReconnect = {}, forceHadronLevelCR = {};
486  int nCount = {}, nShowLHA = {}, nShowInfo = {}, nShowProc = {},
487  nShowEvt = {}, reconnectMode = {};
488 
489  // information for error checkout.
490  int nErrEvent = {};
491  vector<int> iErrId = {}, iErrCol = {}, iErrEpm = {}, iErrNan = {},
492  iErrNanVtx = {};
493 
494  // Setup of beams: flavours, kinematics and PDFs.
495  BeamSetup beamSetup = {};
496 
497  // LHAup object for generating external events.
498  bool doHardProc = false;
499  bool doLHA = false;
500  bool useNewLHA = false;
501  LHAupPtr lhaUpPtr = {};
502 
503  // Pointer to external decay handler and list of particles it handles.
504  DecayHandlerPtr decayHandlePtr = {};
505  vector<int> handledParticles = {};
506 
507  // Pointer to UserHooks object for user interaction with program.
508  UserHooksPtr userHooksPtr = {};
509  bool doVetoProcess = {}, doVetoPartons = {},
510  retryPartonLevel = {}, canVetoHadronization = {};
511 
512  // Pointer to BeamShape object for beam momentum and interaction vertex.
513  BeamShapePtr beamShapePtr = {};
514  bool doVertexSpread = {};
515  double eMinPert = {}, eWidthPert = {};
516 
517  // Pointers to external processes derived from the Pythia base classes.
518  vector<SigmaProcessPtr> sigmaPtrs = {};
519 
520  // Pointers to external phase-space generators derived from Pythia
521  // base classes.
522  vector<PhaseSpacePtr> phaseSpacePtrs = {};
523 
524  // Pointers to external calculation of resonance widths.
525  vector<ResonanceWidthsPtr> resonancePtrs = {};
526 
527  // Pointers to timelike and spacelike showers, including Vincia.
528  // Note, the showerModelPtr must be declared before the
529  // individual shower pointers, partonLevel, and hadronLevel. This
530  // ensures shared pointers that belong to showerModelPtr are
531  // unloaded correctly.
532  ShowerModelPtr showerModelPtr = {};
533  TimeShowerPtr timesDecPtr = {};
534  TimeShowerPtr timesPtr = {};
535  SpaceShowerPtr spacePtr = {};
536 
537  // Pointers to fragmentation models.
538  vector<FragmentationModelPtr> fragPtrs = {};
539 
540  // Pointer to assign space-time vertices during parton evolution.
541  PartonVertexPtr partonVertexPtr;
542 
543  // The main generator class to define the core process of the event.
544  ProcessLevel processLevel = {};
545 
546  // The main generator class to produce the parton level of the event.
547  PartonLevel partonLevel = {};
548 
549  // The main generator class to perform trial showers of the event.
550  PartonLevel trialPartonLevel = {};
551 
552  // Flags for defining the merging scheme.
553  bool doMerging = {};
554 
555  // The StringInteractions class.
556  StringIntPtr stringInteractionsPtr;
557 
558  // The Colour reconnection class.
559  ColRecPtr colourReconnectionPtr = {};
560 
561  // The junction splitting class.
562  JunctionSplitting junctionSplitting = {};
563 
564  // The main generator class to produce the hadron level of the event.
565  HadronLevel hadronLevel = {};
566 
567  // The total cross section classes are used both on process and parton level.
568  SigmaTotal sigmaTot = {};
569  SigmaLowEnergy sigmaLowEnergy;
570  NucleonExcitations nucleonExcitations = {};
571  SigmaCombined sigmaCmb = {};
572 
573  // The fragmentation pointer is used in low energy processes and HadronLevel.
574  LundFragmentationPtr fragPtr{};
575 
576  // The RHadrons class is used both at PartonLevel and HadronLevel.
577  RHadronsPtr rHadronsPtr{};
578 
579  // Flags for handling generation of heavy ion collisons.
580  bool hasHeavyIons = {}, doHeavyIons = {};
581 
582  // Write the Pythia banner, with symbol and version information.
583  void banner();
584 
585  // Check that combinations of settings are allowed; change if not.
586  void checkSettings();
587 
588  // Simplified treatment for low-energy nonperturbative collisions.
589  bool nextNonPert(int procTypeIn = 0);
590 
591  // Perform R-hadron decays.
592  bool doRHadronDecays();
593 
594  // Check that the final event makes sense.
595  bool check();
596 
597  // Initialization of SLHA data.
598  bool initSLHA();
599  stringstream particleDataBuffer;
600 
601  // Keep track of and initialize all pointers to PhysicsBase-derived objects.
602  vector<PhysicsBase*> physicsPtrs = {};
603 
604  // Mutex that can be used by PhysicsBase-derived objects.
605  mutex mainMutex;
606 
607 };
608 
609 //==========================================================================
610 
611 } // end namespace Pythia8
612 
613 #endif // Pythia8_Pythia_H
Event process
The event record for the parton-level central process.
Definition: Pythia.h:374
HeavyIonsPtr getHeavyIonsPtr()
Definition: Pythia.h:257
bool doLowEnergyProcess(int i1, int i2, int procTypeIn)
Do a low-energy collision between two hadrons in the event record.
Definition: Pythia.h:313
bool setRndmEnginePtr(RndmEnginePtr rndmEnginePtrIn)
Possibility to pass in pointer for external random number generation.
Definition: Pythia.h:151
double sigmaPartial(int id1, int id2, double eCM12, double m1, double m2, int type, int mixLoHi=0)
Get partial cross sections.
Definition: SigmaLowEnergy.cc:1563
bool setPDFAPtr(PDFPtr pdfAPtrIn)
Routine to pass in pointers to PDF&#39;s. Usage optional.
Definition: BeamSetup.cc:105
LHAupPtr getLHAupPtr()
Possibility to get the pointer to the LHA accessor.
Definition: Pythia.h:266
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1615
Definition: NucleonExcitations.h:23
PDFPtr getPDFPtr(int idIn, int sequence=1, string beam="A", bool resolved=true)
Return a parton density set among list of possibilities.
Definition: Pythia.h:352
Definition: PhysicsBase.h:26
Definition: Info.h:45
bool rndmEnginePtr(RndmEnginePtr rndmEngPtrIn)
Possibility to pass in pointer for external random number generation.
Definition: Basics.cc:27
Definition: BeamSetup.h:33
Rndm rndm
Random number generator.
Definition: Pythia.h:392
bool LHAeventSkip(int nSkip)
Skip a number of Les Houches events at input.
Definition: Pythia.h:360
bool readString(string line, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in one update for a setting or particle data from a single line.
Definition: Pythia.h:99
ParticleData particleData
ParticleData: the particle data table/database.
Definition: Pythia.h:389
The Event class holds all info on the generated event.
Definition: Event.h:408
Definition: BeamParticle.h:133
bool addResonancePtr(ResonanceWidthsPtr resonancePtrIn)
Possibility to add further pointers to allow for multiple resonances.
Definition: Pythia.h:217
Definition: Weights.h:468
Definition: SigmaLowEnergy.h:135
bool moreDecays(Event &event, bool allowPartons=false)
Special routine to allow more decays if on/off switches changed.
Definition: HadronLevel.cc:341
PDFPtr getPDFPtr(int idIn, int sequence=1, string beam="A", bool resolved=true)
Return a parton density set among list of possibilities.
Definition: BeamSetup.cc:1146
bool setResonancePtr(ResonanceWidthsPtr resonancePtrIn)
Possibility to pass in pointer for external resonance.
Definition: Pythia.h:212
Definition: SLHAinterface.h:27
vector< shared_ptr< UserHooks > > hooks
The vector of user hooks.
Definition: UserHooks.h:766
Definition: SigmaTotal.h:141
bool setDecayPtr(DecayHandlerPtr decayHandlePtrIn, vector< int > handledParticlesIn={})
Possibility to pass in pointer for external handling of some decays.
Definition: Pythia.h:144
bool setPhotonFluxPtr(PDFPtr photonFluxAIn, PDFPtr photonFluxBIn)
Set photon fluxes externally. Used with option "PDF:lepton2gammaSet = 2".
Definition: BeamSetup.h:53
Definition: JunctionSplitting.h:30
bool setFragmentationPtr(FragmentationModelPtr fragmentationPtrIn)
Possibility to pass in pointer for external fragmentation model.
Definition: Pythia.h:231
Class for doing Pythia runs in parallel.
Definition: PythiaParallel.h:18
bool checkVersion()
Check consistency of version numbers (called by constructors).
Definition: Pythia.cc:385
Settings settings
Settings: databases of flags/modes/parms/words to control run.
Definition: Pythia.h:386
bool forceRHadronDecays()
Special routine to force R-hadron decay when not done before.
Definition: Pythia.h:310
void initInfoPtr(Info &infoPtrIn)
This function is called from above for physics objects used in a run.
Definition: PhysicsBase.cc:21
Definition: ProcessLevel.h:36
Gets cross sections for hadron-hadron collisions at low energies.
Definition: SigmaLowEnergy.h:22
bool setLHAupPtr(LHAupPtr lhaUpPtrIn)
Possibility to pass in pointer to external LHA-interfaced generator.
Definition: Pythia.h:138
HIUserHooksPtr hiHooksPtr
Pointer to a HIUserHooks object to modify heavy ion modelling.
Definition: Pythia.h:417
bool setKinematics(double eCMIn)
Switch beam kinematics.
Definition: Pythia.cc:1387
Event event
The event record for the complete event history.
Definition: Pythia.h:377
bool init()
Initialize.
Definition: Pythia.cc:422
bool setPDFPtr(PDFPtr pdfAPtrIn, PDFPtr pdfBPtrIn, PDFPtr pdfHardAPtrIn=nullptr, PDFPtr pdfHardBPtrIn=nullptr, PDFPtr pdfPomAPtrIn=nullptr, PDFPtr pdfPomBPtrIn=nullptr, PDFPtr pdfGamAPtrIn=nullptr, PDFPtr pdfGamBPtrIn=nullptr, PDFPtr pdfHardGamAPtrIn=nullptr, PDFPtr pdfHardGamBPtrIn=nullptr, PDFPtr pdfUnresAPtrIn=nullptr, PDFPtr pdfUnresBPtrIn=nullptr, PDFPtr pdfUnresGamAPtrIn=nullptr, PDFPtr pdfUnresGamBPtrIn=nullptr, PDFPtr pdfVMDAPtrIn=nullptr, PDFPtr pdfVMDBPtrIn=nullptr)
Possibility to pass in pointers to PDF&#39;s.
Definition: Pythia.h:116
bool setMergingPtr(MergingPtr mergingPtrIn)
Possibility to pass in pointer for full merging class.
Definition: Pythia.h:177
bool moreDecays(bool allowPartons=false)
Special routine to allow more decays if on/off switches changed.
Definition: Pythia.h:304
bool setPDFPtr(PDFPtr pdfAPtrIn, PDFPtr pdfBPtrIn, PDFPtr pdfHardAPtrIn=nullptr, PDFPtr pdfHardBPtrIn=nullptr, PDFPtr pdfPomAPtrIn=nullptr, PDFPtr pdfPomBPtrIn=nullptr, PDFPtr pdfGamAPtrIn=nullptr, PDFPtr pdfGamBPtrIn=nullptr, PDFPtr pdfHardGamAPtrIn=nullptr, PDFPtr pdfHardGamBPtrIn=nullptr, PDFPtr pdfUnresAPtrIn=nullptr, PDFPtr pdfUnresBPtrIn=nullptr, PDFPtr pdfUnresGamAPtrIn=nullptr, PDFPtr pdfUnresGamBPtrIn=nullptr, PDFPtr pdfVMDAPtrIn=nullptr, PDFPtr pdfVMDBPtrIn=nullptr)
Possibility to pass in pointers to PDF&#39;s.
Definition: BeamSetup.cc:21
const BeamParticle & beamA
The two incoming beams.
Definition: Pythia.h:423
~Pythia()
Destructor.
Definition: Pythia.h:89
Definition: Logger.h:23
SLHAinterface slhaInterface
SLHA Interface.
Definition: Pythia.h:401
bool insertResonancePtr(int idx, ResonanceWidthsPtr resonancePtrIn)
Possibility to insert further pointers to allow for multiple resonances.
Definition: Pythia.h:221
bool insertUserHooksPtr(int idx, UserHooksPtr userHooksPtrIn)
Possibility to insert a user hook.
Definition: Pythia.h:169
HeavyIonsPtr heavyIonsPtr
Pointer to a HeavyIons object for generating heavy ion collisions.
Definition: Pythia.h:414
bool readString(string line, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in one update from a single line.
Definition: Settings.cc:384
bool setPartonVertexPtr(PartonVertexPtr partonVertexPtrIn)
Possibility to pass in pointer for setting of parton space-time vertices.
Definition: Pythia.h:269
BeamShapePtr getBeamShapePtr()
Possibility to access the pointer to the BeamShape object.
Definition: Pythia.h:260
CoupSUSY coupSUSY
SUSY couplings.
Definition: Pythia.h:398
bool readFile(string fileName, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in updates from a user-defined file.
Definition: Settings.cc:701
bool addFragmentationPtr(FragmentationModelPtr fragmentationPtrIn)
Possibility to allow for multiple external fragmentation models.
Definition: Pythia.h:236
bool setBeamIDs(int idAin, int idBin=0)
Switch to new beam particle identities; for similar hadrons only.
Definition: Pythia.cc:1365
bool decay(int iDec, Event &event, bool allowPartons=true)
Definition: HadronLevel.h:68
CoupSM coupSM
Standard Model couplings, including alphaS and alphaEM.
Definition: Pythia.h:395
BeamParticle beamA
The two incoming beams.
Definition: BeamSetup.h:131
Definition: HadronLevel.h:44
Definition: Basics.h:393
PartonSystems partonSystems
The partonic content of each subcollision system (auxiliary to event).
Definition: Pythia.h:404
ShowerModelPtr getShowerModelPtr()
Possibility to get the pointer to the parton-shower model.
Definition: Pythia.h:263
bool insertFragmentationPtr(int idx, FragmentationModelPtr fragmentationPtrIn)
Possibility to insert external fragmentation model, in specific position.
Definition: Pythia.h:240
bool setBeamShapePtr(BeamShapePtr beamShapePtrIn)
Possibility to pass in pointer for beam shape.
Definition: Pythia.h:185
Definition: PartonLevel.h:45
Status
Enumerate the different status codes the event generation can have.
Definition: PhysicsBase.h:31
Definition: HadronWidths.h:22
bool setMergingHooksPtr(MergingHooksPtr mergingHooksPtrIn)
Possibility to pass in pointer for merging hooks.
Definition: Pythia.h:181
MergingHooksPtr mergingHooksPtr
Definition: Pythia.h:411
void stat()
Main routine to provide final statistics on generation.
Definition: Pythia.cc:1705
bool doLowEnergyProcess(int i1, int i2, int procTypeIn, Event &event)
Special routine to do a low-energy hadron-hadron scattering.
Definition: HadronLevel.h:83
bool addSigmaPtr(SigmaProcessPtr sigmaPtrIn, PhaseSpacePtr phaseSpacePtrIn=nullptr)
Possibility to add further pointers to allow for multiple cross sections.
Definition: Pythia.h:197
bool setHeavyIonsPtr(HeavyIonsPtr heavyIonsPtrIn)
Possibility to pass in pointer for modelling of heavy ion collisions.
Definition: Pythia.h:247
bool setLHAupPtr(LHAupPtr lhaUpPtrIn)
Possibility to pass in pointer to external LHA-interfaced generator.
Definition: BeamSetup.h:59
Definition: StandardModel.h:142
bool forceHadronLevel(bool findJunctions=true)
Generate only the hadronization/decay stage.
Definition: Pythia.cc:1490
bool readFile(string fileName, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in updates for settings or particle data from user-defined file.
Definition: Pythia.h:104
HadronWidths hadronWidths
HadronWidths: the hadron widths data table/database.
Definition: Pythia.h:420
double getSigmaTotal()
Get total cross section for two hadrons in the event record or standalone.
Definition: Pythia.h:319
double m2(const Vec4 &v1)
The squared invariant mass of one or more four-vectors.
Definition: Basics.cc:597
bool setUserHooksPtr(UserHooksPtr userHooksPtrIn)
Possibility to pass in pointer for user hooks.
Definition: Pythia.h:155
bool setPDFBPtr(PDFPtr pdfBPtrIn)
Routine to pass in pointers to PDF&#39;s. Usage optional.
Definition: BeamSetup.cc:129
bool addUserHooksPtr(UserHooksPtr userHooksPtrIn)
Possibility to add further pointers to allow multiple user hooks.
Definition: Pythia.h:159
bool setShowerModelPtr(ShowerModelPtr showerModelPtrIn)
Possibility to pass in pointer for external showers.
Definition: Pythia.h:227
The PartonSystems class describes the whole set of subcollisions.
Definition: PartonSystems.h:42
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The Pythia class contains the top-level routines to generate an event.
Definition: Pythia.h:72
double getSigmaPartial(int procTypeIn)
Definition: Pythia.h:335
UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
Definition: UserHooks.h:319
const Info & info
Public information and statistic on the generation.
Definition: Pythia.h:380
Pythia(string xmlDir="../share/Pythia8/xmldoc", bool printBanner=true)
Constructor. (See Pythia.cc file.)
Definition: Pythia.cc:33
bool next()
Generate the next event.
Definition: Pythia.h:276
BeamShapePtr getBeamShapePtr()
Possibility to access the pointer to the BeamShape object.
Definition: BeamSetup.h:81
double sigmaTotal(int id1, int id2, double eCM12, double m1, double m2, int mixLoHi=0)
Get total cross section.
Definition: SigmaLowEnergy.cc:1526
bool setPhotonFluxPtr(PDFPtr photonFluxAIn, PDFPtr photonFluxBIn)
Set photon fluxes externally. Used with option "PDF:lepton2gammaSet = 2".
Definition: Pythia.h:134
int forceTimeShower(int iBeg, int iEnd, double pTmax, int nBranchMax=0)
Generate only a single timelike shower as in a decay.
Definition: Pythia.h:290
bool insertSigmaPtr(int idx, SigmaProcessPtr sigmaPtrIn, PhaseSpacePtr phaseSpacePtrIn=nullptr)
Definition: Pythia.h:204
void LHAeventList()
List the current Les Houches event.
Definition: Pythia.h:357
MergingPtr mergingPtr
Merging object as wrapper for matrix element merging routines.
Definition: Pythia.h:407
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:423
bool setSigmaPtr(SigmaProcessPtr sigmaPtrIn, PhaseSpacePtr phaseSpacePtrIn=nullptr)
Definition: Pythia.h:190
Logger logger
Logger: for diagnostic messages, errors, statistics, etc.
Definition: Pythia.h:383
bool flag(string key)
Read in settings values: shorthand, not new functionality.
Definition: Pythia.h:368
bool setBeamShapePtr(BeamShapePtr beamShapePtrIn)
Possibility to pass in pointer for beam shape.
Definition: BeamSetup.h:77
Definition: Basics.h:32
bool setHIHooks(HIUserHooksPtr hiHooksPtrIn)
Definition: Pythia.h:252
Definition: SusyCouplings.h:27
Definition: Settings.h:196
void clear()
Reset system list to empty.
Definition: PartonSystems.h:50