1 // Pythia.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2013 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
6 // This file contains the main class for event generation.
7 // Pythia: provide the main user interface to everything else.
9 #ifndef Pythia8_Pythia_H
10 #define Pythia8_Pythia_H
14 #include "BeamParticle.h"
15 #include "BeamShape.h"
17 #include "FragmentationFlavZpT.h"
18 #include "HadronLevel.h"
21 #include "LesHouches.h"
22 #include "PartonLevel.h"
23 #include "ParticleData.h"
24 #include "PartonDistributions.h"
25 #include "PartonSystems.h"
26 #include "ProcessLevel.h"
27 #include "PythiaStdlib.h"
28 #include "ResonanceWidths.h"
31 #include "SigmaTotal.h"
32 #include "SpaceShower.h"
33 #include "StandardModel.h"
34 #include "SusyLesHouches.h"
35 #include "TimeShower.h"
36 #include "UserHooks.h"
37 #include "MergingHooks.h"
42 //==========================================================================
44 // The Pythia class contains the top-level routines to generate an event.
50 // Constructor. (See Pythia.cc file.)
51 Pythia(string xmlDir = "../xmldoc", bool printBanner = true);
53 // Destructor. (See Pythia.cc file.)
56 // Read in one update for a setting or particle data from a single line.
57 bool readString(string, bool warn = true);
59 // Read in updates for settings or particle data from user-defined file.
60 bool readFile(string fileName, bool warn = true,
61 int subrun = SUBRUNDEFAULT);
62 bool readFile(string fileName, int subrun) {
63 return readFile(fileName, true, subrun);}
64 bool readFile(istream& is = cin, bool warn = true,
65 int subrun = SUBRUNDEFAULT);
66 bool readFile(istream& is, int subrun) {
67 return readFile(is, true, subrun);}
69 // Possibility to pass in pointers to PDF's.
70 bool setPDFPtr( PDF* pdfAPtrIn, PDF* pdfBPtrIn, PDF* pdfHardAPtrIn = 0,
71 PDF* pdfHardBPtrIn = 0, PDF* pdfPomAPtrIn = 0, PDF* pdfPomBPtrIn = 0);
73 // Possibility to pass in pointer to external LHA-interfaced generator.
74 bool setLHAupPtr( LHAup* lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn; return true;}
76 // Possibility to pass in pointer for external handling of some decays.
77 bool setDecayPtr( DecayHandler* decayHandlePtrIn,
78 vector<int> handledParticlesIn) {decayHandlePtr = decayHandlePtrIn;
79 handledParticles.resize(0);
80 for(int i = 0; i < int(handledParticlesIn.size()); ++i)
81 handledParticles.push_back( handledParticlesIn[i] ); return true;}
83 // Possibility to pass in pointer for external random number generation.
84 bool setRndmEnginePtr( RndmEngine* rndmEnginePtrIn)
85 { return rndm.rndmEnginePtr( rndmEnginePtrIn);}
87 // Possibility to pass in pointer for user hooks.
88 bool setUserHooksPtr( UserHooks* userHooksPtrIn)
89 { userHooksPtr = userHooksPtrIn; return true;}
91 // Possibility to pass in pointer for merging hooks.
92 bool setMergingHooksPtr( MergingHooks* mergingHooksPtrIn)
93 { mergingHooksPtr = mergingHooksPtrIn; return true;}
95 // Possibility to pass in pointer for beam shape.
96 bool setBeamShapePtr( BeamShape* beamShapePtrIn)
97 { beamShapePtr = beamShapePtrIn; return true;}
99 // Possibility to pass in pointer(s) for external cross section.
100 bool setSigmaPtr( SigmaProcess* sigmaPtrIn)
101 { sigmaPtrs.push_back( sigmaPtrIn); return true;}
103 // Possibility to pass in pointer(s) for external resonance.
104 bool setResonancePtr( ResonanceWidths* resonancePtrIn)
105 { resonancePtrs.push_back( resonancePtrIn); return true;}
107 // Possibility to pass in pointer for external showers.
108 bool setShowerPtr( TimeShower* timesDecPtrIn,
109 TimeShower* timesPtrIn = 0, SpaceShower* spacePtrIn = 0)
110 { timesDecPtr = timesDecPtrIn; timesPtr = timesPtrIn;
111 spacePtr = spacePtrIn; return true;}
113 // Initialization using the Beams variables.
116 // Deprecated: initialization in the CM frame.
117 bool init( int idAin, int idBin, double eCMin);
119 // Deprecated: initialization with two collinear beams, e.g. fixed target.
120 bool init( int idAin, int idBin, double eAin, double eBin);
122 // Deprecated: initialization with two acollinear beams.
123 bool init( int idAin, int idBin, double pxAin, double pyAin,
124 double pzAin, double pxBin, double pyBin, double pzBin);
126 // Deprecated: initialization by a Les Houches Event File.
127 bool init( string LesHouchesEventFile, bool skipInit = false);
129 // Deprecated: initialization according to the Les Houches Accord.
130 bool init( LHAup* lhaUpPtrIn);
132 // Generate the next event.
135 // Generate only a single timelike shower as in a decay.
136 int forceTimeShower( int iBeg, int iEnd, double pTmax, int nBranchMax = 0)
137 { return timesDecPtr->shower( iBeg, iEnd, event, pTmax, nBranchMax); }
139 // Generate only the hadronization/decay stage.
140 bool forceHadronLevel( bool findJunctions = true);
142 // Special routine to allow more decays if on/off switches changed.
143 bool moreDecays() {return hadronLevel.moreDecays(event);}
145 // Special routine to force R-hadron decay when not done before.
146 bool forceRHadronDecays() {return doRHadronDecays();}
148 // List the current Les Houches event.
149 void LHAeventList(ostream& os = cout) {
150 if (lhaUpPtr != 0) lhaUpPtr->listEvent(os);}
152 // Skip a number of Les Houches events at input.
153 bool LHAeventSkip(int nSkip) {
154 if (lhaUpPtr != 0) return lhaUpPtr->skipEvent(nSkip); return false;}
156 // Main routine to provide final statistics on generation.
159 // Deprecated: alternative to provide final statistics on generation.
160 void statistics(bool all = false, bool reset = false);
162 // Read in settings values: shorthand, not new functionality.
163 bool flag(string key) {return settings.flag(key);}
164 int mode(string key) {return settings.mode(key);}
165 double parm(string key) {return settings.parm(key);}
166 string word(string key) {return settings.word(key);}
168 // Auxiliary to set parton densities among list of possibilities.
169 PDF* getPDFPtr(int idIn, int sequence = 1);
171 // The event record for the parton-level central process.
174 // The event record for the complete event history.
177 // Information on the generation: current subprocess and error statistics.
180 // Settings: databases of flags/modes/parms/words to control run.
183 // ParticleData: the particle data table/database.
184 ParticleData particleData;
186 // Random number generator.
189 // Standard Model couplings, including alphaS and alphaEM.
192 Couplings* couplingsPtr;
194 // SusyLesHouches - SLHA object for interface to SUSY spectra.
197 // The partonic content of each subcollision system (auxiliary to event).
198 PartonSystems partonSystems;
200 // Merging object as wrapper for matrix element merging routines.
203 // Pointer to MergingHooks object for user interaction with the merging.
204 // MergingHooks also more generally steers the matrix element merging.
205 MergingHooks* mergingHooksPtr;
209 // Copy and = constructors are made private so they cannot be used.
210 Pythia(const Pythia&);
211 Pythia& operator=(const Pythia&);
213 // Constants: could only be changed in the code itself.
214 static const double VERSIONNUMBERCODE;
215 static const int NTRY, SUBRUNDEFAULT;
217 // Initialization data, extracted from database.
219 bool doProcessLevel, doPartonLevel, doHadronLevel, doDiffraction,
220 doResDec, doFSRinRes, decayRHadrons, abortIfVeto, checkEvent,
223 double epTolErr, epTolWarn;
225 // Initialization data, extracted from init(...) call.
226 bool isConstructed, isInit, isUnresolvedA, isUnresolvedB, showSaV,
228 int idA, idB, frameType, boostType, nCount, nShowLHA, nShowInfo,
230 double mA, mB, pxA, pxB, pyA, pyB, pzA, pzB, eA, eB,
231 pzAcm, pzBcm, eCM, betaZ, gammaZ;
232 Vec4 pAinit, pBinit, pAnow, pBnow;
233 RotBstMatrix MfromCM, MtoCM;
235 // information for error checkout.
237 vector<int> iErrId, iErrCol, iErrNan, iErrNanVtx;
239 // Pointers to the parton distributions of the two incoming beams.
243 // Extra PDF pointers to be used in hard processes only.
247 // Extra Pomeron PDF pointers to be used in diffractive processes only.
251 // Keep track when "new" has been used and needs a "delete" for PDF's.
252 bool useNewPdfA, useNewPdfB, useNewPdfHard, useNewPdfPomA, useNewPdfPomB;
254 // The two incoming beams.
258 // Alternative Pomeron beam-inside-beam.
259 BeamParticle beamPomA;
260 BeamParticle beamPomB;
262 // LHAup object for generating external events.
263 bool doLHA, useNewLHA;
266 // Pointer to external decay handler and list of particles it handles.
267 DecayHandler* decayHandlePtr;
268 vector<int> handledParticles;
270 // Pointer to UserHooks object for user interaction with program.
271 UserHooks* userHooksPtr;
272 bool hasUserHooks, doVetoProcess, doVetoPartons, retryPartonLevel;
274 // Pointer to BeamShape object for beam momentum and interaction vertex.
275 BeamShape* beamShapePtr;
276 bool useNewBeamShape, doMomentumSpread, doVertexSpread;
278 // Pointers to external processes derived from the Pythia base classes.
279 vector<SigmaProcess*> sigmaPtrs;
281 // Pointers to external calculation of resonance widths.
282 vector<ResonanceWidths*> resonancePtrs;
284 // Pointers to timelike and spacelike showers.
285 TimeShower* timesDecPtr;
286 TimeShower* timesPtr;
287 SpaceShower* spacePtr;
288 bool useNewTimes, useNewSpace;
290 // The main generator class to define the core process of the event.
291 ProcessLevel processLevel;
293 // The main generator class to produce the parton level of the event.
294 PartonLevel partonLevel;
296 // The main generator class to perform trial showers of the event.
297 PartonLevel trialPartonLevel;
299 // Flags for defining the merging scheme.
300 bool hasMergingHooks, hasOwnMergingHooks, doMerging;
302 // The main generator class to produce the hadron level of the event.
303 HadronLevel hadronLevel;
305 // The total cross section class is used both on process and parton level.
308 // The RHadrons class is used both at PartonLevel and HadronLevel.
311 // Write the Pythia banner, with symbol and version information.
312 void banner(ostream& os = cout);
314 // Check for lines in file that mark the beginning of new subrun.
315 int readSubrun(string line, bool warn = true, ostream& os = cout);
317 // Check that combinations of settings are allowed; change if not.
318 void checkSettings();
320 // Check that beams and beam combination can be handled.
323 // Calculate kinematics at initialization.
324 bool initKinematics();
326 // Set up pointers to PDFs.
329 // Recalculate kinematics for each event when beam momentum has a spread.
330 void nextKinematics();
332 // Boost from CM frame to lab frame, or inverse. Set production vertex.
333 void boostAndVertex(bool toLab, bool setVertex);
335 // Perform R-hadron decays.
336 bool doRHadronDecays();
338 // Check that the final event makes sense.
339 bool check(ostream& os = cout);
341 // Initialization of SLHA data.
346 //==========================================================================
348 } // end namespace Pythia8
350 #endif // Pythia8_Pythia_H