PYTHIA  8.313
Info.h
1 // Info.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 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 classes that keep track of generic event info.
7 // Info: contains information on the generation process and errors.
8 // Info: user interface to make parts of Info publicly available.
9 
10 #ifndef Pythia8_Info_H
11 #define Pythia8_Info_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/LHEF3.h"
15 #include "Pythia8/Logger.h"
16 #include "Pythia8/PythiaStdlib.h"
17 #include "Pythia8/SharedPointers.h"
18 #include "Pythia8/Weights.h"
19 
20 namespace Pythia8 {
21 
22 // Forward declaration of various classes.
23 class Settings;
24 class ParticleData;
25 class Rndm;
26 class BeamSetup;
27 class CoupSM;
28 class CoupSUSY;
29 class BeamParticle;
30 class PartonSystems;
31 class SigmaTotal;
32 class SigmaCombined;
33 class HadronWidths;
34 
35 // Forward declaration of HIInfo class.
36 class HIInfo;
37 
38 //==========================================================================
39 
40 // The Info class contains a mixed bag of information on the event
41 // generation activity, especially on the current subprocess properties,
42 // and on the number of errors encountered. This is used by the
43 // generation machinery.
44 
45 class Info {
46 
47 public:
48 
49  // Constructors.
50  Info() = default;
51  Info(bool) : Info() {}
52 
53  // Destructor for clean-up.
54  ~Info() {
55  if (hasOwnEventAttributes && eventAttributes != nullptr)
56  delete eventAttributes;
57  }
58 
59  // Assignment operator gives a shallow copy; no objects pointed to
60  // are copied.
61  Info & operator=(const Info &) = default;
62 
63  // Pointers to other class objects, carried piggyback by Info.
64 
65  // Set pointers to other class objects.
66  void setPtrs(Settings* settingsPtrIn, ParticleData* particleDataPtrIn,
67  Logger* loggerPtrIn, Rndm* rndmPtrIn, BeamSetup* beamSetupIn,
68  CoupSM* coupSMPtrIn, CoupSUSY* coupSUSYPtrIn,
69  PartonSystems* partonSystemsPtrIn, SigmaTotal* sigmaTotPtrIn,
70  SigmaCombined* sigmaCmbPtrIn, HadronWidths* hadronWidthsPtrIn,
71  WeightContainer* weightContainerPtrIn) {
72  settingsPtr = settingsPtrIn; particleDataPtr = particleDataPtrIn;
73  loggerPtr = loggerPtrIn; rndmPtr = rndmPtrIn; beamSetupPtr = beamSetupIn;
74  coupSMPtr = coupSMPtrIn; coupSUSYPtr = coupSUSYPtrIn;
75  partonSystemsPtr = partonSystemsPtrIn; sigmaTotPtr = sigmaTotPtrIn;
76  sigmaCmbPtr = sigmaCmbPtrIn; hadronWidthsPtr = hadronWidthsPtrIn;
77  weightContainerPtr = weightContainerPtrIn; }
78 
79  // Pointer to the settings database.
81 
82  // Pointer to the particle data table.
84 
85  // Pointer to the logger.
87 
88  // Pointer to the random number generator.
90 
91  // Pointers to beam configuration.
93 
94  // Pointers to Standard Model and Beyond SM couplings.
96  CoupSUSY* coupSUSYPtr{};
97 
98  // Pointer to information on subcollision parton locations.
100 
101  // Pointers to the total/elastic/diffractive cross sections.
103  SigmaCombined* sigmaCmbPtr{};
104 
105  // Pointer to the hadron widths data table.
107 
108  // Pointer to the UserHooks object set for the run.
109  UserHooksPtr userHooksPtr{};
110 
111  // Pointer to information about a HeavyIons run and the current event.
112  // (Is NULL if HeavyIons object is inactive.)
114 
115  WeightContainer* weightContainerPtr{};
116 
117  // Listing of most available information on current event.
118  void list() const;
119 
120  // Beam particles (in rest frame). CM energy of event.
121  int idA() const {return idASave;}
122  int idB() const {return idBSave;}
123  double pzA() const {return pzASave;}
124  double pzB() const {return pzBSave;}
125  double eA() const {return eASave;}
126  double eB() const {return eBSave;}
127  double mA() const {return mASave;}
128  double mB() const {return mBSave;}
129  double eCM() const {return eCMSave;}
130  double s() const {return sSave;}
131 
132  // Warnings from initialization.
133  bool tooLowPTmin() const {return lowPTmin;}
134 
135  // Process name and code, and the number of final-state particles.
136  string name() const {return nameSave;}
137  int code() const {return codeSave;}
138  int nFinal() const {return nFinalSave;}
139 
140  // Are beam particles resolved, with pdf's? Are they diffractive?
141  bool isResolved() const {return isRes;}
142  bool isDiffractiveA() const {return isDiffA;}
143  bool isDiffractiveB() const {return isDiffB;}
144  bool isDiffractiveC() const {return isDiffC;}
145  bool isNonDiffractive() const {return isND;}
146  bool isElastic() const {return (codeSave == 102);}
147  // Retained for backwards compatibility.
148  bool isMinBias() const {return isND;}
149 
150  // Information for Les Houches Accord and reading files.
151  bool isLHA() const {return isLH;}
152  bool atEndOfFile() const {return atEOF;}
153 
154  // For nondiffractive and Les Houches Accord identify hardest subprocess.
155  bool hasSub(int i = 0) const {return hasSubSave[i];}
156  string nameSub(int i = 0) const {return nameSubSave[i];}
157  int codeSub(int i = 0) const {return codeSubSave[i];}
158  int nFinalSub(int i = 0) const {return nFinalSubSave[i];}
159 
160  // Incoming parton flavours and x values.
161  int id1(int i = 0) const {return id1Save[i];}
162  int id2(int i = 0) const {return id2Save[i];}
163  double x1(int i = 0) const {return x1Save[i];}
164  double x2(int i = 0) const {return x2Save[i];}
165  double y(int i = 0) const {return 0.5 * log( x1Save[i] / x2Save[i]);}
166  double tau(int i = 0) const {return x1Save[i] * x2Save[i];}
167 
168  // Hard process flavours, x values, parton densities, couplings, Q2 scales.
169  int id1pdf(int i = 0) const {return id1pdfSave[i];}
170  int id2pdf(int i = 0) const {return id2pdfSave[i];}
171  double x1pdf(int i = 0) const {return x1pdfSave[i];}
172  double x2pdf(int i = 0) const {return x2pdfSave[i];}
173  double pdf1(int i = 0) const {return pdf1Save[i];}
174  double pdf2(int i = 0) const {return pdf2Save[i];}
175  double QFac(int i = 0) const {return sqrtpos(Q2FacSave[i]);}
176  double Q2Fac(int i = 0) const {return Q2FacSave[i];}
177  bool isValence1() const {return isVal1;}
178  bool isValence2() const {return isVal2;}
179  double alphaS(int i = 0) const {return alphaSSave[i];}
180  double alphaEM(int i = 0) const {return alphaEMSave[i];}
181  double QRen(int i = 0) const {return sqrtpos(Q2RenSave[i]);}
182  double Q2Ren(int i = 0) const {return Q2RenSave[i];}
183  double scalup(int i = 0) const {return scalupSave[i];}
184 
185  // DIS-specific kinematic variables.
186  double Q2DIS() const {return Q2DISSave;}
187  double WDIS() const {return WDISSave;}
188  double xDIS() const {return xDISSave;}
189  double yDIS() const {return yDISSave;}
190 
191  // Kinematics of photons from lepton beams.
192  double xGammaA() const {return x1GammaSave;}
193  double xGammaB() const {return x2GammaSave;}
194  double Q2GammaA() const {return Q2Gamma1Save;}
195  double Q2GammaB() const {return Q2Gamma2Save;}
196  double eCMsub() const {return eCMsubSave;}
197  double thetaScatLepA() const {return thetaLepton1;}
198  double thetaScatLepB() const {return thetaLepton2;}
199  double sHatNew() const {return sHatNewSave;}
200  int photonMode() const {return gammaModeEvent;}
201 
202  // Information on VMD state inside a photon.
203  bool isVMDstateA() const {return isVMDstateAEvent;}
204  bool isVMDstateB() const {return isVMDstateBEvent;}
205  int idVMDA() const {return idVMDASave;}
206  int idVMDB() const {return idVMDBSave;}
207  double mVMDA() const {return mVMDASave;}
208  double mVMDB() const {return mVMDBSave;}
209  double scaleVMDA() const {return scaleVMDASave;}
210  double scaleVMDB() const {return scaleVMDBSave;}
211 
212  // Mandelstam variables (notation as if subcollision).
213  double mHat(int i = 0) const {return sqrt(sH[i]);}
214  double sHat(int i = 0) const {return sH[i];}
215  double tHat(int i = 0) const {return tH[i];}
216  double uHat(int i = 0) const {return uH[i];}
217  double pTHat(int i = 0) const {return pTH[i];}
218  double pT2Hat(int i = 0) const {return pTH[i] * pTH[i];}
219  double m3Hat(int i = 0) const {return m3H[i];}
220  double m4Hat(int i = 0) const {return m4H[i];}
221  double thetaHat(int i = 0) const {return thetaH[i];}
222  double phiHat(int i = 0) const {return phiH[i];}
223 
224  // Weight of current event; normally 1, but used for Les Houches events
225  // or when reweighting phase space selection. Conversion from mb to pb
226  // for LHA strategy +-4. Uncertainty variations can be accessed by
227  // providing an index >= 1 (0 = no variations). Also cumulative sum.
228  double weight(int i = 0) const;
229  double weightSum() const;
230  double lhaStrategy() const {return lhaStrategySave;}
231 
232  // Further access to uncertainty weights: number and labels.
233  int nWeights() const {
234  return weightContainerPtr->weightsShowerPtr->getWeightsSize() +
235  weightContainerPtr->weightsFragmentation.getWeightsSize() - 1;}
236  string weightLabel(int iWgt) const;
237  int nWeightGroups() const {return weightContainerPtr->
238  weightsShowerPtr->nWeightGroups() + weightContainerPtr->
239  weightsFragmentation.nWeightGroups();}
240  string getGroupName(int iGN) const;
241  double getGroupWeight(int iGW) const;
242 
243  // Number of times other steps have been carried out.
244  int nISR() const {return nISRSave;}
245  int nFSRinProc() const {return nFSRinProcSave;}
246  int nFSRinRes() const {return nFSRinResSave;}
247 
248  // Maximum pT scales for MPI, ISR and FSR (in hard process).
249  double pTmaxMPI() const {return pTmaxMPISave;}
250  double pTmaxISR() const {return pTmaxISRSave;}
251  double pTmaxFSR() const {return pTmaxFSRSave;}
252 
253  // Current evolution scale (for UserHooks).
254  double pTnow() const {return pTnowSave;}
255 
256  // Impact parameter picture, global information
257  double a0MPI() const {return a0MPISave;}
258 
259  // Impact parameter picture, as set by hardest interaction.
260  double bMPI() const {return (bIsSet) ? bMPISave : 1.;}
261  double enhanceMPI() const {return (bIsSet) ? enhanceMPISave : 1.;}
262  double enhanceMPIavg() const {return (bIsSet) ? enhanceMPIavgSave : 1.;}
263  double eMPI(int i) const {return (bIsSet) ? eMPISave[i] : 1.;}
264  double bMPIold() const {return (bIsSet) ? bMPIoldSave : 1.;}
265  double enhanceMPIold() const {return (bIsSet) ? enhanceMPIoldSave : 1.;}
266  double enhanceMPIoldavg() const {return (bIsSet)
267  ? enhanceMPIoldavgSave : 1.;}
268 
269  // Number of multiparton interactions, with code and pT for them.
270  int nMPI() const {return nMPISave;}
271  int codeMPI(int i) const {return codeMPISave[i];}
272  double pTMPI(int i) const {return pTMPISave[i];}
273  int iAMPI(int i) const {return iAMPISave[i];}
274  int iBMPI(int i) const {return iBMPISave[i];}
275 
276  // Cross section estimate, optionally process by process.
277  vector<int> codesHard();
278 
279  // Name of the specified process.
280  string nameProc(int i = 0) const {
281  if (i == 0) return "sum";
282  auto itr = procNameM.find(i);
283  if (itr != procNameM.end()) return itr->second;
284  loggerPtr->ERROR_MSG("process code not found", to_string(i));
285  return "unknown process";}
286 
287  // The number of phase-space points tried.
288  long nTried(int i = 0) const {
289  if (i == 0) return nTry;
290  auto itr = nTryM.find(i);
291  if (itr != nTryM.end()) return itr->second;
292  loggerPtr->ERROR_MSG("process code not found", to_string(i));
293  return 0;}
294 
295  // The number of selected hard processes.
296  long nSelected(int i = 0) const {
297  if (i == 0) return nSel;
298  auto itr = nSelM.find(i);
299  if (itr != nSelM.end()) return itr->second;
300  loggerPtr->ERROR_MSG("process code not found", to_string(i));
301  return 0;}
302 
303  // The number of accepted events.
304  long nAccepted(int i = 0) const {
305  if (i == 0) return nAcc;
306  auto itr = nAccM.find(i);
307  if (itr != nAccM.end()) return itr->second;
308  loggerPtr->ERROR_MSG("process code not found", to_string(i));
309  return 0;}
310 
311  // The estimated cross-section in units of mb.
312  double sigmaGen(int i = 0) const {
313  if (i == 0) return sigGen;
314  auto itr = sigGenM.find(i);
315  if (itr != sigGenM.end()) return itr->second;
316  loggerPtr->ERROR_MSG("process code not found", to_string(i));
317  return 0;}
318 
319  // The uncertainty on the estimated cross-section in units of mb.
320  double sigmaErr(int i = 0) const {
321  if (i == 0) return sigErr;
322  auto itr = sigErrM.find(i);
323  if (itr != sigErrM.end()) return itr->second;
324  //loggerPtr->ERROR_MSG("process code not found", to_string(code));
325  return 0;}
326 
327  // Counters for number of loops in various places.
328  int getCounter( int i) const {return counters[i];}
329 
330  // Set or increase the value stored in a counter.
331  void setCounter( int i, int value = 0) {counters[i] = value;}
332  void addCounter( int i, int value = 1) {counters[i] += value;}
333 
334  // Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
335  void setTooLowPTmin(bool lowPTminIn) {lowPTmin = lowPTminIn;}
336 
337  // Set info on valence character of hard collision partons.
338  void setValence( bool isVal1In, bool isVal2In) {isVal1 = isVal1In;
339  isVal2 = isVal2In;}
340 
341  // Set and get some MPI/ISR/FSR properties needed for matching,
342  // i.e. mainly of internal relevance.
343  void hasHistory( bool hasHistoryIn) {hasHistorySave = hasHistoryIn;}
344  bool hasHistory() {return hasHistorySave;}
345  void zNowISR( double zNowIn) {zNowISRSave = zNowIn;}
346  double zNowISR() {return zNowISRSave;}
347  void pT2NowISR( double pT2NowIn) {pT2NowISRSave = pT2NowIn;}
348  double pT2NowISR() {return pT2NowISRSave;}
349 
350  // Return merging weight.
351  double mergingWeight(int i=0) const {
352  return weightContainerPtr->weightsMerging.getWeightsValue(i);}
353 
354  // Return the complete NLO weight.
355  double mergingWeightNLO(int i=0) const {
356  return weightContainerPtr->weightsMerging.getWeightsValue(i); }
357 
358  // Return an LHEF header
359  string header(const string &key) const {
360  if (headers.find(key) == headers.end()) return "";
361  else return headers.at(key);
362  }
363 
364  // Return a list of all header key names
365  vector<string> headerKeys() const;
366 
367  // Return the number of processes in the LHEF.
368  int nProcessesLHEF() const { return int(sigmaLHEFSave.size());}
369  // Return the cross section information read from LHEF.
370  double sigmaLHEF(int iProcess) const { return sigmaLHEFSave[iProcess];}
371 
372  // LHEF3 information: Public for easy access.
374 
375  // Save process information as read from init block of LHEF.
376  vector<double> sigmaLHEFSave;
377 
378  // Contents of the LHAinitrwgt tag
380 
381  // Contents of the LHAgenerator tags.
382  vector<LHAgenerator > *generators{};
383 
384  // A map of the LHAweightgroup tags, indexed by name.
385  map<string,LHAweightgroup > *weightgroups{};
386 
387  // A map of the LHAweight tags, indexed by name.
388  map<string,LHAweight > *init_weights{};
389 
390  // Store current-event Les Houches event tags.
392  map<string, string > *eventAttributes{};
393 
394  // The weights associated with this event, as given by the LHAwgt tags
395  map<string,double > *weights_detailed{};
396 
397  // The weights associated with this event, as given by the LHAweights tags
398  vector<double > *weights_compressed{};
399 
400  // Contents of the LHAscales tag.
402 
403  // Contents of the LHAweights tag (compressed format)
405 
406  // Contents of the LHArwgt tag (detailed format)
408 
409  // Vectorized version of LHArwgt tag, for easy and ordered access.
410  vector<double> weights_detailed_vector;
411 
412  // Value of the unit event weight read from a Les Houches event, necessary
413  // if additional weights in an otherwise unweighted input file are in
414  // relation to this number.
415  double eventWeightLHEF = {};
416 
417  // Set the LHEF3 objects read from the init and header blocks.
418  void setLHEF3InitInfo();
419  void setLHEF3InitInfo( int LHEFversionIn, LHAinitrwgt *initrwgtIn,
420  vector<LHAgenerator> *generatorsIn,
421  map<string,LHAweightgroup> *weightgroupsIn,
422  map<string,LHAweight> *init_weightsIn, string headerBlockIn );
423  // Set the LHEF3 objects read from the event block.
424  void setLHEF3EventInfo();
425  void setLHEF3EventInfo( map<string, string> *eventAttributesIn,
426  map<string, double > *weights_detailedIn,
427  vector<double > *weights_compressedIn,
428  LHAscales *scalesIn, LHAweights *weightsIn,
429  LHArwgt *rwgtIn, vector<double> weights_detailed_vecIn,
430  vector<string> weights_detailed_name_vecIn,
431  string eventCommentsIn, double eventWeightLHEFIn );
432 
433  // Retrieve events tag information.
434  string getEventAttribute(string key, bool doRemoveWhitespace = false) const;
435 
436  // Externally set event tag auxiliary information.
437  void setEventAttribute(string key, string value, bool doOverwrite = true) {
438  if (eventAttributes == nullptr) {
439  eventAttributes = new map<string,string>();
440  hasOwnEventAttributes = true;
441  }
442  map<string, string>::iterator it = eventAttributes->find(key);
443  if ( it != eventAttributes->end() && !doOverwrite ) return;
444  if ( it != eventAttributes->end() ) eventAttributes->erase(it);
445  eventAttributes->insert(make_pair(key,value));
446  }
447 
448  // Retrieve LHEF version
449  int LHEFversion() const { return LHEFversionSave; }
450 
451  // Retrieve initrwgt tag information.
452  unsigned int getInitrwgtSize() const;
453 
454  // Retrieve generator tag information.
455  unsigned int getGeneratorSize() const;
456  string getGeneratorValue(unsigned int n = 0) const;
457  string getGeneratorAttribute( unsigned int n, string key,
458  bool doRemoveWhitespace = false) const;
459 
460  // Retrieve rwgt tag information.
461  unsigned int getWeightsDetailedSize() const;
462  double getWeightsDetailedValue(string n) const;
463  string getWeightsDetailedAttribute(string n, string key,
464  bool doRemoveWhitespace = false) const;
465 
466  // Retrieve weights tag information.
467  unsigned int getWeightsCompressedSize() const;
468  double getWeightsCompressedValue(unsigned int n) const;
469  string getWeightsCompressedAttribute(string key,
470  bool doRemoveWhitespace = false) const;
471 
472  // Retrieve scales tag information.
473  string getScalesValue(bool doRemoveWhitespace = false) const;
474  double getScalesAttribute(string key) const;
475 
476  // Retrieve complete header block and event comments
477  // Retrieve scales tag information.
478  string getHeaderBlock() const { return headerBlock;}
479  string getEventComments() const { return eventComments;}
480 
481  // Set LHEF headers
482  void setHeader(const string &key, const string &val)
483  { headers[key] = val; }
484 
485  // Set abort in parton level.
486  void setAbortPartonLevel(bool abortIn) {abortPartonLevel = abortIn;}
487  bool getAbortPartonLevel() const {return abortPartonLevel;}
488 
489  // Get information on hard diffractive events.
490  bool hasUnresolvedBeams() const {return hasUnresBeams;}
491  bool hasPomPsystem() const {return hasPomPsys;}
492  bool isHardDiffractive() const {return isHardDiffA || isHardDiffB;}
493  bool isHardDiffractiveA() const {return isHardDiffA;}
494  bool isHardDiffractiveB() const {return isHardDiffB;}
495  double xPomeronA() const {return xPomA;}
496  double xPomeronB() const {return xPomB;}
497  double tPomeronA() const {return tPomA;}
498  double tPomeronB() const {return tPomB;}
499 
500  // History information needed to setup the weak shower for 2 -> n.
501  vector<int> getWeakModes() const {return weakModes;}
502  vector<pair<int,int> > getWeakDipoles() const {return weakDipoles;}
503  vector<Vec4> getWeakMomenta() const {return weakMomenta;}
504  vector<int> getWeak2to2lines() const {return weak2to2lines;}
505  void setWeakModes(vector<int> weakModesIn) {weakModes = weakModesIn;}
506  void setWeakDipoles(vector<pair<int,int> > weakDipolesIn)
507  {weakDipoles = weakDipolesIn;}
508  void setWeakMomenta(vector<Vec4> weakMomentaIn)
509  {weakMomenta = weakMomentaIn;}
510  void setWeak2to2lines(vector<int> weak2to2linesIn)
511  {weak2to2lines = weak2to2linesIn;}
512 
513  // Check if onia is included in showers (used for MPI process setup).
514  void setOniumShower(bool oniumShowerIn) {oniumShower = oniumShowerIn;}
515  bool getOniumShower() const {return oniumShower;}
516 
517  // From here on what used to be the private part of the class.
518 
519  // Allow conversion from mb to pb.
520  static const double CONVERTMB2PB;
521 
522  // Store common beam quantities.
523  int idASave{}, idBSave{};
524  double pzASave{}, eASave{}, mASave{}, pzBSave{}, eBSave{}, mBSave{},
525  eCMSave{}, sSave{};
526 
527  // Store initialization information.
528  bool lowPTmin;
529 
530  // Store common integrated cross section quantities.
531  long nTry{}, nSel{}, nAcc{};
532  double sigGen{}, sigErr{}, wtAccSum{};
533  map<int, string> procNameM;
534  map<int, long> nTryM, nSelM, nAccM;
535  map<int, double> sigGenM, sigErrM;
536  int lhaStrategySave{};
537 
538  // Store common MPI information.
539  double a0MPISave;
540 
541  // Store current-event quantities.
542  bool isRes{}, isDiffA{}, isDiffB{}, isDiffC{}, isND{}, isLH{},
543  hasSubSave[4], bIsSet{}, evolIsSet{}, atEOF{}, isVal1{}, isVal2{},
544  hasHistorySave{}, abortPartonLevel{}, isHardDiffA{}, isHardDiffB{},
545  hasUnresBeams{}, hasPomPsys{};
546  int codeSave{}, codeSubSave[4], nFinalSave{}, nFinalSubSave[4], nTotal{},
547  id1Save[4], id2Save[4], id1pdfSave[4], id2pdfSave[4], nMPISave{},
548  nISRSave{}, nFSRinProcSave{}, nFSRinResSave{};
549  double x1Save[4], x2Save[4], x1pdfSave[4], x2pdfSave[4], pdf1Save[4],
550  pdf2Save[4], Q2FacSave[4], alphaEMSave[4], alphaSSave[4],
551  Q2RenSave[4], scalupSave[4], sH[4], tH[4], uH[4], pTH[4], m3H[4],
552  m4H[4], thetaH[4], phiH[4], bMPISave{}, enhanceMPISave{},
553  enhanceMPIavgSave{}, bMPIoldSave{}, enhanceMPIoldSave{},
554  enhanceMPIoldavgSave{}, pTmaxMPISave{}, pTmaxISRSave{},
555  pTmaxFSRSave{}, pTnowSave{}, zNowISRSave{}, pT2NowISRSave{}, xPomA{},
556  xPomB{}, tPomA{}, tPomB{};
557  string nameSave{}, nameSubSave[4];
558  vector<int> codeMPISave, iAMPISave, iBMPISave;
559  vector<double> pTMPISave, eMPISave;
560 
561  // DIS-specific kinematic variables.
562  double Q2DISSave{}, WDISSave{}, xDISSave{}, yDISSave{};
563 
564  // Variables related to photon kinematics.
565  bool isVMDstateAEvent{}, isVMDstateBEvent{};
566  int gammaModeEvent{}, idVMDASave{}, idVMDBSave{};
567  double x1GammaSave{}, x2GammaSave{}, Q2Gamma1Save{}, Q2Gamma2Save{},
568  eCMsubSave{}, thetaLepton1{}, thetaLepton2{}, sHatNewSave{},
569  mVMDASave{}, mVMDBSave{}, scaleVMDASave{}, scaleVMDBSave{};
570 
571  // Vector of various loop counters.
572  int counters[50];
573 
574  // Map for LHEF headers.
575  map<string, string> headers;
576 
577  // Strings for complete header block and event comments.
578  string headerBlock{}, eventComments{};
579 
580  // Set info on the two incoming beams: only from Pythia class.
581  void setBeamIDs( int idAin, int idBin) { idASave = idAin; idBSave = idBin;}
582  void setBeamA( int idAin, double pzAin, double eAin, double mAin) {
583  idASave = idAin; pzASave = pzAin; eASave = eAin; mASave = mAin;}
584  void setBeamB( int idBin, double pzBin, double eBin, double mBin) {
585  idBSave = idBin; pzBSave = pzBin; eBSave = eBin; mBSave = mBin;}
586  void setECM( double eCMin) {eCMSave = eCMin; sSave = eCMSave * eCMSave;}
587 
588  // Set info on DIS-specific kinematic variables.
589  void setDISKinematics(double Q2In, double WIn, double xIn, double yIn) {
590  Q2DISSave = Q2In; WDISSave = WIn; xDISSave = xIn; yDISSave = yIn;}
591  // Set info related to gamma+gamma subcollision.
592  void setX1Gamma( double x1GammaIn) { x1GammaSave = x1GammaIn; }
593  void setX2Gamma( double x2GammaIn) { x2GammaSave = x2GammaIn; }
594  void setQ2Gamma1( double Q2gammaIn) { Q2Gamma1Save = Q2gammaIn; }
595  void setQ2Gamma2( double Q2gammaIn) { Q2Gamma2Save = Q2gammaIn; }
596  void setTheta1( double theta1In) { thetaLepton1 = theta1In; }
597  void setTheta2( double theta2In) { thetaLepton2 = theta2In; }
598  void setECMsub( double eCMsubIn) { eCMsubSave = eCMsubIn; }
599  void setsHatNew( double sHatNewIn) { sHatNewSave = sHatNewIn; }
600  void setGammaMode( double gammaModeIn) { gammaModeEvent = gammaModeIn; }
601  // Set info on VMD state.
602  void setVMDstateA(bool isVMDAIn, int idAIn, double mAIn, double scaleAIn)
603  {isVMDstateAEvent = isVMDAIn; idVMDASave = idAIn; mVMDASave = mAIn;
604  scaleVMDASave = scaleAIn;}
605  void setVMDstateB(bool isVMDBIn, int idBIn, double mBIn, double scaleBIn)
606  {isVMDstateBEvent = isVMDBIn; idVMDBSave = idBIn; mVMDBSave = mBIn;
607  scaleVMDBSave = scaleBIn;}
608 
609  // Reset info for current event: only from Pythia class.
610  void clear() {
611  isRes = isDiffA = isDiffB = isDiffC = isND = isLH = bIsSet
612  = evolIsSet = atEOF = isVal1 = isVal2 = hasHistorySave
613  = isHardDiffA = isHardDiffB = hasUnresBeams = hasPomPsys = false;
614  codeSave = nFinalSave = nTotal = nMPISave = nISRSave = nFSRinProcSave
615  = nFSRinResSave = 0;
616  bMPISave = enhanceMPISave = enhanceMPIavgSave = bMPIoldSave
617  = enhanceMPIoldSave = enhanceMPIoldavgSave = 1.;
618  pTmaxMPISave = pTmaxISRSave = pTmaxFSRSave = pTnowSave = zNowISRSave
619  = pT2NowISRSave = 0.;
620  nameSave = " ";
621  for (int i = 0; i < 4; ++i) {
622  hasSubSave[i] = false;
623  codeSubSave[i] = nFinalSubSave[i] = id1pdfSave[i] = id2pdfSave[i]
624  = id1Save[i] = id2Save[i] = 0;
625  x1pdfSave[i] = x2pdfSave[i] = pdf1Save[i] = pdf2Save[i]
626  = Q2FacSave[i] = alphaEMSave[i] = alphaSSave[i] = Q2RenSave[i]
627  = scalupSave[i] = x1Save[i] = x2Save[i] = sH[i] = tH[i] = uH[i]
628  = pTH[i] = m3H[i] = m4H[i] = thetaH[i] = phiH[i] = 0.;
629  nameSubSave[i] = " ";
630  }
631  codeMPISave.resize(0); iAMPISave.resize(0); iBMPISave.resize(0);
632  pTMPISave.resize(0); eMPISave.resize(0); setHardDiff();
633  weightContainerPtr->clear();
634  }
635 
636  // Reset info arrays only for parton and hadron level.
637  int sizeMPIarrays() const { return codeMPISave.size(); }
638  void resizeMPIarrays(int newSize) { codeMPISave.resize(newSize);
639  iAMPISave.resize(newSize); iBMPISave.resize(newSize);
640  pTMPISave.resize(newSize); eMPISave.resize(newSize); }
641 
642  // Set info on the (sub)process: from ProcessLevel, ProcessContainer or
643  // MultipartonInteractions classes.
644  void setType( string nameIn, int codeIn, int nFinalIn,
645  bool isNonDiffIn = false, bool isResolvedIn = true,
646  bool isDiffractiveAin = false, bool isDiffractiveBin = false,
647  bool isDiffractiveCin = false, bool isLHAin = false) {
648  nameSave = nameIn; codeSave = codeIn; nFinalSave = nFinalIn;
649  isND = isNonDiffIn; isRes = isResolvedIn; isDiffA = isDiffractiveAin;
650  isDiffB = isDiffractiveBin; isDiffC = isDiffractiveCin; isLH = isLHAin;
651  nTotal = 2 + nFinalSave; bIsSet = false; hasSubSave[0] = false;
652  nameSubSave[0] = " "; codeSubSave[0] = 0; nFinalSubSave[0] = 0;
653  evolIsSet = false;}
654  void setSubType( int iDS, string nameSubIn, int codeSubIn,
655  int nFinalSubIn) { hasSubSave[iDS] = true; nameSubSave[iDS] = nameSubIn;
656  codeSubSave[iDS] = codeSubIn; nFinalSubSave[iDS] = nFinalSubIn;}
657  void setPDFalpha( int iDS, int id1pdfIn, int id2pdfIn, double x1pdfIn,
658  double x2pdfIn, double pdf1In, double pdf2In, double Q2FacIn,
659  double alphaEMIn, double alphaSIn, double Q2RenIn, double scalupIn) {
660  id1pdfSave[iDS] = id1pdfIn; id2pdfSave[iDS] = id2pdfIn;
661  x1pdfSave[iDS] = x1pdfIn; x2pdfSave[iDS] = x2pdfIn;
662  pdf1Save[iDS] = pdf1In; pdf2Save[iDS] = pdf2In; Q2FacSave[iDS] = Q2FacIn;
663  alphaEMSave[iDS] = alphaEMIn; alphaSSave[iDS] = alphaSIn;
664  Q2RenSave[iDS] = Q2RenIn; scalupSave[iDS] = scalupIn;}
665  void setScalup( int iDS, double scalupIn) {scalupSave[iDS] = scalupIn;}
666  void setKin( int iDS, int id1In, int id2In, double x1In, double x2In,
667  double sHatIn, double tHatIn, double uHatIn, double pTHatIn,
668  double m3HatIn, double m4HatIn, double thetaHatIn, double phiHatIn) {
669  id1Save[iDS] = id1In; id2Save[iDS] = id2In; x1Save[iDS] = x1In;
670  x2Save[iDS] = x2In; sH[iDS] = sHatIn; tH[iDS] = tHatIn;
671  uH[iDS] = uHatIn; pTH[iDS] = pTHatIn; m3H[iDS] = m3HatIn;
672  m4H[iDS] = m4HatIn; thetaH[iDS] = thetaHatIn; phiH[iDS] = phiHatIn;}
673  void setTypeMPI( int codeMPIIn, double pTMPIIn, int iAMPIIn = 0,
674  int iBMPIIn = 0, double eMPIIn = 1.) {codeMPISave.push_back(codeMPIIn);
675  pTMPISave.push_back(pTMPIIn); iAMPISave.push_back(iAMPIIn);
676  iBMPISave.push_back(iBMPIIn); eMPISave.push_back(eMPIIn); }
677 
678  // Set info on cross section: from ProcessLevel (reset in Pythia).
679  void sigmaReset() { nTry = nSel = nAcc = 0; sigGen = sigErr = wtAccSum = 0.;
680  procNameM.clear(); nTryM.clear(); nSelM.clear(); nAccM.clear();
681  sigGenM.clear(); sigErrM.clear();}
682  void setSigma( int i, string procNameIn, long nTryIn, long nSelIn,
683  long nAccIn, double sigGenIn, double sigErrIn, double wtAccSumIn) {
684  if (i == 0) {nTry = nTryIn; nSel = nSelIn; nAcc = nAccIn;
685  sigGen = sigGenIn; sigErr = sigErrIn; wtAccSum = wtAccSumIn;}
686  else {procNameM[i] = procNameIn; nTryM[i] = nTryIn; nSelM[i] = nSelIn;
687  nAccM[i] = nAccIn; sigGenM[i] = sigGenIn; sigErrM[i] = sigErrIn;} }
688  void addSigma( int i, long nTryIn, long nSelIn, long nAccIn, double sigGenIn,
689  double sigErrIn) { nTryM[i] += nTryIn; nSelM[i] += nSelIn;
690  nAccM[i] += nAccIn; sigGenM[i] += sigGenIn;
691  sigErrM[i] = sqrtpos(sigErrM[i]*sigErrM[i] + sigErrIn*sigErrIn); }
692 
693  // Set info on impact parameter: from PartonLevel.
694  void setImpact( double bMPIIn, double enhanceMPIIn, double enhanceMPIavgIn,
695  bool bIsSetIn = true, bool pushBack = false) {
696  if (pushBack) {bMPIoldSave = bMPISave; enhanceMPIoldSave = enhanceMPISave;
697  enhanceMPIoldavgSave = enhanceMPIavgSave;}
698  bMPISave = bMPIIn; enhanceMPISave = eMPISave[0] = enhanceMPIIn,
699  enhanceMPIavgSave = enhanceMPIavgIn, bIsSet = bIsSetIn;}
700 
701  // Set info on pTmax scales and number of evolution steps: from PartonLevel.
702  void setPartEvolved( int nMPIIn, int nISRIn) {
703  nMPISave = nMPIIn; nISRSave = nISRIn;}
704  void setEvolution( double pTmaxMPIIn, double pTmaxISRIn, double pTmaxFSRIn,
705  int nMPIIn, int nISRIn, int nFSRinProcIn, int nFSRinResIn) {
706  pTmaxMPISave = pTmaxMPIIn; pTmaxISRSave = pTmaxISRIn;
707  pTmaxFSRSave = pTmaxFSRIn; nMPISave = nMPIIn; nISRSave = nISRIn;
708  nFSRinProcSave = nFSRinProcIn; nFSRinResSave = nFSRinResIn;
709  evolIsSet = true;}
710 
711  // Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
712  void setPTnow( double pTnowIn) {pTnowSave = pTnowIn;}
713 
714  // Set a0 from MultipartonInteractions.
715  void seta0MPI(double a0MPIin) {a0MPISave = a0MPIin;}
716 
717  // Set info whether reading of Les Houches Accord file at end.
718  void setEndOfFile( bool atEOFin) {atEOF = atEOFin;}
719 
720  // Set event weight, either for LHEF3 or for uncertainty bands. If weight
721  // has units (lhaStrategy 4 or -4), input in mb
722  void setWeight( double weightIn, int lhaStrategyIn) {
723  for (int i = 0; i < nWeights(); ++i)
724  weightContainerPtr->weightsShowerPtr->setValueByIndex(i,1.);
725  // Nominal weight in weightContainer saved in pb for lhaStrategy +-4
726  weightContainerPtr->setWeightNominal(
727  abs(lhaStrategyIn) == 4 ? CONVERTMB2PB*weightIn : weightIn);
728  lhaStrategySave = lhaStrategyIn;}
729  //
730  // Set info on resolved processes.
731  void setIsResolved(bool isResIn) {isRes = isResIn;}
732 
733  // Set info on hard diffraction.
734  void setHardDiff(bool hasUnresBeamsIn = false, bool hasPomPsysIn = false,
735  bool isHardDiffAIn = false, bool isHardDiffBIn = false,
736  double xPomAIn = 0., double xPomBIn = 0., double tPomAIn = 0.,
737  double tPomBIn = 0.) { hasUnresBeams = hasUnresBeamsIn;
738  hasPomPsys = hasPomPsysIn; isHardDiffA = isHardDiffAIn;
739  isHardDiffB = isHardDiffBIn;
740  xPomA = xPomAIn; xPomB = xPomBIn;
741  tPomA = tPomAIn; tPomB = tPomBIn;}
742 
743  // Move process information to an another diffractive system.
744  void reassignDiffSystem(int iDSold, int iDSnew);
745 
746  // Set information in hard diffractive events.
747  void setHasUnresolvedBeams(bool hasUnresBeamsIn)
748  {hasUnresBeams = hasUnresBeamsIn;}
749  void setHasPomPsystem(bool hasPomPsysIn) {hasPomPsys = hasPomPsysIn;}
750 
751  // Variables for weak and onia shower setup.
752  vector<int> weakModes, weak2to2lines;
753  vector<Vec4> weakMomenta;
754  vector<pair<int, int> > weakDipoles;
755  bool oniumShower{false};
756 
757  int numberOfWeights() const {
758  return weightContainerPtr->numberOfWeights(); }
759  double weightValueByIndex(int key=0) const {
760  return weightContainerPtr->weightValueByIndex(key); }
761  string weightNameByIndex(int key=0) const {
762  return weightContainerPtr->weightNameByIndex(key); }
763  vector<double> weightValueVector() const {
764  return weightContainerPtr->weightValueVector(); }
765  vector<string> weightNameVector() const {
766  return weightContainerPtr->weightNameVector(); }
767 
768 };
769 
770 //==========================================================================
771 
772 } // end namespace Pythia8
773 
774 #endif // Pythia8_Info_H
Info()=default
Constructors.
SigmaTotal * sigmaTotPtr
Pointers to the total/elastic/diffractive cross sections.
Definition: Info.h:102
void list() const
Listing of most available information on current event.
Definition: Info.cc:32
string name() const
Process name and code, and the number of final-state particles.
Definition: Info.h:136
void setPartEvolved(int nMPIIn, int nISRIn)
Set info on pTmax scales and number of evolution steps: from PartonLevel.
Definition: Info.h:702
BeamSetup * beamSetupPtr
Pointers to beam configuration.
Definition: Info.h:92
Rndm * rndmPtr
Pointer to the random number generator.
Definition: Info.h:89
bool hasUnresolvedBeams() const
Get information on hard diffractive events.
Definition: Info.h:490
double Q2DISSave
DIS-specific kinematic variables.
Definition: Info.h:562
UserHooksPtr userHooksPtr
Pointer to the UserHooks object set for the run.
Definition: Info.h:109
bool isVMDstateAEvent
Variables related to photon kinematics.
Definition: Info.h:565
string getScalesValue(bool doRemoveWhitespace=false) const
Retrieve scales tag information.
Definition: Info.cc:423
Definition: Info.h:45
double bMPI() const
Impact parameter picture, as set by hardest interaction.
Definition: Info.h:260
Definition: BeamSetup.h:33
bool hasOwnEventAttributes
Store current-event Les Houches event tags.
Definition: Info.h:391
void seta0MPI(double a0MPIin)
Set a0 from MultipartonInteractions.
Definition: Info.h:715
Definition: Weights.h:433
WeightsFragmentation weightsFragmentation
Fragmentation weights.
Definition: Weights.h:460
Logger * loggerPtr
Pointer to the logger.
Definition: Info.h:86
Definition: SigmaLowEnergy.h:135
vector< double > * weights_compressed
The weights associated with this event, as given by the LHAweights tags.
Definition: Info.h:398
double pTnow() const
Current evolution scale (for UserHooks).
Definition: Info.h:254
double pTmaxMPI() const
Maximum pT scales for MPI, ISR and FSR (in hard process).
Definition: Info.h:249
int counters[50]
Vector of various loop counters.
Definition: Info.h:572
Info & operator=(const Info &)=default
int idASave
Store common beam quantities.
Definition: Info.h:523
void setCounter(int i, int value=0)
Set or increase the value stored in a counter.
Definition: Info.h:331
CoupSM * coupSMPtr
Pointers to Standard Model and Beyond SM couplings.
Definition: Info.h:95
Definition: SigmaTotal.h:141
int nMPI() const
Number of multiparton interactions, with code and pT for them.
Definition: Info.h:270
vector< double > weights_detailed_vector
Vectorized version of LHArwgt tag, for easy and ordered access.
Definition: Info.h:410
bool isVMDstateA() const
Information on VMD state inside a photon.
Definition: Info.h:203
double mergingWeightNLO(int i=0) const
Return the complete NLO weight.
Definition: Info.h:355
map< string, string > headers
Map for LHEF headers.
Definition: Info.h:575
vector< int > codesHard()
Cross section estimate, optionally process by process.
Definition: Info.cc:212
double getWeightsValue(int iPos) const override
Modified weight getter to include first order weight.
Definition: Weights.h:258
double mHat(int i=0) const
Mandelstam variables (notation as if subcollision).
Definition: Info.h:213
LHAweights * weights
Contents of the LHAweights tag (compressed format)
Definition: Info.h:404
long nAccepted(int i=0) const
The number of accepted events.
Definition: Info.h:304
double sigmaErr(int i=0) const
The uncertainty on the estimated cross-section in units of mb.
Definition: Info.h:320
void reassignDiffSystem(int iDSold, int iDSnew)
Move process information to an another diffractive system.
Definition: Info.cc:453
Collect different scales relevant for an event.
Definition: LHEF3.h:281
map< string, double > * weights_detailed
The weights associated with this event, as given by the LHAwgt tags.
Definition: Info.h:395
The LHAweights struct represents the information in a weights tag.
Definition: LHEF3.h:245
Definition: Logger.h:23
string header(const string &key) const
Return an LHEF header.
Definition: Info.h:359
vector< LHAgenerator > * generators
Contents of the LHAgenerator tags.
Definition: Info.h:382
unsigned int getGeneratorSize() const
Retrieve generator tag information.
Definition: Info.cc:331
unsigned int getInitrwgtSize() const
Retrieve initrwgt tag information.
Definition: Info.cc:322
int id1pdf(int i=0) const
Hard process flavours, x values, parton densities, couplings, Q2 scales.
Definition: Info.h:169
bool hasSub(int i=0) const
For nondiffractive and Les Houches Accord identify hardest subprocess.
Definition: Info.h:155
void setWeight(double weightIn, int lhaStrategyIn)
Definition: Info.h:722
int LHEFversion() const
Retrieve LHEF version.
Definition: Info.h:449
The LHArwgt assigns a group-name to a set of LHAwgt objects.
Definition: LHEF3.h:473
void setOniumShower(bool oniumShowerIn)
Check if onia is included in showers (used for MPI process setup).
Definition: Info.h:514
long nSelected(int i=0) const
The number of selected hard processes.
Definition: Info.h:296
The LHAinitrwgt assigns a group-name to a set of LHAweightgroup objects.
Definition: LHEF3.h:511
void setEventAttribute(string key, string value, bool doOverwrite=true)
Externally set event tag auxiliary information.
Definition: Info.h:437
unsigned int getWeightsDetailedSize() const
Retrieve rwgt tag information.
Definition: Info.cc:362
int id1(int i=0) const
Incoming parton flavours and x values.
Definition: Info.h:161
LHAscales * scales
Contents of the LHAscales tag.
Definition: Info.h:401
Definition: Basics.h:388
bool isLHA() const
Information for Les Houches Accord and reading files.
Definition: Info.h:151
double a0MPI() const
Impact parameter picture, global information.
Definition: Info.h:257
void setVMDstateA(bool isVMDAIn, int idAIn, double mAIn, double scaleAIn)
Set info on VMD state.
Definition: Info.h:602
double mergingWeight(int i=0) const
Return merging weight.
Definition: Info.h:351
void setLHEF3EventInfo()
Set the LHEF3 objects read from the event block.
Definition: Info.cc:266
void setEndOfFile(bool atEOFin)
Set info whether reading of Les Houches Accord file at end.
Definition: Info.h:718
int nISR() const
Number of times other steps have been carried out.
Definition: Info.h:244
double a0MPISave
Store common MPI information.
Definition: Info.h:539
void setBeamIDs(int idAin, int idBin)
Set info on the two incoming beams: only from Pythia class.
Definition: Info.h:581
void setHasUnresolvedBeams(bool hasUnresBeamsIn)
Set information in hard diffractive events.
Definition: Info.h:747
long nTried(int i=0) const
The number of phase-space points tried.
Definition: Info.h:288
int numberOfWeights()
Functions to retrieve information stored in the subcategory members.
Definition: Weights.cc:983
HIInfo * hiInfo
Definition: Info.h:113
Definition: HadronWidths.h:22
void setWeightNominal(double weightNow)
The WeightContainer class.
Definition: Weights.cc:964
void clear()
Reset info for current event: only from Pythia class.
Definition: Info.h:610
double xGammaA() const
Kinematics of photons from lepton beams.
Definition: Info.h:192
string headerBlock
Strings for complete header block and event comments.
Definition: Info.h:578
void hasHistory(bool hasHistoryIn)
Definition: Info.h:343
map< string, LHAweight > * init_weights
A map of the LHAweight tags, indexed by name.
Definition: Info.h:388
map< string, LHAweightgroup > * weightgroups
A map of the LHAweightgroup tags, indexed by name.
Definition: Info.h:385
vector< string > weightNameVector()
Definition: Weights.cc:1027
vector< int > getWeakModes() const
History information needed to setup the weak shower for 2 -> n.
Definition: Info.h:501
bool lowPTmin
Store initialization information.
Definition: Info.h:528
void setHardDiff(bool hasUnresBeamsIn=false, bool hasPomPsysIn=false, bool isHardDiffAIn=false, bool isHardDiffBIn=false, double xPomAIn=0., double xPomBIn=0., double tPomAIn=0., double tPomBIn=0.)
Set info on hard diffraction.
Definition: Info.h:734
vector< string > headerKeys() const
Return a list of all header key names.
Definition: Info.cc:224
~Info()
Destructor for clean-up.
Definition: Info.h:54
void setValueByIndex(int iPos, double val)
Functions to set values of weights.
Definition: Weights.h:96
int nProcessesLHEF() const
Return the number of processes in the LHEF.
Definition: Info.h:368
double sigmaGen(int i=0) const
The estimated cross-section in units of mb.
Definition: Info.h:312
Definition: StandardModel.h:135
void setX1Gamma(double x1GammaIn)
Set info related to gamma+gamma subcollision.
Definition: Info.h:592
unsigned int getWeightsCompressedSize() const
Retrieve weights tag information.
Definition: Info.cc:394
void clear()
Reset all members to default stage.
Definition: Weights.cc:1051
bool isRes
Store current-event quantities.
Definition: Info.h:542
void setValence(bool isVal1In, bool isVal2In)
Set info on valence character of hard collision partons.
Definition: Info.h:338
long nTry
Store common integrated cross section quantities.
Definition: Info.h:531
int nWeights() const
Further access to uncertainty weights: number and labels.
Definition: Info.h:233
double sigmaLHEF(int iProcess) const
Return the cross section information read from LHEF.
Definition: Info.h:370
bool isResolved() const
Are beam particles resolved, with pdf&#39;s? Are they diffractive?
Definition: Info.h:141
vector< int > weakModes
Variables for weak and onia shower setup.
Definition: Info.h:752
ParticleData * particleDataPtr
Pointer to the particle data table.
Definition: Info.h:83
PartonSystems * partonSystemsPtr
Pointer to information on subcollision parton locations.
Definition: Info.h:99
void setIsResolved(bool isResIn)
Set info on resolved processes.
Definition: Info.h:731
string getHeaderBlock() const
Definition: Info.h:478
double eventWeightLHEF
Definition: Info.h:415
void setDISKinematics(double Q2In, double WIn, double xIn, double yIn)
Set info on DIS-specific kinematic variables.
Definition: Info.h:589
LHArwgt * rwgt
Contents of the LHArwgt tag (detailed format)
Definition: Info.h:407
HadronWidths * hadronWidthsPtr
Pointer to the hadron widths data table.
Definition: Info.h:106
void setImpact(double bMPIIn, double enhanceMPIIn, double enhanceMPIavgIn, bool bIsSetIn=true, bool pushBack=false)
Set info on impact parameter: from PartonLevel.
Definition: Info.h:694
static const double CONVERTMB2PB
From here on what used to be the private part of the class.
Definition: Info.h:520
The PartonSystems class describes the whole set of subcollisions.
Definition: PartonSystems.h:42
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
void setTooLowPTmin(bool lowPTminIn)
Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
Definition: Info.h:335
double Q2DIS() const
DIS-specific kinematic variables.
Definition: Info.h:186
LHAinitrwgt * initrwgt
Contents of the LHAinitrwgt tag.
Definition: Info.h:379
void setAbortPartonLevel(bool abortIn)
Set abort in parton level.
Definition: Info.h:486
vector< double > weightValueVector()
Definition: Weights.cc:1000
bool isMinBias() const
Retained for backwards compatibility.
Definition: Info.h:148
WeightsMerging weightsMerging
Merging weights.
Definition: Weights.h:457
void sigmaReset()
Set info on cross section: from ProcessLevel (reset in Pythia).
Definition: Info.h:679
void setPTnow(double pTnowIn)
Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
Definition: Info.h:712
void setPtrs(Settings *settingsPtrIn, ParticleData *particleDataPtrIn, Logger *loggerPtrIn, Rndm *rndmPtrIn, BeamSetup *beamSetupIn, CoupSM *coupSMPtrIn, CoupSUSY *coupSUSYPtrIn, PartonSystems *partonSystemsPtrIn, SigmaTotal *sigmaTotPtrIn, SigmaCombined *sigmaCmbPtrIn, HadronWidths *hadronWidthsPtrIn, WeightContainer *weightContainerPtrIn)
Pointers to other class objects, carried piggyback by Info.
Definition: Info.h:66
string getEventAttribute(string key, bool doRemoveWhitespace=false) const
Retrieve events tag information.
Definition: Info.cc:307
string nameProc(int i=0) const
Name of the specified process.
Definition: Info.h:280
int sizeMPIarrays() const
Reset info arrays only for parton and hadron level.
Definition: Info.h:637
bool tooLowPTmin() const
Warnings from initialization.
Definition: Info.h:133
vector< double > sigmaLHEFSave
Save process information as read from init block of LHEF.
Definition: Info.h:376
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
void setLHEF3InitInfo()
Set the LHEF3 objects read from the init and header blocks.
Definition: Info.cc:235
int getCounter(int i) const
Counters for number of loops in various places.
Definition: Info.h:328
Settings * settingsPtr
Pointer to the settings database.
Definition: Info.h:80
WeightsShower * weightsShowerPtr
Definition: Weights.h:453
void setHeader(const string &key, const string &val)
Set LHEF headers.
Definition: Info.h:482
int LHEFversionSave
LHEF3 information: Public for easy access.
Definition: Info.h:373
Definition: HIInfo.h:27
void setType(string nameIn, int codeIn, int nFinalIn, bool isNonDiffIn=false, bool isResolvedIn=true, bool isDiffractiveAin=false, bool isDiffractiveBin=false, bool isDiffractiveCin=false, bool isLHAin=false)
Definition: Info.h:644
Definition: SusyCouplings.h:27
double weight(int i=0) const
Event weights and accumulated weight.
Definition: Info.cc:162
Definition: Settings.h:196
double sqrtpos(const double &x)
Avoid problem with negative square root argument (from roundoff).
Definition: PythiaStdlib.h:191
int idA() const
Beam particles (in rest frame). CM energy of event.
Definition: Info.h:121