1 // Pythia.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2010 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"
20 #include "LesHouches.h"
21 #include "PartonLevel.h"
22 #include "ParticleData.h"
23 #include "PartonDistributions.h"
24 #include "PartonSystems.h"
25 #include "ProcessLevel.h"
26 #include "PythiaStdlib.h"
27 #include "ResonanceWidths.h"
29 #include "SigmaTotal.h"
30 #include "SpaceShower.h"
31 #include "StandardModel.h"
32 #include "SusyLesHouches.h"
33 #include "TimeShower.h"
34 #include "UserHooks.h"
38 //==========================================================================
40 // The Pythia class contains the top-level routines to generate an event.
46 // Constructor. (See Pythia.cc file.)
47 Pythia(string xmlDir = "../xmldoc");
49 // Destructor. (See Pythia.cc file.)
52 // Read in one update for a setting or particle data from a single line.
53 bool readString(string, bool warn = true);
55 // Read in updates for settings or particle data from user-defined file.
56 bool readFile(string fileName, bool warn = true,
57 int subrun = SUBRUNDEFAULT);
58 bool readFile(string fileName, int subrun) {
59 return readFile(fileName, true, subrun);}
60 bool readFile(istream& is = cin, bool warn = true,
61 int subrun = SUBRUNDEFAULT);
62 bool readFile(istream& is, int subrun) {
63 return readFile(is, true, subrun);}
65 // Possibility to pass in pointers to PDF's.
66 bool setPDFPtr( PDF* pdfAPtrIn, PDF* pdfBPtrIn, PDF* pdfHardAPtrIn = 0,
67 PDF* pdfHardBPtrIn = 0, PDF* pdfPomAPtrIn = 0, PDF* pdfPomBPtrIn = 0);
69 // Possibility to pass in pointer for external handling of some decays.
70 bool setDecayPtr( DecayHandler* decayHandlePtrIn,
71 vector<int> handledParticlesIn) {decayHandlePtr = decayHandlePtrIn;
72 handledParticles.resize(0);
73 for(int i = 0; i < int(handledParticlesIn.size()); ++i)
74 handledParticles.push_back( handledParticlesIn[i] ); return true;}
76 // Possibility to pass in pointer for external random number generation.
77 bool setRndmEnginePtr( RndmEngine* rndmEnginePtrIn)
78 { return rndm.rndmEnginePtr( rndmEnginePtrIn);}
80 // Possibility to pass in pointer for user hooks.
81 bool setUserHooksPtr( UserHooks* userHooksPtrIn)
82 { userHooksPtr = userHooksPtrIn; return true;}
84 // Possibility to pass in pointer for beam shape.
85 bool setBeamShapePtr( BeamShape* beamShapePtrIn)
86 { beamShapePtr = beamShapePtrIn; return true;}
88 // Possibility to pass in pointer(s) for external cross section.
89 bool setSigmaPtr( SigmaProcess* sigmaPtrIn)
90 { sigmaPtrs.push_back( sigmaPtrIn); return true;}
92 // Possibility to pass in pointer(s) for external resonance.
93 bool setResonancePtr( ResonanceWidths* resonancePtrIn)
94 { resonancePtrs.push_back( resonancePtrIn); return true;}
96 // Possibility to pass in pointer for external showers.
97 bool setShowerPtr( TimeShower* timesDecPtrIn,
98 TimeShower* timesPtrIn = 0, SpaceShower* spacePtrIn = 0)
99 { timesDecPtr = timesDecPtrIn; timesPtr = timesPtrIn;
100 spacePtr = spacePtrIn; return true;}
102 // Initialization in the CM frame.
103 bool init( int idAin, int idBin, double eCMin);
105 // Initialization with two collinear beams, including fixed target.
106 bool init( int idAin, int idBin, double eAin, double eBin);
108 // Initialization with two acollinear beams.
109 bool init( int idAin, int idBin, double pxAin, double pyAin,
110 double pzAin, double pxBin, double pyBin, double pzBin);
112 // Initialization by a Les Houches Event File.
113 bool init( string LesHouchesEventFile, bool skipInit = false);
115 // Initialization using the Main beam variables.
118 // Initialization according to the Les Houches Accord.
119 bool init( LHAup* lhaUpPtrIn);
121 // Initialization of SLHA data
124 // Generate the next event.
127 // Generate only the hadronization/decay stage.
128 bool forceHadronLevel();
130 // Special routine to allow more decays if on/off switches changed.
131 bool moreDecays() {return hadronLevel.moreDecays(event);}
133 // List the current Les Houches event.
134 void LHAeventList(ostream& os = cout) {
135 if (lhaUpPtr > 0) lhaUpPtr->listEvent(os);}
137 // Skip a number of Les Houches events at input.
138 bool LHAeventSkip(int nSkip) {
139 if (lhaUpPtr > 0) return lhaUpPtr->skipEvent(nSkip); return false;}
141 // Main routine to provide final statistics on generation.
142 void statistics(bool all = false, bool reset = false);
144 // Read in settings values: shorthand, not new functionality.
145 bool flag(string key) {return settings.flag(key);}
146 int mode(string key) {return settings.mode(key);}
147 double parm(string key) {return settings.parm(key);}
148 string word(string key) {return settings.word(key);}
150 // The event record for the parton-level central process.
153 // The event record for the complete event history.
156 // Information on the generation: current subprocess and error statistics.
159 // Settings: databases of flags/modes/parms/words to control run.
162 // ParticleData: the particle data table/database.
163 ParticleData particleData;
165 // Random number generator.
168 // Standard Model couplings, including alphaS and alphaEM.
169 Couplings couplings, *couplingsPtr;
172 // SusyLesHouches - SLHA object for interface to SUSY spectra.
175 // The partonic content of each subcollision system (auxiliary to event).
176 PartonSystems partonSystems;
180 // Constants: could only be changed in the code itself.
181 static const int NTRY, SUBRUNDEFAULT;
183 // Initialization data, extracted from database.
185 bool doProcessLevel, doPartonLevel, doHadronLevel, checkEvent,
188 double epTolErr, epTolWarn;
190 // Initialization data, extracted from init(...) call.
191 bool isConstructed, isInit;
192 int idA, idB, frameType;
193 double mA, mB, pxA, pxB, pyA, pyB, pzA, pzB, eA, eB,
194 pzAcm, pzBcm, eCM, betaZ, gammaZ;
195 Vec4 pAinit, pBinit, pAnow, pBnow;
196 RotBstMatrix MfromCM, MtoCM;
198 // information for error checkout.
200 vector<int> iErrId, iErrCol, iErrNan, iErrNanVtx;
202 // Pointers to the parton distributions of the two incoming beams.
206 // Extra PDF pointers to be used in hard processes only.
210 // Extra Pomeron PDF pointers to be used in diffractive processes only.
214 // Keep track when "new" has been used and needs a "delete" for PDF's.
215 bool useNewPdfA, useNewPdfB, useNewPdfHard, useNewPdfPomA, useNewPdfPomB;
217 // The two incoming beams.
221 // Alternative Pomeron beam-inside-beam.
222 BeamParticle beamPomA;
223 BeamParticle beamPomB;
225 // LHAup object for generating external events.
226 bool doLHA, useNewLHA;
229 // Pointer to external decay handler and list of particles it handles.
230 DecayHandler* decayHandlePtr;
231 vector<int> handledParticles;
233 // Pointer to UserHooks object for user interaction with program.
234 UserHooks* userHooksPtr;
235 bool hasUserHooks, doVetoProcess, doVetoPartons;
237 // Pointer to BeamShape object for beam momentum and interaction vertex.
238 BeamShape* beamShapePtr;
239 bool useNewBeamShape, doMomentumSpread, doVertexSpread;
241 // Pointers to external processes derived from the Pythia base classes.
242 vector<SigmaProcess*> sigmaPtrs;
244 // Pointers to external calculation of resonance widths.
245 vector<ResonanceWidths*> resonancePtrs;
247 // Pointers to timelike and spacelike showers.
248 TimeShower* timesDecPtr;
249 TimeShower* timesPtr;
250 SpaceShower* spacePtr;
251 bool useNewTimes, useNewSpace;
253 // The main generator class to define the core process of the event.
254 ProcessLevel processLevel;
256 // The main generator class to produce the parton level of the event.
257 PartonLevel partonLevel;
259 // The main generator class to produce the hadron level of the event.
260 HadronLevel hadronLevel;
262 // The total cross section class is used both on process and parton level.
265 // Write the Pythia banner, with symbol and version information.
266 void banner(ostream& os = cout);
268 // Check for lines in file that mark the beginning of new subrun.
269 int readSubrun(string line, bool warn = true, ostream& os = cout);
271 // Initialization routine to set up the whole generation machinery.
274 // Check that combinations of settings are allowed; change if not.
275 void checkSettings();
277 // Check that beams and beam combination can be handled.
280 // Calculate kinematics at initialization.
281 bool initKinematics();
283 // Set up pointers to PDFs.
286 // Recalculate kinematics for each event when beam momentum has a spread.
287 void nextKinematics();
289 // Boost from CM frame to lab frame, or inverse. Set production vertex.
290 void boostAndVertex(bool toLab, bool setVertex);
292 // Check that the final event makes sense.
293 bool check(ostream& os = cout);
295 // Auxiliary to set parton densities among list of possibilities.
296 PDF* getPDFPtr(int idIn, int sequence = 1);
300 //==========================================================================
302 } // end namespace Pythia8
304 #endif // Pythia8_Pythia_H