StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Info.h
1 // Info.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 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.
5 
6 // This file contains a class that keep track of generic event info.
7 // Info: contains information on the generation process and errors.
8 
9 #ifndef Pythia8_Info_H
10 #define Pythia8_Info_H
11 
12 #include "PythiaStdlib.h"
13 
14 namespace Pythia8 {
15 
16 //==========================================================================
17 
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.
24 
25 class Info {
26 
27 public:
28 
29  // Constructor.
30  Info() : lowPTmin(false), a0MPISave(0.), mergingWeightSave(1.) {
31  for (int i = 0; i < 40; ++i) counters[i] = 0;}
32 
33  // Listing of most available information on current event.
34  void list(ostream& os = cout) const;
35 
36  // Beam particles (in rest frame). CM energy of event.
37  int idA() const {return idASave;}
38  int idB() const {return idBSave;}
39  double pzA() const {return pzASave;}
40  double pzB() const {return pzBSave;}
41  double eA() const {return eASave;}
42  double eB() const {return eBSave;}
43  double mA() const {return mASave;}
44  double mB() const {return mBSave;}
45  double eCM() const {return eCMSave;}
46  double s() const {return sSave;}
47 
48  // Warnings from initialization.
49  bool tooLowPTmin() const {return lowPTmin;}
50 
51  // Process name and code, and the number of final-state particles.
52  string name() const {return nameSave;}
53  int code() const {return codeSave;}
54  int nFinal() const {return nFinalSave;}
55 
56  // Are beam particles resolved, with pdf's? Are they diffractive?
57  bool isResolved() const {return isRes;}
58  bool isDiffractiveA() const {return isDiffA;}
59  bool isDiffractiveB() const {return isDiffB;}
60  bool isMinBias() const {return isMB;}
61 
62  // Information for Les Houches Accord and reading files.
63  bool isLHA() const {return isLH;}
64  bool atEndOfFile() const {return atEOF;}
65 
66  // For minbias and Les Houches Accord identify hardest subprocess.
67  bool hasSub() const {return hasSubSave;}
68  string nameSub() const {return nameSubSave;}
69  int codeSub() const {return codeSubSave;}
70  int nFinalSub() const {return nFinalSubSave;}
71 
72  // Incoming parton flavours and x values.
73  int id1() const {return id1Save;}
74  int id2() const {return id2Save;}
75  double x1() const {return x1Save;}
76  double x2() const {return x2Save;}
77  double y() const {return 0.5 * log( x1Save / x2Save );}
78  double tau() const {return x1Save * x2Save;}
79 
80  // Incoming parton densities, hard process couplings, Q2 scales.
81  double pdf1() const {return pdf1Save;}
82  double pdf2() const {return pdf2Save;}
83  double QFac() const {return sqrtpos(Q2FacSave);}
84  double Q2Fac() const {return Q2FacSave;}
85  bool isValence1() const {return isVal1;}
86  bool isValence2() const {return isVal2;}
87  double alphaS() const {return alphaSSave;}
88  double alphaEM() const {return alphaEMSave;}
89  double QRen() const {return sqrtpos(Q2RenSave);}
90  double Q2Ren() const {return Q2RenSave;}
91 
92  // Mandelstam variables (notation as if subcollision).
93  double mHat() const {return sqrt(sH);}
94  double sHat() const {return sH;}
95  double tHat() const {return tH;}
96  double uHat() const {return uH;}
97  double pTHat() const {return pTH;}
98  double pT2Hat() const {return pTH*pTH;}
99  double m3Hat() const {return m3H;}
100  double m4Hat() const {return m4H;}
101  double thetaHat() const {return thetaH;}
102  double phiHat() const {return phiH;}
103 
104  // Weight of current event; normally 1, but used for Les Houches events
105  // or when reweighting phase space selection. Conversion from mb to pb
106  // for LHA strategy +-4. Also cumulative sum.
107  double weight() const;
108  double weightSum() const;
109  double lhaStrategy() const {return lhaStrategySave;}
110 
111  // Number of times other steps have been carried out.
112  int nISR() const {return nISRSave;}
113  int nFSRinProc() const {return nFSRinProcSave;}
114  int nFSRinRes() const {return nFSRinResSave;}
115 
116  // Maximum pT scales for MPI, ISR and FSR (in hard process).
117  double pTmaxMPI() const {return pTmaxMPISave;}
118  double pTmaxISR() const {return pTmaxISRSave;}
119  double pTmaxFSR() const {return pTmaxFSRSave;}
120 
121  // Current evolution scale (for UserHooks).
122  double pTnow() const {return pTnowSave;}
123 
124  // Impact parameter picture, global information
125  double a0MPI() const {return a0MPISave;}
126 
127  // Impact parameter picture, as set by hardest interaction.
128  double bMPI() const {return (bIsSet) ? bMPISave : 1.;}
129  double enhanceMPI() const {return (bIsSet) ? enhanceMPISave : 1.;}
130  double eMPI(int i) const {return (bIsSet) ? eMPISave[i] : 1.;}
131 
132  // Number of multiparton interactions, with code and pT for them.
133  int nMPI() const {return nMPISave;}
134  int codeMPI(int i) const {return codeMPISave[i];}
135  double pTMPI(int i) const {return pTMPISave[i];}
136  int iAMPI(int i) const {return iAMPISave[i];}
137  int iBMPI(int i) const {return iBMPISave[i];}
138 
139  // Cross section estimate.
140  long nTried() const {return nTry;}
141  long nSelected() const {return nSel;}
142  long nAccepted() const {return nAcc;}
143  double sigmaGen() const {return sigGen;}
144  double sigmaErr() const {return sigErr;}
145 
146  // Counters for number of loops in various places.
147  int getCounter( int i) const {return counters[i];}
148 
149  // Set or increase the value stored in a counter.
150  void setCounter( int i, int value = 0) {counters[i] = value;}
151  void addCounter( int i, int value = 1) {counters[i] += value;}
152 
153  // Reset to empty map of error messages.
154  void errorReset() {messages.clear();}
155 
156  // Print a message the first few times. Insert in database.
157  void errorMsg(string messageIn, string extraIn = " ",
158  bool showAlways = false, ostream& os = cout);
159 
160  // Provide total number of errors/aborts/warnings experienced to date.
161  int errorTotalNumber();
162 
163  // Print statistics on errors/aborts/warnings.
164  void errorStatistics(ostream& os = cout);
165 
166  // Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
167  void setTooLowPTmin(bool lowPTminIn) {lowPTmin = lowPTminIn;}
168 
169  // Set info on valence character of hard collision partons.
170  void setValence( bool isVal1In, bool isVal2In) {isVal1 = isVal1In;
171  isVal2 = isVal2In;}
172 
173  // Set and get some MPI/ISR/FSR properties needed for matching,
174  // i.e. mainly of internal relevance.
175  void hasHistory( bool hasHistoryIn) {hasHistorySave = hasHistoryIn;}
176  bool hasHistory() {return hasHistorySave;}
177  void zNowISR( double zNowIn) {zNowISRSave = zNowIn;}
178  double zNowISR() {return zNowISRSave;}
179  void pT2NowISR( double pT2NowIn) {pT2NowISRSave = pT2NowIn;}
180  double pT2NowISR() {return pT2NowISRSave;}
181 
182  // Return merging weight
183  double mergingWeight() const { return mergingWeightSave;}
184 
185 private:
186 
187  // Number of times the same error message is repeated, unless overridden.
188  static const int TIMESTOPRINT;
189 
190  // Allow conversion from mb to pb.
191  static const double CONVERTMB2PB;
192 
193  // Store common beam quantities.
194  int idASave, idBSave;
195  double pzASave, eASave,mASave, pzBSave, eBSave, mBSave, eCMSave, sSave;
196 
197  // Store initialization information.
198  bool lowPTmin;
199 
200  // Store common integrated cross section quantities.
201  long nTry, nSel, nAcc;
202  double sigGen, sigErr, wtAccSum;
203  int lhaStrategySave;
204 
205  // Store common MPI information
206  double a0MPISave;
207 
208  // Store current-event quantities.
209  bool isRes, isDiffA, isDiffB, isMB, isLH, hasSubSave, bIsSet, evolIsSet,
210  atEOF, isVal1, isVal2, hasHistorySave;
211  int codeSave, codeSubSave, nFinalSave, nFinalSubSave, nTotal,
212  id1Save, id2Save, nMPISave, nISRSave, nFSRinProcSave, nFSRinResSave;
213  double x1Save, x2Save, pdf1Save, pdf2Save, Q2FacSave, alphaEMSave,
214  alphaSSave, Q2RenSave, sH, tH, uH, pTH, m3H, m4H, thetaH, phiH,
215  weightSave, bMPISave, enhanceMPISave, pTmaxMPISave, pTmaxISRSave,
216  pTmaxFSRSave, pTnowSave, zNowISRSave, pT2NowISRSave;
217  string nameSave, nameSubSave;
218  vector<int> codeMPISave, iAMPISave, iBMPISave;
219  vector<double> pTMPISave, eMPISave;
220 
221  // Vector of various loop counters.
222  int counters[50];
223 
224  // Map for all error messages.
225  map<string, int> messages;
226 
227  // Friend classes allowed to set info.
228  friend class Pythia;
229  friend class ProcessLevel;
230  friend class ProcessContainer;
231  friend class PartonLevel;
232  friend class MultipartonInteractions;
233 
234  // Set info on the two incoming beams: only from Pythia class.
235  void setBeamA( int idAin, double pzAin, double eAin, double mAin) {
236  idASave = idAin; pzASave = pzAin; eASave = eAin; mASave = mAin;}
237  void setBeamB( int idBin, double pzBin, double eBin, double mBin) {
238  idBSave = idBin; pzBSave = pzBin; eBSave = eBin; mBSave = mBin;}
239  void setECM( double eCMin) {eCMSave = eCMin; sSave = eCMSave * eCMSave;}
240 
241  // Reset info for current event: only from Pythia class.
242  void clear() { isRes = isDiffA = isDiffB = isMB = isLH = hasSubSave
243  = bIsSet = evolIsSet = atEOF = isVal1 =isVal2 = hasHistorySave = false;
244  codeSave = codeSubSave = nFinalSave = nFinalSubSave = nTotal = id1Save
245  = id2Save = nMPISave = nISRSave = nFSRinProcSave = nFSRinResSave = 0;
246  x1Save = x2Save = pdf1Save = pdf2Save = Q2FacSave = alphaEMSave
247  = alphaSSave = Q2RenSave = sH = tH = uH = pTH = m3H = m4H = thetaH
248  = phiH = 0.; weightSave = bMPISave = enhanceMPISave = mergingWeightSave
249  = 1.; pTmaxMPISave = pTmaxISRSave = pTmaxFSRSave = pTnowSave
250  = zNowISRSave = pT2NowISRSave = 0.; nameSave = nameSubSave = " ";
251  codeMPISave.resize(0); iAMPISave.resize(0); iBMPISave.resize(0);
252  pTMPISave.resize(0); eMPISave.resize(0); }
253 
254  // Set info on the (sub)process: from ProcessLevel, ProcessContainer or
255  // MultipartonInteractions classes.
256  void setType( string nameIn, int codeIn, int nFinalIn,
257  bool isMinBiasIn = false, bool isResolvedIn = true,
258  bool isDiffractiveAin = false, bool isDiffractiveBin = false,
259  bool isLHAin = false) {nameSave = nameIn; codeSave = codeIn;
260  nFinalSave = nFinalIn; isMB = isMinBiasIn; isRes = isResolvedIn;
261  isDiffA = isDiffractiveAin; isDiffB = isDiffractiveBin; isLH = isLHAin;
262  nTotal = 2 + nFinalSave; bIsSet = false; hasSubSave = false;
263  nameSubSave = " "; codeSubSave = 0; nFinalSubSave = 0; evolIsSet = false;}
264  void setSubType( string nameSubIn, int codeSubIn, int nFinalSubIn) {
265  hasSubSave = true; nameSubSave = nameSubIn; codeSubSave = codeSubIn;
266  nFinalSubSave = nFinalSubIn;}
267  void setPDFalpha( int id1In, int id2In, double pdf1In, double pdf2In,
268  double Q2FacIn, double alphaEMIn, double alphaSIn, double Q2RenIn)
269  {id1Save = id1In; id2Save = id2In; pdf1Save = pdf1In; pdf2Save = pdf2In;
270  Q2FacSave = Q2FacIn; alphaEMSave = alphaEMIn; alphaSSave = alphaSIn;
271  Q2RenSave = Q2RenIn;}
272  void setKin( double x1In, double x2In, double sHatIn, double tHatIn,
273  double uHatIn, double pTHatIn, double m3HatIn, double m4HatIn,
274  double thetaHatIn, double phiHatIn) {x1Save = x1In; x2Save = x2In;
275  sH = sHatIn; tH = tHatIn; uH = uHatIn; pTH = pTHatIn; m3H = m3HatIn;
276  m4H = m4HatIn; thetaH = thetaHatIn; phiH = phiHatIn;}
277  void setTypeMPI( int codeMPIIn, double pTMPIIn, int iAMPIIn = 0,
278  int iBMPIIn = 0, double eMPIIn = 1.) {codeMPISave.push_back(codeMPIIn);
279  pTMPISave.push_back(pTMPIIn); iAMPISave.push_back(iAMPIIn);
280  iBMPISave.push_back(iBMPIIn); eMPISave.push_back(eMPIIn); }
281 
282  // Set info on cross section: from ProcessLevel.
283  void setSigma( long nTryIn, long nSelIn, long nAccIn, double sigGenIn,
284  double sigErrIn, double wtAccSumIn) { nTry = nTryIn; nSel = nSelIn;
285  nAcc = nAccIn; sigGen = sigGenIn; sigErr = sigErrIn;
286  wtAccSum = wtAccSumIn; }
287 
288  // Set info on impact parameter: from PartonLevel.
289  void setImpact( double bMPIIn, double enhanceMPIIn) {bMPISave = bMPIIn;
290  enhanceMPISave = eMPISave[0] = enhanceMPIIn, bIsSet = true;}
291 
292  // Set info on pTmax scales and number of evolution steps: from PartonLevel.
293  void setPartEvolved( int nMPIIn, int nISRIn) {
294  nMPISave = nMPIIn; nISRSave = nISRIn;}
295  void setEvolution( double pTmaxMPIIn, double pTmaxISRIn, double pTmaxFSRIn,
296  int nMPIIn, int nISRIn, int nFSRinProcIn, int nFSRinResIn) {
297  pTmaxMPISave = pTmaxMPIIn; pTmaxISRSave = pTmaxISRIn;
298  pTmaxFSRSave = pTmaxFSRIn; nMPISave = nMPIIn; nISRSave = nISRIn;
299  nFSRinProcSave = nFSRinProcIn; nFSRinResSave = nFSRinResIn;
300  evolIsSet = true;}
301 
302  // Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
303  void setPTnow( double pTnowIn) {pTnowSave = pTnowIn;}
304 
305  // Set a0 from MultipartonInteractions.
306  void a0MPI(double a0MPIin) {a0MPISave = a0MPIin;}
307 
308  // Set info whether reading of Les Houches Accord file at end.
309  void setEndOfFile( bool atEOFin) {atEOF = atEOFin;}
310 
311  // Set event weight; currently only for Les Houches description.
312  void setWeight( double weightIn, int lhaStrategyIn) {
313  weightSave = weightIn; lhaStrategySave = lhaStrategyIn; }
314 
315  // Set merging weight for event
316  double mergingWeightSave;
317  void setMergingWeight( double weightIn) { mergingWeightSave = weightIn;}
318 
319 };
320 
321 //==========================================================================
322 
323 } // end namespace Pythia8
324 
325 #endif // Pythia8_Info_H