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) 2018 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 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 "Pythia8/PythiaStdlib.h"
13 #include "Pythia8/LHEF3.h"
14 #include "Pythia8/Basics.h"
15 
16 namespace Pythia8 {
17 
18 // Forard declaration of HIInfo class.
19 class HIInfo;
20 
21 //==========================================================================
22 
23 // The Info class contains a mixed bag of information on the event
24 // generation activity, especially on the current subprocess properties,
25 // and on the number of errors encountered. This is used by the
26 // generation machinery, but can also be read by the user.
27 // Note: some methods that maybe should not be accessible to the user
28 // are still public, to work also for user-written FSR/ISR classes.
29 
30 class Info {
31 
32 public:
33 
34  // Constructor.
35  Info() : LHEFversionSave(0), initrwgt(NULL), generators(NULL),
36  weightgroups(NULL), init_weights(NULL), eventAttributes(NULL),
37  weights_detailed(NULL), weights_compressed(NULL), scales(NULL),
38  weights(NULL), rwgt(NULL), hiinfo(0), eCMSave(0.),
39  lowPTmin(false), a0MPISave(0.),
40  abortPartonLevel(false), weightCKKWLSave(1.), weightFIRSTSave(0.) {
41  for (int i = 0; i < 40; ++i) counters[i] = 0;
42  setNWeights(1);}
43 
44  // Listing of most available information on current event.
45  void list() const;
46 
47  // Beam particles (in rest frame). CM energy of event.
48  int idA() const {return idASave;}
49  int idB() const {return idBSave;}
50  double pzA() const {return pzASave;}
51  double pzB() const {return pzBSave;}
52  double eA() const {return eASave;}
53  double eB() const {return eBSave;}
54  double mA() const {return mASave;}
55  double mB() const {return mBSave;}
56  double eCM() const {return eCMSave;}
57  double s() const {return sSave;}
58 
59  // Warnings from initialization.
60  bool tooLowPTmin() const {return lowPTmin;}
61 
62  // Process name and code, and the number of final-state particles.
63  string name() const {return nameSave;}
64  int code() const {return codeSave;}
65  int nFinal() const {return nFinalSave;}
66 
67  // Are beam particles resolved, with pdf's? Are they diffractive?
68  bool isResolved() const {return isRes;}
69  bool isDiffractiveA() const {return isDiffA;}
70  bool isDiffractiveB() const {return isDiffB;}
71  bool isDiffractiveC() const {return isDiffC;}
72  bool isNonDiffractive() const {return isND;}
73  // Retained for backwards compatibility.
74  bool isMinBias() const {return isND;}
75 
76  // Information for Les Houches Accord and reading files.
77  bool isLHA() const {return isLH;}
78  bool atEndOfFile() const {return atEOF;}
79 
80  // For nondiffractive and Les Houches Accord identify hardest subprocess.
81  bool hasSub(int i = 0) const {return hasSubSave[i];}
82  string nameSub(int i = 0) const {return nameSubSave[i];}
83  int codeSub(int i = 0) const {return codeSubSave[i];}
84  int nFinalSub(int i = 0) const {return nFinalSubSave[i];}
85 
86  // Incoming parton flavours and x values.
87  int id1(int i = 0) const {return id1Save[i];}
88  int id2(int i = 0) const {return id2Save[i];}
89  double x1(int i = 0) const {return x1Save[i];}
90  double x2(int i = 0) const {return x2Save[i];}
91  double y(int i = 0) const {return 0.5 * log( x1Save[i] / x2Save[i]);}
92  double tau(int i = 0) const {return x1Save[i] * x2Save[i];}
93 
94  // Hard process flavours, x values, parton densities, couplings, Q2 scales.
95  int id1pdf(int i = 0) const {return id1pdfSave[i];}
96  int id2pdf(int i = 0) const {return id2pdfSave[i];}
97  double x1pdf(int i = 0) const {return x1pdfSave[i];}
98  double x2pdf(int i = 0) const {return x2pdfSave[i];}
99  double pdf1(int i = 0) const {return pdf1Save[i];}
100  double pdf2(int i = 0) const {return pdf2Save[i];}
101  double QFac(int i = 0) const {return sqrtpos(Q2FacSave[i]);}
102  double Q2Fac(int i = 0) const {return Q2FacSave[i];}
103  bool isValence1() const {return isVal1;}
104  bool isValence2() const {return isVal2;}
105  double alphaS(int i = 0) const {return alphaSSave[i];}
106  double alphaEM(int i = 0) const {return alphaEMSave[i];}
107  double QRen(int i = 0) const {return sqrtpos(Q2RenSave[i]);}
108  double Q2Ren(int i = 0) const {return Q2RenSave[i];}
109  double scalup(int i = 0) const {return scalupSave[i];}
110 
111  // Kinematics of photons from lepton beams.
112  double xGammaA() const {return x1GammaSave;}
113  double xGammaB() const {return x2GammaSave;}
114  double Q2GammaA() const {return Q2Gamma1Save;}
115  double Q2GammaB() const {return Q2Gamma2Save;}
116  double eCMsub() const {return eCMsubSave;}
117  double thetaScatLepA() const {return thetaLepton1;}
118  double thetaScatLepB() const {return thetaLepton2;}
119  double sHatNew() const {return sHatNewSave;}
120  int photonMode() const {return gammaModeEvent;}
121 
122  // Information on VMD state inside a photon.
123  bool isVMDstateA() const {return isVMDstateAEvent;}
124  bool isVMDstateB() const {return isVMDstateBEvent;}
125  int idVMDA() const {return idVMDASave;}
126  int idVMDB() const {return idVMDBSave;}
127  double mVMDA() const {return mVMDASave;}
128  double mVMDB() const {return mVMDBSave;}
129  double scaleVMDA() const {return scaleVMDASave;}
130  double scaleVMDB() const {return scaleVMDBSave;}
131 
132  // Mandelstam variables (notation as if subcollision).
133  double mHat(int i = 0) const {return sqrt(sH[i]);}
134  double sHat(int i = 0) const {return sH[i];}
135  double tHat(int i = 0) const {return tH[i];}
136  double uHat(int i = 0) const {return uH[i];}
137  double pTHat(int i = 0) const {return pTH[i];}
138  double pT2Hat(int i = 0) const {return pTH[i] * pTH[i];}
139  double m3Hat(int i = 0) const {return m3H[i];}
140  double m4Hat(int i = 0) const {return m4H[i];}
141  double thetaHat(int i = 0) const {return thetaH[i];}
142  double phiHat(int i = 0) const {return phiH[i];}
143 
144  // Weight of current event; normally 1, but used for Les Houches events
145  // or when reweighting phase space selection. Conversion from mb to pb
146  // for LHA strategy +-4. Uncertainty variations can be accessed by
147  // providing an index >= 1 (0 = no variations). Also cumulative sum.
148  double weight(int i=0) const;
149  double weightSum() const;
150  double lhaStrategy() const {return lhaStrategySave;}
151 
152  // Further access to uncertainty weights: number and labels
153  int nWeights() { return weightSave.size(); }
154  string weightLabel(int iWeight) {
155  if (iWeight >= 0 && iWeight < (int)weightLabelSave.size())
156  return weightLabelSave[iWeight];
157  else return "";}
158 
159  void initUncertainties(vector<string>*, bool = false);
160  int nVariationGroups() { return externalVariationsSize; }
161  string getGroupName(int iGN) {
162  string tmpString("Null");
163  if( iGN < 0 || iGN >= externalVariationsSize ) return tmpString;
164  return externalGroupNames[iGN];
165  }
166  double getGroupWeight(int iGW) {
167  double tempWeight(1.0);
168  if( iGW < 0 || iGW >= externalVariationsSize ) return tempWeight;
169  for( vector<int>::const_iterator cit = externalMap[iGW].begin();
170  cit < externalMap[iGW].end(); ++cit ) tempWeight *= weightSave[*cit];
171  return tempWeight;
172  }
173  string getInitialName(int iG) { return initialNameSave[iG]; }
174  // Variations that must be known by TimeShower and Spaceshower
175  map<int,double> varPDFplus, varPDFminus, varPDFmember;
176 
177  // Number of times other steps have been carried out.
178  int nISR() const {return nISRSave;}
179  int nFSRinProc() const {return nFSRinProcSave;}
180  int nFSRinRes() const {return nFSRinResSave;}
181 
182  // Maximum pT scales for MPI, ISR and FSR (in hard process).
183  double pTmaxMPI() const {return pTmaxMPISave;}
184  double pTmaxISR() const {return pTmaxISRSave;}
185  double pTmaxFSR() const {return pTmaxFSRSave;}
186 
187  // Current evolution scale (for UserHooks).
188  double pTnow() const {return pTnowSave;}
189 
190  // Impact parameter picture, global information
191  double a0MPI() const {return a0MPISave;}
192 
193  // Impact parameter picture, as set by hardest interaction.
194  double bMPI() const {return (bIsSet) ? bMPISave : 1.;}
195  double enhanceMPI() const {return (bIsSet) ? enhanceMPISave : 1.;}
196  double enhanceMPIavg() const {return (bIsSet) ? enhanceMPIavgSave : 1.;}
197  double eMPI(int i) const {return (bIsSet) ? eMPISave[i] : 1.;}
198  double bMPIold() const {return (bIsSet) ? bMPIoldSave : 1.;}
199  double enhanceMPIold() const {return (bIsSet) ? enhanceMPIoldSave : 1.;}
200  double enhanceMPIoldavg() const {return (bIsSet)
201  ? enhanceMPIoldavgSave : 1.;}
202 
203  // Number of multiparton interactions, with code and pT for them.
204  int nMPI() const {return nMPISave;}
205  int codeMPI(int i) const {return codeMPISave[i];}
206  double pTMPI(int i) const {return pTMPISave[i];}
207  int iAMPI(int i) const {return iAMPISave[i];}
208  int iBMPI(int i) const {return iBMPISave[i];}
209 
210  // Cross section estimate, optionally process by process.
211  vector<int> codesHard();
212  string nameProc(int i = 0) {return (i == 0) ? "sum"
213  : ( (procNameM[i] == "") ? "unknown process" : procNameM[i] );}
214  long nTried(int i = 0) {return (i == 0) ? nTry : nTryM[i];}
215  long nSelected(int i = 0) {return (i == 0) ? nSel : nSelM[i];}
216  long nAccepted(int i = 0) {return (i == 0) ? nAcc : nAccM[i];}
217  double sigmaGen(int i = 0) {return (i == 0) ? sigGen : sigGenM[i];}
218  double sigmaErr(int i = 0) {return (i == 0) ? sigErr : sigErrM[i];}
219 
220  // Counters for number of loops in various places.
221  int getCounter( int i) const {return counters[i];}
222 
223  // Set or increase the value stored in a counter.
224  void setCounter( int i, int value = 0) {counters[i] = value;}
225  void addCounter( int i, int value = 1) {counters[i] += value;}
226 
227  // Reset to empty map of error messages.
228  void errorReset() {messages.clear();}
229 
230  // Print a message the first few times. Insert in database.
231  void errorMsg(string messageIn, string extraIn = " ",
232  bool showAlways = false);
233 
234  // Provide total number of errors/aborts/warnings experienced to date.
235  int errorTotalNumber();
236 
237  // Print statistics on errors/aborts/warnings.
238  void errorStatistics();
239 
240  // Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
241  void setTooLowPTmin(bool lowPTminIn) {lowPTmin = lowPTminIn;}
242 
243  // Set info on valence character of hard collision partons.
244  void setValence( bool isVal1In, bool isVal2In) {isVal1 = isVal1In;
245  isVal2 = isVal2In;}
246 
247  // Set and get some MPI/ISR/FSR properties needed for matching,
248  // i.e. mainly of internal relevance.
249  void hasHistory( bool hasHistoryIn) {hasHistorySave = hasHistoryIn;}
250  bool hasHistory() {return hasHistorySave;}
251  void zNowISR( double zNowIn) {zNowISRSave = zNowIn;}
252  double zNowISR() {return zNowISRSave;}
253  void pT2NowISR( double pT2NowIn) {pT2NowISRSave = pT2NowIn;}
254  double pT2NowISR() {return pT2NowISRSave;}
255 
256  // Update a particular event weight, first entry by default.
257  void updateWeight( double weightIn, int i=0) { weightSave[i] = weightIn;}
258 
259  // Return CKKW-L weight.
260  double getWeightCKKWL() const { return weightCKKWLSave;}
261  // Set CKKW-L weight.
262  void setWeightCKKWL( double weightIn) { weightCKKWLSave = weightIn;}
263  // Return merging weight.
264  double mergingWeight() const { return weightCKKWLSave;}
265 
266  // Return the complete NLO weight.
267  double mergingWeightNLO() const {
268  return (weightCKKWLSave - weightFIRSTSave); }
269  // Return the O(\alpha_s)-term of the CKKW-L weight.
270  double getWeightFIRST() const { return weightFIRSTSave;}
271  // Set the O(\alpha_s)-term of the CKKW-L weight.
272  void setWeightFIRST( double weightIn) { weightFIRSTSave = weightIn;}
273 
274  // Return an LHEF header
275  string header(const string &key) {
276  if (headers.find(key) == headers.end()) return "";
277  else return headers[key];
278  }
279 
280  // Return a list of all header key names
281  vector < string > headerKeys();
282 
283  // Return the number of processes in the LHEF.
284  int nProcessesLHEF() { return int(sigmaLHEFSave.size());}
285  // Return the cross section information read from LHEF.
286  double sigmaLHEF(int iProcess) { return sigmaLHEFSave[iProcess];}
287 
288  // LHEF3 information: Public for easy access.
289  int LHEFversionSave;
290 
291  // Save process information as read from init block of LHEF.
292  vector<double> sigmaLHEFSave;
293 
294  // Contents of the LHAinitrwgt tag
295  LHAinitrwgt *initrwgt;
296 
297  // Contents of the LHAgenerator tags.
298  vector<LHAgenerator > *generators;
299 
300  // A map of the LHAweightgroup tags, indexed by name.
301  map<string,LHAweightgroup > *weightgroups;
302 
303  // A map of the LHAweight tags, indexed by name.
304  map<string,LHAweight > *init_weights;
305 
306  // Store current-event Les Houches event tags.
307  map<string, string > *eventAttributes;
308 
309  // The weights associated with this event, as given by the LHAwgt tags
310  map<string,double > *weights_detailed;
311 
312  // The weights associated with this event, as given by the LHAweights tags
313  vector<double > *weights_compressed;
314 
315  // Contents of the LHAscales tag.
316  LHAscales *scales;
317 
318  // Contents of the LHAweights tag (compressed format)
319  LHAweights *weights;
320 
321  // Contents of the LHArwgt tag (detailed format)
322  LHArwgt *rwgt;
323 
324  // Vectorized version of LHArwgt tag, for easy and ordered access.
325  vector<double> weights_detailed_vector;
326 
327  // Value of the unit event weight read from a Les Houches event, necessary
328  // if additional weights in an otherwise unweighted input file are in
329  // relation to this number.
330  double eventWeightLHEF;
331 
332  // Set the LHEF3 objects read from the init and header blocks.
333  void setLHEF3InitInfo();
334  void setLHEF3InitInfo( int LHEFversionIn, LHAinitrwgt *initrwgtIn,
335  vector<LHAgenerator> *generatorsIn,
336  map<string,LHAweightgroup> *weightgroupsIn,
337  map<string,LHAweight> *init_weightsIn, string headerBlockIn );
338  // Set the LHEF3 objects read from the event block.
339  void setLHEF3EventInfo();
340  void setLHEF3EventInfo( map<string, string> *eventAttributesIn,
341  map<string, double > *weights_detailedIn,
342  vector<double > *weights_compressedIn,
343  LHAscales *scalesIn, LHAweights *weightsIn,
344  LHArwgt *rwgtIn, vector<double> weights_detailed_vecIn,
345  string eventCommentsIn, double eventWeightLHEFIn );
346 
347  // Retrieve events tag information.
348  string getEventAttribute(string key, bool doRemoveWhitespace = false);
349 
350  // Retrieve LHEF version
351  int LHEFversion();
352 
353  // Retrieve initrwgt tag information.
354  unsigned int getInitrwgtSize();
355 
356  // Retrieve generator tag information.
357  unsigned int getGeneratorSize();
358  string getGeneratorValue(unsigned int n = 0);
359  string getGeneratorAttribute( unsigned int n, string key,
360  bool doRemoveWhitespace = false);
361 
362  // Retrieve rwgt tag information.
363  unsigned int getWeightsDetailedSize();
364  double getWeightsDetailedValue(string n);
365  string getWeightsDetailedAttribute(string n, string key,
366  bool doRemoveWhitespace = false);
367 
368  // Retrieve weights tag information.
369  unsigned int getWeightsCompressedSize();
370  double getWeightsCompressedValue(unsigned int n);
371  string getWeightsCompressedAttribute(string key,
372  bool doRemoveWhitespace = false);
373 
374  // Retrieve scales tag information.
375  string getScalesValue(bool doRemoveWhitespace = false);
376  double getScalesAttribute(string key);
377 
378  // Retrieve complete header block and event comments
379  // Retrieve scales tag information.
380  string getHeaderBlock() { return headerBlock;}
381  string getEventComments() { return eventComments;}
382 
383  // Set LHEF headers
384  void setHeader(const string &key, const string &val)
385  { headers[key] = val; }
386 
387  // Set abort in parton level.
388  void setAbortPartonLevel(bool abortIn) {abortPartonLevel = abortIn;}
389  bool getAbortPartonLevel() {return abortPartonLevel;}
390 
391  // Get information on hard diffractive events.
392  bool hasUnresolvedBeams() const {return hasUnresBeams;}
393  bool hasPomPsystem() const {return hasPomPsys;}
394  bool isHardDiffractive() const {return isHardDiffA || isHardDiffB;}
395  bool isHardDiffractiveA() const {return isHardDiffA;}
396  bool isHardDiffractiveB() const {return isHardDiffB;}
397  double xPomeronA() const {return xPomA;}
398  double xPomeronB() const {return xPomB;}
399  double tPomeronA() const {return tPomA;}
400  double tPomeronB() const {return tPomB;}
401 
402  // History information needed to setup the weak shower for 2 -> n.
403  vector<int> getWeakModes() {return weakModes;}
404  vector<pair<int,int> > getWeakDipoles() {return weakDipoles;}
405  vector<Vec4> getWeakMomenta() {return weakMomenta;}
406  vector<int> getWeak2to2lines() {return weak2to2lines;}
407  void setWeakModes(vector<int> weakModesIn) {weakModes = weakModesIn;}
408  void setWeakDipoles(vector<pair<int,int> > weakDipolesIn)
409  {weakDipoles = weakDipolesIn;}
410  void setWeakMomenta(vector<Vec4> weakMomentaIn)
411  {weakMomenta = weakMomentaIn;}
412  void setWeak2to2lines(vector<int> weak2to2linesIn)
413  {weak2to2lines = weak2to2linesIn;}
414 
415  // Access to a HIInfo object containing information about a
416  // HeavyIons run and the current Event. (Is NULL if HeavyIons object
417  // is inactive.
418  HIInfo * hiinfo;
419 
420  private:
421 
422  // Number of times the same error message is repeated, unless overridden.
423  static const int TIMESTOPRINT;
424 
425  // Allow conversion from mb to pb.
426  static const double CONVERTMB2PB;
427 
428  // Store common beam quantities.
429  int idASave, idBSave;
430  double pzASave, eASave,mASave, pzBSave, eBSave, mBSave, eCMSave, sSave;
431 
432  // Store initialization information.
433  bool lowPTmin;
434 
435  // Store common integrated cross section quantities.
436  long nTry, nSel, nAcc;
437  double sigGen, sigErr, wtAccSum;
438  map<int, string> procNameM;
439  map<int, long> nTryM, nSelM, nAccM;
440  map<int, double> sigGenM, sigErrM;
441  int lhaStrategySave;
442 
443  // Store common MPI information.
444  double a0MPISave;
445 
446  // Store current-event quantities.
447  bool isRes, isDiffA, isDiffB, isDiffC, isND, isLH, hasSubSave[4],
448  bIsSet, evolIsSet, atEOF, isVal1, isVal2, hasHistorySave,
449  abortPartonLevel, isHardDiffA, isHardDiffB, hasUnresBeams,
450  hasPomPsys;
451  int codeSave, codeSubSave[4], nFinalSave, nFinalSubSave[4], nTotal,
452  id1Save[4], id2Save[4], id1pdfSave[4], id2pdfSave[4], nMPISave,
453  nISRSave, nFSRinProcSave, nFSRinResSave;
454  double x1Save[4], x2Save[4], x1pdfSave[4], x2pdfSave[4], pdf1Save[4],
455  pdf2Save[4], Q2FacSave[4], alphaEMSave[4], alphaSSave[4],
456  Q2RenSave[4], scalupSave[4], sH[4], tH[4], uH[4], pTH[4], m3H[4],
457  m4H[4], thetaH[4], phiH[4], bMPISave, enhanceMPISave,
458  enhanceMPIavgSave, bMPIoldSave, enhanceMPIoldSave,
459  enhanceMPIoldavgSave, pTmaxMPISave, pTmaxISRSave, pTmaxFSRSave,
460  pTnowSave, zNowISRSave, pT2NowISRSave, xPomA, xPomB, tPomA, tPomB;
461  string nameSave, nameSubSave[4];
462  vector<int> codeMPISave, iAMPISave, iBMPISave;
463  vector<double> pTMPISave, eMPISave, weightSave;
464  vector<string> weightLabelSave;
465  vector<string> externalVariations;
466  vector<vector<string> > externalVarNames;
467  vector<string> externalGroupNames;
468  vector<string> initialNameSave;
469  vector<vector<int> > externalMap;
470  int externalVariationsSize;
471 
472  // Variables related to photon kinematics.
473  bool isVMDstateAEvent, isVMDstateBEvent;
474  int gammaModeEvent, idVMDASave, idVMDBSave;
475  double x1GammaSave, x2GammaSave, Q2Gamma1Save, Q2Gamma2Save, eCMsubSave,
476  thetaLepton1, thetaLepton2, sHatNewSave, mVMDASave, mVMDBSave,
477  scaleVMDASave, scaleVMDBSave;
478 
479  // Vector of various loop counters.
480  int counters[50];
481 
482  // Map for all error messages.
483  map<string, int> messages;
484 
485  // Map for LHEF headers.
486  map<string, string> headers;
487 
488  // Strings for complete header block and event comments.
489  string headerBlock, eventComments;
490 
491  // Map for plugin libraries.
492  map<string, pair<void*, int> > plugins;
493 
494  // Friend classes allowed to set info.
495  friend class Pythia;
496  friend class ProcessLevel;
497  friend class ProcessContainer;
498  friend class PartonLevel;
499  friend class MultipartonInteractions;
500  friend class LHAup;
501  friend class LHAPDF;
502  friend class Diffraction;
503  friend class Settings;
504  friend class TimeShower;
505  friend class SpaceShower;
506  friend class GammaKinematics;
507  friend class HeavyIons;
508 
509  // Set info on the two incoming beams: only from Pythia class.
510  void setBeamA( int idAin, double pzAin, double eAin, double mAin) {
511  idASave = idAin; pzASave = pzAin; eASave = eAin; mASave = mAin;}
512  void setBeamB( int idBin, double pzBin, double eBin, double mBin) {
513  idBSave = idBin; pzBSave = pzBin; eBSave = eBin; mBSave = mBin;}
514  void setECM( double eCMin) {eCMSave = eCMin; sSave = eCMSave * eCMSave;}
515 
516  // Set info related to gamma+gamma subcollision.
517  void setX1Gamma( double x1GammaIn) { x1GammaSave = x1GammaIn; }
518  void setX2Gamma( double x2GammaIn) { x2GammaSave = x2GammaIn; }
519  void setQ2Gamma1( double Q2gammaIn) { Q2Gamma1Save = Q2gammaIn; }
520  void setQ2Gamma2( double Q2gammaIn) { Q2Gamma2Save = Q2gammaIn; }
521  void setTheta1( double theta1In) { thetaLepton1 = theta1In; }
522  void setTheta2( double theta2In) { thetaLepton2 = theta2In; }
523  void setECMsub( double eCMsubIn) { eCMsubSave = eCMsubIn; }
524  void setsHatNew( double sHatNewIn) { sHatNewSave = sHatNewIn; }
525  void setGammaMode( double gammaModeIn) { gammaModeEvent = gammaModeIn; }
526  // Set info on VMD state.
527  void setVMDstateA(bool isVMDAIn, int idAIn, double mAIn, double scaleAIn)
528  {isVMDstateAEvent = isVMDAIn; idVMDASave = idAIn; mVMDASave = mAIn;
529  scaleVMDASave = scaleAIn;}
530  void setVMDstateB(bool isVMDBIn, int idBIn, double mBIn, double scaleBIn)
531  {isVMDstateBEvent = isVMDBIn; idVMDBSave = idBIn; mVMDBSave = mBIn;
532  scaleVMDBSave = scaleBIn;}
533 
534  // Reset info for current event: only from Pythia class.
535  void clear() {
536  isRes = isDiffA = isDiffB = isDiffC = isND = isLH = bIsSet
537  = evolIsSet = atEOF = isVal1 = isVal2 = hasHistorySave
538  = isHardDiffA = isHardDiffB = hasUnresBeams = hasPomPsys = false;
539  codeSave = nFinalSave = nTotal = nMPISave = nISRSave = nFSRinProcSave
540  = nFSRinResSave = 0;
541  bMPISave = enhanceMPISave = enhanceMPIavgSave = bMPIoldSave
542  = enhanceMPIoldSave = enhanceMPIoldavgSave = weightCKKWLSave = 1.;
543  pTmaxMPISave = pTmaxISRSave = pTmaxFSRSave = pTnowSave = zNowISRSave
544  = pT2NowISRSave = weightFIRSTSave = 0.;
545  nameSave = " ";
546  for (int i = 0; i < 4; ++i) {
547  hasSubSave[i] = false;
548  codeSubSave[i] = nFinalSubSave[i] = id1pdfSave[i] = id2pdfSave[i]
549  = id1Save[i] = id2Save[i] = 0;
550  x1pdfSave[i] = x2pdfSave[i] = pdf1Save[i] = pdf2Save[i]
551  = Q2FacSave[i] = alphaEMSave[i] = alphaSSave[i] = Q2RenSave[i]
552  = scalupSave[i] = x1Save[i] = x2Save[i] = sH[i] = tH[i] = uH[i]
553  = pTH[i] = m3H[i] = m4H[i] = thetaH[i] = phiH[i] = 0.;
554  nameSubSave[i] = " ";
555  }
556  codeMPISave.resize(0); iAMPISave.resize(0); iBMPISave.resize(0);
557  pTMPISave.resize(0); eMPISave.resize(0); setHardDiff();
558  for (int i = 0; i < nWeights(); ++i) weightSave[i]=1.;}
559 
560  // Reset info arrays only for parton and hadron level.
561  int sizeMPIarrays() const { return codeMPISave.size(); }
562  void resizeMPIarrays(int newSize) { codeMPISave.resize(newSize);
563  iAMPISave.resize(newSize); iBMPISave.resize(newSize);
564  pTMPISave.resize(newSize); eMPISave.resize(newSize); }
565 
566  // Set info on the (sub)process: from ProcessLevel, ProcessContainer or
567  // MultipartonInteractions classes.
568  void setType( string nameIn, int codeIn, int nFinalIn,
569  bool isNonDiffIn = false, bool isResolvedIn = true,
570  bool isDiffractiveAin = false, bool isDiffractiveBin = false,
571  bool isDiffractiveCin = false, bool isLHAin = false) {
572  nameSave = nameIn; codeSave = codeIn; nFinalSave = nFinalIn;
573  isND = isNonDiffIn; isRes = isResolvedIn; isDiffA = isDiffractiveAin;
574  isDiffB = isDiffractiveBin; isDiffC = isDiffractiveCin; isLH = isLHAin;
575  nTotal = 2 + nFinalSave; bIsSet = false; hasSubSave[0] = false;
576  nameSubSave[0] = " "; codeSubSave[0] = 0; nFinalSubSave[0] = 0;
577  evolIsSet = false;}
578  void setSubType( int iDS, string nameSubIn, int codeSubIn,
579  int nFinalSubIn) { hasSubSave[iDS] = true; nameSubSave[iDS] = nameSubIn;
580  codeSubSave[iDS] = codeSubIn; nFinalSubSave[iDS] = nFinalSubIn;}
581  void setPDFalpha( int iDS, int id1pdfIn, int id2pdfIn, double x1pdfIn,
582  double x2pdfIn, double pdf1In, double pdf2In, double Q2FacIn,
583  double alphaEMIn, double alphaSIn, double Q2RenIn, double scalupIn) {
584  id1pdfSave[iDS] = id1pdfIn; id2pdfSave[iDS] = id2pdfIn;
585  x1pdfSave[iDS] = x1pdfIn; x2pdfSave[iDS] = x2pdfIn;
586  pdf1Save[iDS] = pdf1In; pdf2Save[iDS] = pdf2In; Q2FacSave[iDS] = Q2FacIn;
587  alphaEMSave[iDS] = alphaEMIn; alphaSSave[iDS] = alphaSIn;
588  Q2RenSave[iDS] = Q2RenIn; scalupSave[iDS] = scalupIn;}
589  void setScalup( int iDS, double scalupIn) {scalupSave[iDS] = scalupIn;}
590  void setKin( int iDS, int id1In, int id2In, double x1In, double x2In,
591  double sHatIn, double tHatIn, double uHatIn, double pTHatIn,
592  double m3HatIn, double m4HatIn, double thetaHatIn, double phiHatIn) {
593  id1Save[iDS] = id1In; id2Save[iDS] = id2In; x1Save[iDS] = x1In;
594  x2Save[iDS] = x2In; sH[iDS] = sHatIn; tH[iDS] = tHatIn;
595  uH[iDS] = uHatIn; pTH[iDS] = pTHatIn; m3H[iDS] = m3HatIn;
596  m4H[iDS] = m4HatIn; thetaH[iDS] = thetaHatIn; phiH[iDS] = phiHatIn;}
597  void setTypeMPI( int codeMPIIn, double pTMPIIn, int iAMPIIn = 0,
598  int iBMPIIn = 0, double eMPIIn = 1.) {codeMPISave.push_back(codeMPIIn);
599  pTMPISave.push_back(pTMPIIn); iAMPISave.push_back(iAMPIIn);
600  iBMPISave.push_back(iBMPIIn); eMPISave.push_back(eMPIIn); }
601 
602  // Set info on cross section: from ProcessLevel (reset in Pythia).
603  void sigmaReset() { nTry = nSel = nAcc = 0; sigGen = sigErr = wtAccSum = 0.;
604  procNameM.clear(); nTryM.clear(); nSelM.clear(); nAccM.clear();
605  sigGenM.clear(); sigErrM.clear();}
606  void setSigma( int i, string procNameIn, long nTryIn, long nSelIn,
607  long nAccIn, double sigGenIn, double sigErrIn, double wtAccSumIn) {
608  if (i == 0) {nTry = nTryIn; nSel = nSelIn; nAcc = nAccIn;
609  sigGen = sigGenIn; sigErr = sigErrIn; wtAccSum = wtAccSumIn;}
610  else {procNameM[i] = procNameIn; nTryM[i] = nTryIn; nSelM[i] = nSelIn;
611  nAccM[i] = nAccIn; sigGenM[i] = sigGenIn; sigErrM[i] = sigErrIn;} }
612  void addSigma( int i, long nTryIn, long nSelIn, long nAccIn, double sigGenIn,
613  double sigErrIn) { nTryM[i] += nTryIn; nSelM[i] += nSelIn;
614  nAccM[i] += nAccIn; sigGenM[i] += sigGenIn;
615  sigErrM[i] = sqrtpos(sigErrM[i]*sigErrM[i] + sigErrIn*sigErrIn); }
616 
617  // Set info on impact parameter: from PartonLevel.
618  void setImpact( double bMPIIn, double enhanceMPIIn, double enhanceMPIavgIn,
619  bool bIsSetIn = true, bool pushBack = false) {
620  if (pushBack) {bMPIoldSave = bMPISave; enhanceMPIoldSave = enhanceMPISave;
621  enhanceMPIoldavgSave = enhanceMPIavgSave;}
622  bMPISave = bMPIIn; enhanceMPISave = eMPISave[0] = enhanceMPIIn,
623  enhanceMPIavgSave = enhanceMPIavgIn, bIsSet = bIsSetIn;}
624 
625  // Set info on pTmax scales and number of evolution steps: from PartonLevel.
626  void setPartEvolved( int nMPIIn, int nISRIn) {
627  nMPISave = nMPIIn; nISRSave = nISRIn;}
628  void setEvolution( double pTmaxMPIIn, double pTmaxISRIn, double pTmaxFSRIn,
629  int nMPIIn, int nISRIn, int nFSRinProcIn, int nFSRinResIn) {
630  pTmaxMPISave = pTmaxMPIIn; pTmaxISRSave = pTmaxISRIn;
631  pTmaxFSRSave = pTmaxFSRIn; nMPISave = nMPIIn; nISRSave = nISRIn;
632  nFSRinProcSave = nFSRinProcIn; nFSRinResSave = nFSRinResIn;
633  evolIsSet = true;}
634 
635  // Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
636  void setPTnow( double pTnowIn) {pTnowSave = pTnowIn;}
637 
638  // Set a0 from MultipartonInteractions.
639  void seta0MPI(double a0MPIin) {a0MPISave = a0MPIin;}
640 
641  // Set info whether reading of Les Houches Accord file at end.
642  void setEndOfFile( bool atEOFin) {atEOF = atEOFin;}
643 
644  // Set number and labels of weights (for uncertainty evaluations).
645  void setNWeights(int mWeights) {
646  mWeights = max(1,mWeights);
647  int lWeights = weightSave.size();
648  weightSave.resize(mWeights);
649  weightLabelSave.resize(mWeights);
650  for (int i=lWeights; i<mWeights; ++i) weightLabelSave[i]="";
651  }
652  void setWeightLabel(int iWeight, string labelIn) {
653  if (iWeight >= 0 && iWeight < (int)weightLabelSave.size())
654  weightLabelSave[iWeight] = labelIn;
655  }
656 
657  // Set event weight, either for LHEF3 or for uncertainty bands.
658  void setWeight( double weightIn, int lhaStrategyIn) {
659  for (int i=0; i<nWeights(); ++i) weightSave[i] = weightIn;
660  if (nWeights() == 0) weightSave.push_back(weightIn);
661  lhaStrategySave = lhaStrategyIn;}
662 
663  // Apply weight modification (used for automated uncertainty variations).
664  void reWeight( int iWeight, double rwIn) {
665  if (iWeight >= 0 || iWeight < nWeights()) weightSave[iWeight] *= rwIn;}
666 
667  // Save merging weight (i.e. CKKW-L-type weight, summed O(\alpha_s) weight).
668  double weightCKKWLSave, weightFIRSTSave;
669 
670  // Set info on resolved processes.
671  void setIsResolved(bool isResIn) {isRes = isResIn;}
672 
673  // Set info on hard diffraction.
674  void setHardDiff( bool hasUnresBeamsIn = false, bool hasPomPsysIn = false,
675  bool isHardDiffAIn = false, bool isHardDiffBIn = false,
676  double xPomAIn = 0., double xPomBIn = 0., double tPomAIn = 0.,
677  double tPomBIn = 0.) { hasUnresBeams = hasUnresBeamsIn;
678  hasPomPsys = hasPomPsysIn; isHardDiffA = isHardDiffAIn;
679  isHardDiffB = isHardDiffBIn;
680  xPomA = xPomAIn; xPomB = xPomBIn;
681  tPomA = tPomAIn; tPomB = tPomBIn;}
682 
683  // Set information in hard diffractive events.
684  void setHasUnresolvedBeams(bool hasUnresBeamsIn)
685  {hasUnresBeams = hasUnresBeamsIn;}
686  void setHasPomPsystem(bool hasPomPsysIn) {hasPomPsys = hasPomPsysIn;}
687 
688  // Variables for weak shower setup.
689  vector<int> weakModes, weak2to2lines;
690  vector<Vec4> weakMomenta;
691  vector<pair<int, int> > weakDipoles;
692 
693 };
694 
695 //==========================================================================
696 
697 } // end namespace Pythia8
698 
699 #endif // Pythia8_Info_H