1 // Info.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 a class that keep track of generic event info.
7 // Info: contains information on the generation process and errors.
10 #define Pythia8_Info_H
12 #include "PythiaStdlib.h"
16 //==========================================================================
18 // The Info class contains a mixed bag of information on the event
19 // generation activity, especially on the current subprocess properties,
20 // and on the number of errors encountered. This is used by the
21 // generation machinery, but can also be read by the user.
22 // Note: some methods that maybe should not be accessible to the user
23 // are still public, to work also for user-written FSR/ISR classes.
30 Info() : eCMSave(0.), lowPTmin(false), a0MPISave(0.), weightCKKWLSave(1.),
32 for (int i = 0; i < 40; ++i) counters[i] = 0;}
34 // Listing of most available information on current event.
35 void list(ostream& os = cout) const;
37 // Beam particles (in rest frame). CM energy of event.
38 int idA() const {return idASave;}
39 int idB() const {return idBSave;}
40 double pzA() const {return pzASave;}
41 double pzB() const {return pzBSave;}
42 double eA() const {return eASave;}
43 double eB() const {return eBSave;}
44 double mA() const {return mASave;}
45 double mB() const {return mBSave;}
46 double eCM() const {return eCMSave;}
47 double s() const {return sSave;}
49 // Warnings from initialization.
50 bool tooLowPTmin() const {return lowPTmin;}
52 // Process name and code, and the number of final-state particles.
53 string name() const {return nameSave;}
54 int code() const {return codeSave;}
55 int nFinal() const {return nFinalSave;}
57 // Are beam particles resolved, with pdf's? Are they diffractive?
58 bool isResolved() const {return isRes;}
59 bool isDiffractiveA() const {return isDiffA;}
60 bool isDiffractiveB() const {return isDiffB;}
61 bool isDiffractiveC() const {return isDiffC;}
62 bool isMinBias() const {return isMB;}
64 // Information for Les Houches Accord and reading files.
65 bool isLHA() const {return isLH;}
66 bool atEndOfFile() const {return atEOF;}
68 // For minbias and Les Houches Accord identify hardest subprocess.
69 bool hasSub(int i = 0) const {return hasSubSave[i];}
70 string nameSub(int i = 0) const {return nameSubSave[i];}
71 int codeSub(int i = 0) const {return codeSubSave[i];}
72 int nFinalSub(int i = 0) const {return nFinalSubSave[i];}
74 // Incoming parton flavours and x values.
75 int id1(int i = 0) const {return id1Save[i];}
76 int id2(int i = 0) const {return id2Save[i];}
77 double x1(int i = 0) const {return x1Save[i];}
78 double x2(int i = 0) const {return x2Save[i];}
79 double y(int i = 0) const {return 0.5 * log( x1Save[i] / x2Save[i]);}
80 double tau(int i = 0) const {return x1Save[i] * x2Save[i];}
82 // Hard process flavours, x values, parton densities, couplings, Q2 scales.
83 int id1pdf(int i = 0) const {return id1pdfSave[i];}
84 int id2pdf(int i = 0) const {return id2pdfSave[i];}
85 double x1pdf(int i = 0) const {return x1pdfSave[i];}
86 double x2pdf(int i = 0) const {return x2pdfSave[i];}
87 double pdf1(int i = 0) const {return pdf1Save[i];}
88 double pdf2(int i = 0) const {return pdf2Save[i];}
89 double QFac(int i = 0) const {return sqrtpos(Q2FacSave[i]);}
90 double Q2Fac(int i = 0) const {return Q2FacSave[i];}
91 bool isValence1() const {return isVal1;}
92 bool isValence2() const {return isVal2;}
93 double alphaS(int i = 0) const {return alphaSSave[i];}
94 double alphaEM(int i = 0) const {return alphaEMSave[i];}
95 double QRen(int i = 0) const {return sqrtpos(Q2RenSave[i]);}
96 double Q2Ren(int i = 0) const {return Q2RenSave[i];}
98 // Mandelstam variables (notation as if subcollision).
99 double mHat(int i = 0) const {return sqrt(sH[i]);}
100 double sHat(int i = 0) const {return sH[i];}
101 double tHat(int i = 0) const {return tH[i];}
102 double uHat(int i = 0) const {return uH[i];}
103 double pTHat(int i = 0) const {return pTH[i];}
104 double pT2Hat(int i = 0) const {return pTH[i] * pTH[i];}
105 double m3Hat(int i = 0) const {return m3H[i];}
106 double m4Hat(int i = 0) const {return m4H[i];}
107 double thetaHat(int i = 0) const {return thetaH[i];}
108 double phiHat(int i = 0) const {return phiH[i];}
110 // Weight of current event; normally 1, but used for Les Houches events
111 // or when reweighting phase space selection. Conversion from mb to pb
112 // for LHA strategy +-4. Also cumulative sum.
113 double weight() const;
114 double weightSum() const;
115 double lhaStrategy() const {return lhaStrategySave;}
117 // Number of times other steps have been carried out.
118 int nISR() const {return nISRSave;}
119 int nFSRinProc() const {return nFSRinProcSave;}
120 int nFSRinRes() const {return nFSRinResSave;}
122 // Maximum pT scales for MPI, ISR and FSR (in hard process).
123 double pTmaxMPI() const {return pTmaxMPISave;}
124 double pTmaxISR() const {return pTmaxISRSave;}
125 double pTmaxFSR() const {return pTmaxFSRSave;}
127 // Current evolution scale (for UserHooks).
128 double pTnow() const {return pTnowSave;}
130 // Impact parameter picture, global information
131 double a0MPI() const {return a0MPISave;}
133 // Impact parameter picture, as set by hardest interaction.
134 double bMPI() const {return (bIsSet) ? bMPISave : 1.;}
135 double enhanceMPI() const {return (bIsSet) ? enhanceMPISave : 1.;}
136 double eMPI(int i) const {return (bIsSet) ? eMPISave[i] : 1.;}
138 // Number of multiparton interactions, with code and pT for them.
139 int nMPI() const {return nMPISave;}
140 int codeMPI(int i) const {return codeMPISave[i];}
141 double pTMPI(int i) const {return pTMPISave[i];}
142 int iAMPI(int i) const {return iAMPISave[i];}
143 int iBMPI(int i) const {return iBMPISave[i];}
145 // Cross section estimate.
146 long nTried(int i = 0) {return (i == 0) ? nTry : nTryM[i];}
147 long nSelected(int i = 0) {return (i == 0) ? nSel : nSelM[i];}
148 long nAccepted(int i = 0) {return (i == 0) ? nAcc : nAccM[i];}
149 double sigmaGen(int i = 0) {return (i == 0) ? sigGen : sigGenM[i];}
150 double sigmaErr(int i = 0) {return (i == 0) ? sigErr : sigErrM[i];}
152 // Counters for number of loops in various places.
153 int getCounter( int i) const {return counters[i];}
155 // Set or increase the value stored in a counter.
156 void setCounter( int i, int value = 0) {counters[i] = value;}
157 void addCounter( int i, int value = 1) {counters[i] += value;}
159 // Reset to empty map of error messages.
160 void errorReset() {messages.clear();}
162 // Print a message the first few times. Insert in database.
163 void errorMsg(string messageIn, string extraIn = " ",
164 bool showAlways = false, ostream& os = cout);
166 // Provide total number of errors/aborts/warnings experienced to date.
167 int errorTotalNumber();
169 // Print statistics on errors/aborts/warnings.
170 void errorStatistics(ostream& os = cout);
172 // Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
173 void setTooLowPTmin(bool lowPTminIn) {lowPTmin = lowPTminIn;}
175 // Set info on valence character of hard collision partons.
176 void setValence( bool isVal1In, bool isVal2In) {isVal1 = isVal1In;
179 // Set and get some MPI/ISR/FSR properties needed for matching,
180 // i.e. mainly of internal relevance.
181 void hasHistory( bool hasHistoryIn) {hasHistorySave = hasHistoryIn;}
182 bool hasHistory() {return hasHistorySave;}
183 void zNowISR( double zNowIn) {zNowISRSave = zNowIn;}
184 double zNowISR() {return zNowISRSave;}
185 void pT2NowISR( double pT2NowIn) {pT2NowISRSave = pT2NowIn;}
186 double pT2NowISR() {return pT2NowISRSave;}
188 // Return CKKW-L weight.
189 double getWeightCKKWL() const { return weightCKKWLSave;}
190 // Set CKKW-L weight.
191 void setWeightCKKWL( double weightIn) { weightCKKWLSave = weightIn;}
192 // Return merging weight.
193 double mergingWeight() const { return weightCKKWLSave;}
195 // Return the complete NLO weight.
196 double mergingWeightNLO() const {
197 return (weightCKKWLSave - weightFIRSTSave); }
198 // Return the O(\alpha_s)-term of the CKKW-L weight.
199 double getWeightFIRST() const { return weightFIRSTSave;}
200 // Set the O(\alpha_s)-term of the CKKW-L weight.
201 void setWeightFIRST( double weightIn) { weightFIRSTSave = weightIn;}
203 // Return an LHEF header
204 string header(const string &key) {
205 if (headers.find(key) == headers.end()) return "";
206 else return headers[key];
209 // Return a list of all header key names
210 vector < string > headerKeys() {
211 vector < string > keys;
212 for (map < string, string >::iterator it = headers.begin();
213 it != headers.end(); it++)
214 keys.push_back(it->first);
220 // Number of times the same error message is repeated, unless overridden.
221 static const int TIMESTOPRINT;
223 // Allow conversion from mb to pb.
224 static const double CONVERTMB2PB;
226 // Store common beam quantities.
227 int idASave, idBSave;
228 double pzASave, eASave,mASave, pzBSave, eBSave, mBSave, eCMSave, sSave;
230 // Store initialization information.
233 // Store common integrated cross section quantities.
234 long nTry, nSel, nAcc;
235 double sigGen, sigErr, wtAccSum;
236 map<int, long> nTryM, nSelM, nAccM;
237 map<int, double> sigGenM, sigErrM;
240 // Store common MPI information.
243 // Store current-event quantities.
244 bool isRes, isDiffA, isDiffB, isDiffC, isMB, isLH, hasSubSave[4],
245 bIsSet, evolIsSet, atEOF, isVal1, isVal2, hasHistorySave;
246 int codeSave, codeSubSave[4], nFinalSave, nFinalSubSave[4], nTotal,
247 id1Save[4], id2Save[4], id1pdfSave[4], id2pdfSave[4], nMPISave,
248 nISRSave, nFSRinProcSave, nFSRinResSave;
249 double x1Save[4], x2Save[4], x1pdfSave[4], x2pdfSave[4], pdf1Save[4],
250 pdf2Save[4], Q2FacSave[4], alphaEMSave[4], alphaSSave[4],
251 Q2RenSave[4], sH[4], tH[4], uH[4], pTH[4], m3H[4], m4H[4],
252 thetaH[4], phiH[4], weightSave, bMPISave, enhanceMPISave,
253 pTmaxMPISave, pTmaxISRSave, pTmaxFSRSave, pTnowSave,
254 zNowISRSave, pT2NowISRSave;
255 string nameSave, nameSubSave[4];
256 vector<int> codeMPISave, iAMPISave, iBMPISave;
257 vector<double> pTMPISave, eMPISave;
259 // Vector of various loop counters.
262 // Map for all error messages.
263 map<string, int> messages;
265 // Map for LHEF headers
266 map<string, string> headers;
268 // Friend classes allowed to set info.
270 friend class ProcessLevel;
271 friend class ProcessContainer;
272 friend class PartonLevel;
273 friend class MultipartonInteractions;
276 // Set info on the two incoming beams: only from Pythia class.
277 void setBeamA( int idAin, double pzAin, double eAin, double mAin) {
278 idASave = idAin; pzASave = pzAin; eASave = eAin; mASave = mAin;}
279 void setBeamB( int idBin, double pzBin, double eBin, double mBin) {
280 idBSave = idBin; pzBSave = pzBin; eBSave = eBin; mBSave = mBin;}
281 void setECM( double eCMin) {eCMSave = eCMin; sSave = eCMSave * eCMSave;}
283 // Reset info for current event: only from Pythia class.
285 isRes = isDiffA = isDiffB = isDiffC = isMB = isLH = bIsSet
286 = evolIsSet = atEOF = isVal1 = isVal2 = hasHistorySave = false;
287 codeSave = nFinalSave = nTotal = nMPISave = nISRSave = nFSRinProcSave
289 weightSave = bMPISave = enhanceMPISave = weightCKKWLSave = 1.;
290 pTmaxMPISave = pTmaxISRSave = pTmaxFSRSave = pTnowSave = zNowISRSave
291 = pT2NowISRSave = weightFIRSTSave = 0.;
293 for (int i = 0; i < 4; ++i) {
294 hasSubSave[i] = false;
295 codeSubSave[i] = nFinalSubSave[i] = id1pdfSave[i] = id2pdfSave[i]
296 = id1Save[i] = id2Save[i] = 0;
297 x1pdfSave[i] = x2pdfSave[i] = pdf1Save[i] = pdf2Save[i]
298 = Q2FacSave[i] = alphaEMSave[i] = alphaSSave[i] = Q2RenSave[i]
299 = x1Save[i] = x2Save[i] = sH[i] = tH[i] = uH[i] = pTH[i]
300 = m3H[i] = m4H[i] = thetaH[i] = phiH[i] = 0.;
301 nameSubSave[i] = " ";
303 codeMPISave.resize(0); iAMPISave.resize(0); iBMPISave.resize(0);
304 pTMPISave.resize(0); eMPISave.resize(0); }
306 // Reset info arrays only for parton and hadron level.
307 int sizeMPIarrays() const { return codeMPISave.size(); }
308 void resizeMPIarrays(int newSize) { codeMPISave.resize(newSize);
309 iAMPISave.resize(newSize); iBMPISave.resize(newSize);
310 pTMPISave.resize(newSize); eMPISave.resize(newSize); }
312 // Set info on the (sub)process: from ProcessLevel, ProcessContainer or
313 // MultipartonInteractions classes.
314 void setType( string nameIn, int codeIn, int nFinalIn,
315 bool isMinBiasIn = false, bool isResolvedIn = true,
316 bool isDiffractiveAin = false, bool isDiffractiveBin = false,
317 bool isDiffractiveCin = false, bool isLHAin = false) {
318 nameSave = nameIn; codeSave = codeIn; nFinalSave = nFinalIn;
319 isMB = isMinBiasIn; isRes = isResolvedIn; isDiffA = isDiffractiveAin;
320 isDiffB = isDiffractiveBin; isDiffC = isDiffractiveCin; isLH = isLHAin;
321 nTotal = 2 + nFinalSave; bIsSet = false; hasSubSave[0] = false;
322 nameSubSave[0] = " "; codeSubSave[0] = 0; nFinalSubSave[0] = 0;
324 void setSubType( int iDS, string nameSubIn, int codeSubIn,
325 int nFinalSubIn) { hasSubSave[iDS] = true; nameSubSave[iDS] = nameSubIn;
326 codeSubSave[iDS] = codeSubIn; nFinalSubSave[iDS] = nFinalSubIn;}
327 void setPDFalpha( int iDS, int id1pdfIn, int id2pdfIn, double x1pdfIn,
328 double x2pdfIn, double pdf1In, double pdf2In, double Q2FacIn,
329 double alphaEMIn, double alphaSIn, double Q2RenIn) {
330 id1pdfSave[iDS] = id1pdfIn; id2pdfSave[iDS] = id2pdfIn;
331 x1pdfSave[iDS] = x1pdfIn; x2pdfSave[iDS] = x2pdfIn;
332 pdf1Save[iDS] = pdf1In; pdf2Save[iDS] = pdf2In; Q2FacSave[iDS] = Q2FacIn;
333 alphaEMSave[iDS] = alphaEMIn; alphaSSave[iDS] = alphaSIn;
334 Q2RenSave[iDS] = Q2RenIn;}
335 void setKin( int iDS, int id1In, int id2In, double x1In, double x2In,
336 double sHatIn, double tHatIn, double uHatIn, double pTHatIn,
337 double m3HatIn, double m4HatIn, double thetaHatIn, double phiHatIn) {
338 id1Save[iDS] = id1In; id2Save[iDS] = id2In; x1Save[iDS] = x1In;
339 x2Save[iDS] = x2In; sH[iDS] = sHatIn; tH[iDS] = tHatIn;
340 uH[iDS] = uHatIn; pTH[iDS] = pTHatIn; m3H[iDS] = m3HatIn;
341 m4H[iDS] = m4HatIn; thetaH[iDS] = thetaHatIn; phiH[iDS] = phiHatIn;}
342 void setTypeMPI( int codeMPIIn, double pTMPIIn, int iAMPIIn = 0,
343 int iBMPIIn = 0, double eMPIIn = 1.) {codeMPISave.push_back(codeMPIIn);
344 pTMPISave.push_back(pTMPIIn); iAMPISave.push_back(iAMPIIn);
345 iBMPISave.push_back(iBMPIIn); eMPISave.push_back(eMPIIn); }
347 // Set info on cross section: from ProcessLevel (reset in Pythia).
348 void sigmaReset() { nTry = nSel = nAcc = 0; sigGen = sigErr = wtAccSum = 0.;
349 nTryM.clear(); nSelM.clear(); nAccM.clear(); sigGenM.clear();
351 void setSigma( int i, long nTryIn, long nSelIn, long nAccIn,
352 double sigGenIn, double sigErrIn, double wtAccSumIn) {
353 if (i == 0) {nTry = nTryIn; nSel = nSelIn; nAcc = nAccIn;
354 sigGen = sigGenIn; sigErr = sigErrIn; wtAccSum = wtAccSumIn;}
355 else { nTryM[i] = nTryIn; nSelM[i] = nSelIn; nAccM[i] = nAccIn;
356 sigGenM[i] = sigGenIn; sigErrM[i] = sigErrIn;} }
358 // Set info on impact parameter: from PartonLevel.
359 void setImpact( double bMPIIn, double enhanceMPIIn, bool bIsSetIn = true) {
360 bMPISave = bMPIIn; enhanceMPISave = eMPISave[0] = enhanceMPIIn,
363 // Set info on pTmax scales and number of evolution steps: from PartonLevel.
364 void setPartEvolved( int nMPIIn, int nISRIn) {
365 nMPISave = nMPIIn; nISRSave = nISRIn;}
366 void setEvolution( double pTmaxMPIIn, double pTmaxISRIn, double pTmaxFSRIn,
367 int nMPIIn, int nISRIn, int nFSRinProcIn, int nFSRinResIn) {
368 pTmaxMPISave = pTmaxMPIIn; pTmaxISRSave = pTmaxISRIn;
369 pTmaxFSRSave = pTmaxFSRIn; nMPISave = nMPIIn; nISRSave = nISRIn;
370 nFSRinProcSave = nFSRinProcIn; nFSRinResSave = nFSRinResIn;
373 // Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
374 void setPTnow( double pTnowIn) {pTnowSave = pTnowIn;}
376 // Set a0 from MultipartonInteractions.
377 void seta0MPI(double a0MPIin) {a0MPISave = a0MPIin;}
379 // Set info whether reading of Les Houches Accord file at end.
380 void setEndOfFile( bool atEOFin) {atEOF = atEOFin;}
382 // Set event weight; currently only for Les Houches description.
383 void setWeight( double weightIn, int lhaStrategyIn) {
384 weightSave = weightIn; lhaStrategySave = lhaStrategyIn; }
386 // Save merging weight (i.e. CKKW-L-type weight, summed O(\alpha_s) weight)
387 double weightCKKWLSave, weightFIRSTSave;
390 void setHeader(const string &key, const string &val)
391 { headers[key] = val; }
395 //==========================================================================
397 } // end namespace Pythia8
399 #endif // Pythia8_Info_H