Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-08 08:30:41

0001 // Info.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2025 Torbjorn Sjostrand.
0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
0004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
0005 
0006 // This file contains classes that keep track of generic event info.
0007 // Info: contains information on the generation process and errors.
0008 // Info: user interface to make parts of Info publicly available.
0009 
0010 #ifndef Pythia8_Info_H
0011 #define Pythia8_Info_H
0012 
0013 #include "Pythia8/Basics.h"
0014 #include "Pythia8/LHEF3.h"
0015 #include "Pythia8/Logger.h"
0016 #include "Pythia8/PythiaStdlib.h"
0017 #include "Pythia8/SharedPointers.h"
0018 #include "Pythia8/Weights.h"
0019 
0020 namespace Pythia8 {
0021 
0022 // Forward declaration of various classes.
0023 class Settings;
0024 class ParticleData;
0025 class Rndm;
0026 class BeamSetup;
0027 class CoupSM;
0028 class CoupSUSY;
0029 class BeamParticle;
0030 class PartonSystems;
0031 class SigmaTotal;
0032 class SigmaCombined;
0033 class HadronWidths;
0034 
0035 // Forward declaration of HIInfo class.
0036 class HIInfo;
0037 
0038 //==========================================================================
0039 
0040 // The Info class contains a mixed bag of information on the event
0041 // generation activity, especially on the current subprocess properties,
0042 // and on the number of errors encountered. This is used by the
0043 // generation machinery.
0044 
0045 class Info {
0046 
0047 public:
0048 
0049   // Constructors.
0050   Info() = default;
0051   Info(bool) : Info() {}
0052 
0053   // Destructor for clean-up.
0054   ~Info() {
0055     if (hasOwnEventAttributes && eventAttributes != nullptr)
0056       delete eventAttributes;
0057   }
0058 
0059   // Assignment operator gives a shallow copy; no objects pointed to
0060   // are copied.
0061   Info & operator=(const Info &) = default;
0062 
0063   // Pointers to other class objects, carried piggyback by Info.
0064 
0065   // Set pointers to other class objects.
0066   void setPtrs(Settings* settingsPtrIn, ParticleData* particleDataPtrIn,
0067     Logger* loggerPtrIn, Rndm* rndmPtrIn, BeamSetup* beamSetupIn,
0068     CoupSM* coupSMPtrIn, CoupSUSY* coupSUSYPtrIn,
0069     PartonSystems* partonSystemsPtrIn, SigmaTotal* sigmaTotPtrIn,
0070     SigmaCombined* sigmaCmbPtrIn, HadronWidths* hadronWidthsPtrIn,
0071     WeightContainer* weightContainerPtrIn) {
0072     settingsPtr = settingsPtrIn; particleDataPtr = particleDataPtrIn;
0073     loggerPtr = loggerPtrIn; rndmPtr = rndmPtrIn; beamSetupPtr = beamSetupIn;
0074     coupSMPtr = coupSMPtrIn; coupSUSYPtr = coupSUSYPtrIn;
0075     partonSystemsPtr = partonSystemsPtrIn; sigmaTotPtr = sigmaTotPtrIn;
0076     sigmaCmbPtr = sigmaCmbPtrIn; hadronWidthsPtr = hadronWidthsPtrIn;
0077     weightContainerPtr = weightContainerPtrIn; }
0078 
0079   // Pointer to the settings database.
0080   Settings*      settingsPtr{};
0081 
0082   // Pointer to the particle data table.
0083   ParticleData*  particleDataPtr{};
0084 
0085   // Pointer to the logger.
0086   Logger*        loggerPtr{};
0087 
0088   // Pointer to the random number generator.
0089   Rndm*          rndmPtr{};
0090 
0091   // Pointers to beam configuration.
0092   BeamSetup*     beamSetupPtr{};
0093 
0094   // Pointers to Standard Model and Beyond SM couplings.
0095   CoupSM*        coupSMPtr{};
0096   CoupSUSY*      coupSUSYPtr{};
0097 
0098   // Pointer to information on subcollision parton locations.
0099   PartonSystems* partonSystemsPtr{};
0100 
0101   // Pointers to the total/elastic/diffractive cross sections.
0102   SigmaTotal*    sigmaTotPtr{};
0103   SigmaCombined* sigmaCmbPtr{};
0104 
0105   // Pointer to the hadron widths data table.
0106   HadronWidths*  hadronWidthsPtr;
0107 
0108   // Pointer to the UserHooks object set for the run.
0109   UserHooksPtr   userHooksPtr{};
0110 
0111   // Mutex that should be locked for thread-unsafe code.
0112   mutex*         mutexPtr{};
0113 
0114   // Pointer to information about a HeavyIons run and the current event.
0115   // (Is nullptr if HeavyIons object is inactive.)
0116   HIInfo*        hiInfo{};
0117 
0118   WeightContainer* weightContainerPtr{};
0119 
0120   // Listing of most available information on current event.
0121   void   list() const;
0122 
0123   // Beam particles (in rest frame). CM energy of event.
0124   int    idA()                const {return idASave;}
0125   int    idB()                const {return idBSave;}
0126   double pzA()                const {return pzASave;}
0127   double pzB()                const {return pzBSave;}
0128   double eA()                 const {return eASave;}
0129   double eB()                 const {return eBSave;}
0130   double mA()                 const {return mASave;}
0131   double mB()                 const {return mBSave;}
0132   double eCM()                const {return eCMSave;}
0133   double s()                  const {return sSave;}
0134 
0135   // Warnings from initialization.
0136   bool   tooLowPTmin()        const {return lowPTmin;}
0137 
0138   // Process name and code, and the number of final-state particles.
0139   string name()               const {return nameSave;}
0140   int    code()               const {return codeSave;}
0141   int    nFinal()             const {return nFinalSave;}
0142 
0143   // Are beam particles resolved, with pdf's? Are they diffractive?
0144   bool   isResolved()         const {return isRes;}
0145   bool   isDiffractiveA()     const {return isDiffA;}
0146   bool   isDiffractiveB()     const {return isDiffB;}
0147   bool   isDiffractiveC()     const {return isDiffC;}
0148   bool   isNonDiffractive()   const {return isND;}
0149   bool   isElastic()          const {return (codeSave == 102);}
0150   // Retained for backwards compatibility.
0151   bool   isMinBias()          const {return isND;}
0152 
0153   // Information for Les Houches Accord and reading files.
0154   bool   isLHA()              const {return isLH;}
0155   bool   atEndOfFile()        const {return atEOF;}
0156 
0157   // For nondiffractive and Les Houches Accord identify hardest subprocess.
0158   bool   hasSub(int i = 0)    const {return hasSubSave[i];}
0159   string nameSub(int i = 0)   const {return nameSubSave[i];}
0160   int    codeSub(int i = 0)   const {return codeSubSave[i];}
0161   int    nFinalSub(int i = 0) const {return nFinalSubSave[i];}
0162 
0163   // Incoming parton flavours and x values.
0164   int    id1(int i = 0)       const {return id1Save[i];}
0165   int    id2(int i = 0)       const {return id2Save[i];}
0166   double x1(int i = 0)        const {return x1Save[i];}
0167   double x2(int i = 0)        const {return x2Save[i];}
0168   double y(int i = 0)         const {return 0.5 * log( x1Save[i] / x2Save[i]);}
0169   double tau(int i = 0)       const {return x1Save[i] * x2Save[i];}
0170 
0171   // Hard process flavours, x values, parton densities, couplings, Q2 scales.
0172   int    id1pdf(int i = 0)    const {return id1pdfSave[i];}
0173   int    id2pdf(int i = 0)    const {return id2pdfSave[i];}
0174   double x1pdf(int i = 0)     const {return x1pdfSave[i];}
0175   double x2pdf(int i = 0)     const {return x2pdfSave[i];}
0176   double pdf1(int i = 0)      const {return pdf1Save[i];}
0177   double pdf2(int i = 0)      const {return pdf2Save[i];}
0178   double QFac(int i = 0)      const {return sqrtpos(Q2FacSave[i]);}
0179   double Q2Fac(int i = 0)     const {return Q2FacSave[i];}
0180   bool   isValence1()         const {return isVal1;}
0181   bool   isValence2()         const {return isVal2;}
0182   double alphaS(int i = 0)    const {return alphaSSave[i];}
0183   double alphaEM(int i = 0)   const {return alphaEMSave[i];}
0184   double QRen(int i = 0)      const {return sqrtpos(Q2RenSave[i]);}
0185   double Q2Ren(int i = 0)     const {return Q2RenSave[i];}
0186   double scalup(int i = 0)    const {return scalupSave[i];}
0187 
0188   // DIS-specific kinematic variables.
0189   double Q2DIS()              const {return Q2DISSave;}
0190   double WDIS()               const {return WDISSave;}
0191   double xDIS()               const {return xDISSave;}
0192   double yDIS()               const {return yDISSave;}
0193   bool   isDIS()              const {return isDISSave;}
0194 
0195   // Kinematics of photons from lepton beams.
0196   double xGammaA()            const {return x1GammaSave;}
0197   double xGammaB()            const {return x2GammaSave;}
0198   double Q2GammaA()           const {return Q2Gamma1Save;}
0199   double Q2GammaB()           const {return Q2Gamma2Save;}
0200   double eCMsub()             const {return eCMsubSave;}
0201   double thetaScatLepA()      const {return thetaLepton1;}
0202   double thetaScatLepB()      const {return thetaLepton2;}
0203   double sHatNew()            const {return sHatNewSave;}
0204   int    photonMode()         const {return gammaModeEvent;}
0205 
0206   // Information on VMD state inside a photon.
0207   bool   isVMDstateA()        const {return isVMDstateAEvent;}
0208   bool   isVMDstateB()        const {return isVMDstateBEvent;}
0209   int    idVMDA()             const {return idVMDASave;}
0210   int    idVMDB()             const {return idVMDBSave;}
0211   double mVMDA()              const {return mVMDASave;}
0212   double mVMDB()              const {return mVMDBSave;}
0213   double scaleVMDA()          const {return scaleVMDASave;}
0214   double scaleVMDB()          const {return scaleVMDBSave;}
0215 
0216   // Mandelstam variables (notation as if subcollision).
0217   double mHat(int i = 0)      const {return sqrt(sH[i]);}
0218   double sHat(int i = 0)      const {return sH[i];}
0219   double tHat(int i = 0)      const {return tH[i];}
0220   double uHat(int i = 0)      const {return uH[i];}
0221   double pTHat(int i = 0)     const {return pTH[i];}
0222   double pT2Hat(int i = 0)    const {return pTH[i] * pTH[i];}
0223   double m3Hat(int i = 0)     const {return m3H[i];}
0224   double m4Hat(int i = 0)     const {return m4H[i];}
0225   double thetaHat(int i = 0)  const {return thetaH[i];}
0226   double phiHat(int i = 0)    const {return phiH[i];}
0227 
0228   // Weight of current event; normally 1, but used for Les Houches events
0229   // or when reweighting phase space selection. Conversion from mb to pb
0230   // for LHA strategy +-4. Uncertainty variations can be accessed by
0231   // providing an index >= 1 (0 = no variations). Also cumulative sum.
0232   double weight(int i = 0)    const;
0233   double weightSum()          const;
0234   double lhaStrategy()        const {return lhaStrategySave;}
0235 
0236   // Further access to uncertainty weights: number and labels.
0237   int    nWeights() const {
0238     return weightContainerPtr->weightsShowerPtr->getWeightsSize() +
0239       weightContainerPtr->weightsFragmentation.getWeightsSize() - 1;}
0240   string weightLabel(int iWgt) const;
0241   int    nWeightGroups() const {return weightContainerPtr->
0242       weightsShowerPtr->nWeightGroups() + weightContainerPtr->
0243       weightsFragmentation.nWeightGroups();}
0244   string getGroupName(int iGN) const;
0245   double getGroupWeight(int iGW) const;
0246 
0247   // Number of times other steps have been carried out.
0248   int    nISR()               const {return nISRSave;}
0249   int    nFSRinProc()         const {return nFSRinProcSave;}
0250   int    nFSRinRes()          const {return nFSRinResSave;}
0251 
0252   // Maximum pT scales for MPI, ISR and FSR (in hard process).
0253   double pTmaxMPI()           const {return pTmaxMPISave;}
0254   double pTmaxISR()           const {return pTmaxISRSave;}
0255   double pTmaxFSR()           const {return pTmaxFSRSave;}
0256 
0257   // Current evolution scale (for UserHooks).
0258   double pTnow()              const {return pTnowSave;}
0259 
0260   // Impact parameter picture, global information
0261   double a0MPI()              const {return a0MPISave;}
0262 
0263   // Impact parameter picture, as set by hardest interaction.
0264   double bMPI()               const {return (bIsSet) ? bMPISave : 1.;}
0265   double enhanceMPI()         const {return (bIsSet) ? enhanceMPISave : 1.;}
0266   double enhanceMPIavg()      const {return (bIsSet) ? enhanceMPIavgSave : 1.;}
0267   double eMPI(int i)          const {return (bIsSet) ? eMPISave[i] : 1.;}
0268   double bMPIold()            const {return (bIsSet) ? bMPIoldSave : 1.;}
0269   double enhanceMPIold()      const {return (bIsSet) ? enhanceMPIoldSave : 1.;}
0270   double enhanceMPIoldavg()   const {return (bIsSet)
0271                                      ? enhanceMPIoldavgSave : 1.;}
0272 
0273   // Number of multiparton interactions, with code and pT for them.
0274   int    nMPI()               const {return nMPISave;}
0275   int    codeMPI(int i)       const {return codeMPISave[i];}
0276   double pTMPI(int i)         const {return pTMPISave[i];}
0277   int    iAMPI(int i)         const {return iAMPISave[i];}
0278   int    iBMPI(int i)         const {return iBMPISave[i];}
0279 
0280   // Cross section estimate, optionally process by process.
0281   vector<int> codesHard();
0282 
0283   // Name of the specified process.
0284   string nameProc(int i = 0)  const {
0285     if (i == 0) return "sum";
0286     auto itr = procNameM.find(i);
0287     if (itr != procNameM.end()) return itr->second;
0288     loggerPtr->ERROR_MSG("process code not found", to_string(i));
0289     return "unknown process";}
0290 
0291   // The number of phase-space points tried.
0292   long   nTried(int i = 0) const {
0293     if (i == 0) return nTry;
0294     auto itr = nTryM.find(i);
0295     if (itr != nTryM.end()) return itr->second;
0296     loggerPtr->ERROR_MSG("process code not found", to_string(i));
0297     return 0;}
0298 
0299   // The number of selected hard processes.
0300   long   nSelected(int i = 0) const {
0301     if (i == 0) return nSel;
0302     auto itr = nSelM.find(i);
0303     if (itr != nSelM.end()) return itr->second;
0304     loggerPtr->ERROR_MSG("process code not found", to_string(i));
0305     return 0;}
0306 
0307   // The number of accepted events.
0308   long   nAccepted(int i = 0) const {
0309     if (i == 0) return nAcc;
0310     auto itr = nAccM.find(i);
0311     if (itr != nAccM.end()) return itr->second;
0312     loggerPtr->ERROR_MSG("process code not found", to_string(i));
0313     return 0;}
0314 
0315   // The estimated cross-section in units of mb.
0316   double sigmaGen(int i = 0)  const {
0317     if (i == 0) return sigGen;
0318     auto itr = sigGenM.find(i);
0319     if (itr != sigGenM.end()) return itr->second;
0320     loggerPtr->ERROR_MSG("process code not found", to_string(i));
0321     return 0;}
0322 
0323   // The uncertainty on the estimated cross-section in units of mb.
0324   double sigmaErr(int i = 0)  const {
0325     if (i == 0) return sigErr;
0326     auto itr = sigErrM.find(i);
0327     if (itr != sigErrM.end()) return itr->second;
0328     //loggerPtr->ERROR_MSG("process code not found", to_string(code));
0329     return 0;}
0330 
0331   // Counters for number of loops in various places.
0332   int    getCounter( int i)   const {return counters[i];}
0333 
0334   // Set or increase the value stored in a counter.
0335   void   setCounter( int i, int value = 0) {counters[i]  = value;}
0336   void   addCounter( int i, int value = 1) {counters[i] += value;}
0337 
0338   // Set initialization warning flag when too low pTmin in ISR/FSR/MPI.
0339   void   setTooLowPTmin(bool lowPTminIn) {lowPTmin = lowPTminIn;}
0340 
0341   // Set info on valence character of hard collision partons.
0342   void setValence( bool isVal1In, bool isVal2In) {isVal1 = isVal1In;
0343     isVal2 = isVal2In;}
0344 
0345   // Set and get some MPI/ISR/FSR properties needed for matching,
0346   // i.e. mainly of internal relevance.
0347   void   hasHistory( bool hasHistoryIn) {hasHistorySave = hasHistoryIn;}
0348   bool   hasHistory() {return hasHistorySave;}
0349   void   zNowISR( double zNowIn) {zNowISRSave = zNowIn;}
0350   double zNowISR() {return zNowISRSave;}
0351   void   pT2NowISR( double pT2NowIn) {pT2NowISRSave = pT2NowIn;}
0352   double pT2NowISR() {return pT2NowISRSave;}
0353 
0354   // Return merging weight.
0355   double mergingWeight(int i=0) const {
0356     return weightContainerPtr->weightsMerging.getWeightsValue(i);}
0357 
0358   // Return the complete NLO weight.
0359   double mergingWeightNLO(int i=0) const {
0360     return weightContainerPtr->weightsMerging.getWeightsValue(i); }
0361 
0362   // Return an LHEF header
0363   string header(const string &key) const {
0364     if (headers.find(key) == headers.end()) return "";
0365     else return headers.at(key);
0366   }
0367 
0368   // Return a list of all header key names
0369   vector<string> headerKeys() const;
0370 
0371   // Return the number of processes in the LHEF.
0372   int nProcessesLHEF() const { return int(sigmaLHEFSave.size());}
0373   // Return the cross section information read from LHEF.
0374   double sigmaLHEF(int iProcess) const { return sigmaLHEFSave[iProcess];}
0375 
0376   // LHEF3 information: Public for easy access.
0377   int LHEFversionSave;
0378 
0379   // Save process information as read from init block of LHEF.
0380   vector<double> sigmaLHEFSave;
0381 
0382   // Contents of the LHAinitrwgt tag
0383   LHAinitrwgt *initrwgt{};
0384 
0385   // Contents of the LHAgenerator tags.
0386   vector<LHAgenerator > *generators{};
0387 
0388   // A map of the LHAweightgroup tags, indexed by name.
0389   map<string,LHAweightgroup > *weightgroups{};
0390 
0391   // A map of the LHAweight tags, indexed by name.
0392   map<string,LHAweight > *init_weights{};
0393 
0394   // Store current-event Les Houches event tags.
0395   bool hasOwnEventAttributes{};
0396   map<string, string > *eventAttributes{};
0397 
0398   // The weights associated with this event, as given by the LHAwgt tags
0399   map<string,double > *weights_detailed{};
0400 
0401   // The weights associated with this event, as given by the LHAweights tags
0402   vector<double > *weights_compressed{};
0403 
0404   // Contents of the LHAscales tag.
0405   LHAscales *scales{};
0406 
0407   // Contents of the LHAweights tag (compressed format)
0408   LHAweights *weights{};
0409 
0410   // Contents of the LHArwgt tag (detailed format)
0411   LHArwgt *rwgt{};
0412 
0413   // Vectorized version of LHArwgt tag, for easy and ordered access.
0414   vector<double> weights_detailed_vector;
0415 
0416   // Value of the unit event weight read from a Les Houches event, necessary
0417   // if additional weights in an otherwise unweighted input file are in
0418   // relation to this number.
0419   double eventWeightLHEF = {};
0420 
0421   // Set the LHEF3 objects read from the init and header blocks.
0422   void setLHEF3InitInfo();
0423   void setLHEF3InitInfo( int LHEFversionIn, LHAinitrwgt *initrwgtIn,
0424     vector<LHAgenerator> *generatorsIn,
0425     map<string,LHAweightgroup> *weightgroupsIn,
0426     map<string,LHAweight> *init_weightsIn, string headerBlockIn );
0427   // Set the LHEF3 objects read from the event block.
0428   void setLHEF3EventInfo();
0429   void setLHEF3EventInfo( map<string, string> *eventAttributesIn,
0430     map<string, double > *weights_detailedIn,
0431     vector<double > *weights_compressedIn,
0432     LHAscales *scalesIn, LHAweights *weightsIn,
0433     LHArwgt *rwgtIn, vector<double> weights_detailed_vecIn,
0434     vector<string> weights_detailed_name_vecIn,
0435     string eventCommentsIn, double eventWeightLHEFIn );
0436 
0437   // Retrieve events tag information.
0438   string getEventAttribute(string key, bool doRemoveWhitespace = false) const;
0439 
0440   // Externally set event tag auxiliary information.
0441   void setEventAttribute(string key, string value, bool doOverwrite = true) {
0442     if (eventAttributes == nullptr) {
0443       eventAttributes = new map<string,string>();
0444       hasOwnEventAttributes = true;
0445     }
0446     map<string, string>::iterator it = eventAttributes->find(key);
0447     if ( it != eventAttributes->end() && !doOverwrite ) return;
0448     if ( it != eventAttributes->end() ) eventAttributes->erase(it);
0449     eventAttributes->insert(make_pair(key,value));
0450   }
0451 
0452   // Retrieve LHEF version
0453   int LHEFversion() const { return LHEFversionSave; }
0454 
0455   // Retrieve initrwgt tag information.
0456   unsigned int getInitrwgtSize() const;
0457 
0458   // Retrieve generator tag information.
0459   unsigned int getGeneratorSize() const;
0460   string getGeneratorValue(unsigned int n = 0) const;
0461   string getGeneratorAttribute( unsigned int n, string key,
0462     bool doRemoveWhitespace = false) const;
0463 
0464   // Retrieve rwgt tag information.
0465   unsigned int getWeightsDetailedSize() const;
0466   double getWeightsDetailedValue(string n) const;
0467   string getWeightsDetailedAttribute(string n, string key,
0468     bool doRemoveWhitespace = false) const;
0469 
0470   // Retrieve weights tag information.
0471   unsigned int getWeightsCompressedSize() const;
0472   double getWeightsCompressedValue(unsigned int n) const;
0473   string getWeightsCompressedAttribute(string key,
0474     bool doRemoveWhitespace = false) const;
0475 
0476   // Retrieve scales tag information.
0477   string getScalesValue(bool doRemoveWhitespace = false) const;
0478   double getScalesAttribute(string key) const;
0479 
0480   // Retrieve complete header block and event comments
0481   // Retrieve scales tag information.
0482   string getHeaderBlock() const { return headerBlock;}
0483   string getEventComments() const { return eventComments;}
0484 
0485   // Set LHEF headers
0486   void setHeader(const string &key, const string &val)
0487     { headers[key] = val; }
0488 
0489   // Set abort in parton level.
0490   void setAbortPartonLevel(bool abortIn) {abortPartonLevel = abortIn;}
0491   bool getAbortPartonLevel() const {return abortPartonLevel;}
0492 
0493   // Get information on hard diffractive events.
0494   bool   hasUnresolvedBeams() const {return hasUnresBeams;}
0495   bool   hasPomPsystem()      const {return hasPomPsys;}
0496   bool   isHardDiffractive()  const {return isHardDiffA || isHardDiffB;}
0497   bool   isHardDiffractiveA() const {return isHardDiffA;}
0498   bool   isHardDiffractiveB() const {return isHardDiffB;}
0499   double xPomeronA()          const {return xPomA;}
0500   double xPomeronB()          const {return xPomB;}
0501   double tPomeronA()          const {return tPomA;}
0502   double tPomeronB()          const {return tPomB;}
0503 
0504   // History information needed to setup the weak shower for 2 -> n.
0505   vector<int> getWeakModes() const {return weakModes;}
0506   vector<pair<int,int> > getWeakDipoles() const {return weakDipoles;}
0507   vector<Vec4> getWeakMomenta() const {return weakMomenta;}
0508   vector<int> getWeak2to2lines() const {return weak2to2lines;}
0509   void setWeakModes(vector<int> weakModesIn) {weakModes = weakModesIn;}
0510   void setWeakDipoles(vector<pair<int,int> > weakDipolesIn)
0511     {weakDipoles = weakDipolesIn;}
0512   void setWeakMomenta(vector<Vec4> weakMomentaIn)
0513     {weakMomenta = weakMomentaIn;}
0514   void setWeak2to2lines(vector<int> weak2to2linesIn)
0515     {weak2to2lines = weak2to2linesIn;}
0516 
0517   // Check if onia is included in showers (used for MPI process setup).
0518   void setOniumShower(bool oniumShowerIn) {oniumShower = oniumShowerIn;}
0519   bool getOniumShower() const {return oniumShower;}
0520 
0521   // Check whether in initialization stage (init) or generating events (next).
0522   void setInInit(bool inInitIn) {inInitSave = inInitIn;}
0523   bool getInInit() const {return inInitSave;}
0524 
0525   // From here on what used to be the private part of the class.
0526 
0527   // Allow conversion from mb to pb.
0528   static const double CONVERTMB2PB;
0529 
0530   // Store common beam quantities.
0531   int    idASave{}, idBSave{};
0532   double pzASave{}, eASave{}, mASave{}, pzBSave{},  eBSave{}, mBSave{},
0533          eCMSave{}, sSave{};
0534 
0535   // Store initialization information.
0536   bool   lowPTmin, inInitSave{};
0537 
0538   // Store common integrated cross section quantities.
0539   long   nTry{}, nSel{}, nAcc{};
0540   double sigGen{}, sigErr{}, wtAccSum{};
0541   map<int, string> procNameM;
0542   map<int, long> nTryM, nSelM, nAccM;
0543   map<int, double> sigGenM, sigErrM;
0544   int    lhaStrategySave{};
0545 
0546   // Store common MPI information.
0547   double a0MPISave;
0548 
0549   // Store current-event quantities.
0550   bool   isRes{}, isDiffA{}, isDiffB{}, isDiffC{}, isND{}, isLH{},
0551          hasSubSave[4], bIsSet{}, evolIsSet{}, atEOF{}, isVal1{}, isVal2{},
0552          hasHistorySave{}, abortPartonLevel{}, isHardDiffA{}, isHardDiffB{},
0553          hasUnresBeams{}, hasPomPsys{};
0554   int    codeSave{}, codeSubSave[4], nFinalSave{}, nFinalSubSave[4], nTotal{},
0555          id1Save[4], id2Save[4], id1pdfSave[4], id2pdfSave[4], nMPISave{},
0556          nISRSave{}, nFSRinProcSave{}, nFSRinResSave{};
0557   double x1Save[4], x2Save[4], x1pdfSave[4], x2pdfSave[4], pdf1Save[4],
0558          pdf2Save[4], Q2FacSave[4], alphaEMSave[4], alphaSSave[4],
0559          Q2RenSave[4], scalupSave[4], sH[4], tH[4], uH[4], pTH[4], m3H[4],
0560          m4H[4], thetaH[4], phiH[4], bMPISave{}, enhanceMPISave{},
0561          enhanceMPIavgSave{}, bMPIoldSave{}, enhanceMPIoldSave{},
0562          enhanceMPIoldavgSave{}, pTmaxMPISave{}, pTmaxISRSave{},
0563          pTmaxFSRSave{}, pTnowSave{}, zNowISRSave{}, pT2NowISRSave{}, xPomA{},
0564          xPomB{}, tPomA{}, tPomB{};
0565   string nameSave{}, nameSubSave[4];
0566   vector<int>    codeMPISave, iAMPISave, iBMPISave;
0567   vector<double> pTMPISave, eMPISave;
0568 
0569   // DIS-specific kinematic variables.
0570   double Q2DISSave{}, WDISSave{}, xDISSave{}, yDISSave{};
0571   bool   isDISSave{};
0572 
0573   // Variables related to photon kinematics.
0574   bool   isVMDstateAEvent{}, isVMDstateBEvent{};
0575   int    gammaModeEvent{}, idVMDASave{}, idVMDBSave{};
0576   double x1GammaSave{}, x2GammaSave{}, Q2Gamma1Save{}, Q2Gamma2Save{},
0577          eCMsubSave{}, thetaLepton1{}, thetaLepton2{}, sHatNewSave{},
0578          mVMDASave{}, mVMDBSave{}, scaleVMDASave{}, scaleVMDBSave{};
0579 
0580   // Vector of various loop counters.
0581   int    counters[50];
0582 
0583   // Map for LHEF headers.
0584   map<string, string> headers;
0585 
0586   // Strings for complete header block and event comments.
0587   string headerBlock{}, eventComments{};
0588 
0589   // Set info on the two incoming beams: only from Pythia class.
0590   void setBeamIDs( int idAin, int idBin) { idASave = idAin; idBSave = idBin;}
0591   void setBeamA( int idAin, double pzAin, double eAin, double mAin) {
0592     idASave = idAin; pzASave = pzAin; eASave = eAin; mASave = mAin;}
0593   void setBeamB( int idBin, double pzBin, double eBin, double mBin) {
0594     idBSave = idBin; pzBSave = pzBin; eBSave = eBin; mBSave = mBin;}
0595   void setECM( double eCMin) {eCMSave = eCMin; sSave = eCMSave * eCMSave;}
0596 
0597   // Set info on DIS-specific kinematic variables.
0598   void setDISKinematics(double Q2In, double WIn, double xIn, double yIn) {
0599     Q2DISSave = Q2In; WDISSave = WIn; xDISSave = xIn; yDISSave = yIn;
0600     isDISSave = true;}
0601   // Set info related to gamma+gamma subcollision.
0602   void setX1Gamma( double x1GammaIn)     { x1GammaSave    = x1GammaIn;   }
0603   void setX2Gamma( double x2GammaIn)     { x2GammaSave    = x2GammaIn;   }
0604   void setQ2Gamma1( double Q2gammaIn)    { Q2Gamma1Save   = Q2gammaIn;   }
0605   void setQ2Gamma2( double Q2gammaIn)    { Q2Gamma2Save   = Q2gammaIn;   }
0606   void setTheta1( double theta1In)       { thetaLepton1   = theta1In;    }
0607   void setTheta2( double theta2In)       { thetaLepton2   = theta2In;    }
0608   void setECMsub( double eCMsubIn)       { eCMsubSave     = eCMsubIn;    }
0609   void setsHatNew( double sHatNewIn)     { sHatNewSave    = sHatNewIn;   }
0610   void setGammaMode( double gammaModeIn) { gammaModeEvent = gammaModeIn; }
0611   // Set info on VMD state.
0612   void setVMDstateA(bool isVMDAIn, int idAIn, double mAIn, double scaleAIn)
0613     {isVMDstateAEvent = isVMDAIn; idVMDASave = idAIn; mVMDASave = mAIn;
0614     scaleVMDASave = scaleAIn;}
0615   void setVMDstateB(bool isVMDBIn, int idBIn, double mBIn, double scaleBIn)
0616     {isVMDstateBEvent = isVMDBIn; idVMDBSave = idBIn; mVMDBSave = mBIn;
0617     scaleVMDBSave = scaleBIn;}
0618 
0619   // Reset info for current event: only from Pythia class.
0620   void clear() {
0621     isRes = isDiffA = isDiffB = isDiffC = isND = isLH = bIsSet
0622       = evolIsSet = atEOF = isVal1 = isVal2 = hasHistorySave
0623       = isHardDiffA = isHardDiffB = hasUnresBeams = hasPomPsys = false;
0624     codeSave = nFinalSave = nTotal = nMPISave = nISRSave = nFSRinProcSave
0625       = nFSRinResSave = 0;
0626     bMPISave = enhanceMPISave = enhanceMPIavgSave = bMPIoldSave
0627       = enhanceMPIoldSave = enhanceMPIoldavgSave = 1.;
0628     pTmaxMPISave = pTmaxISRSave = pTmaxFSRSave = pTnowSave = zNowISRSave
0629       = pT2NowISRSave = 0.;
0630     nameSave = " ";
0631     for (int i = 0; i < 4; ++i) {
0632       hasSubSave[i] = false;
0633       codeSubSave[i] = nFinalSubSave[i] = id1pdfSave[i] = id2pdfSave[i]
0634         = id1Save[i] = id2Save[i] = 0;
0635       x1pdfSave[i] = x2pdfSave[i] = pdf1Save[i] = pdf2Save[i]
0636         = Q2FacSave[i] = alphaEMSave[i] = alphaSSave[i] = Q2RenSave[i]
0637         = scalupSave[i] = x1Save[i] = x2Save[i] = sH[i] = tH[i] = uH[i]
0638         = pTH[i] = m3H[i] = m4H[i] = thetaH[i] = phiH[i] = 0.;
0639       nameSubSave[i] = " ";
0640     }
0641     codeMPISave.resize(0); iAMPISave.resize(0); iBMPISave.resize(0);
0642     pTMPISave.resize(0); eMPISave.resize(0); setHardDiff();
0643     weightContainerPtr->clear();
0644   }
0645 
0646   // Reset info arrays only for parton and hadron level.
0647   int sizeMPIarrays() const { return codeMPISave.size(); }
0648   void resizeMPIarrays(int newSize) { codeMPISave.resize(newSize);
0649     iAMPISave.resize(newSize); iBMPISave.resize(newSize);
0650     pTMPISave.resize(newSize); eMPISave.resize(newSize); }
0651 
0652   // Set info on the (sub)process: from ProcessLevel, ProcessContainer or
0653   // MultipartonInteractions classes.
0654   void setType( string nameIn, int codeIn, int nFinalIn,
0655     bool isNonDiffIn = false, bool isResolvedIn = true,
0656     bool isDiffractiveAin = false, bool isDiffractiveBin = false,
0657     bool isDiffractiveCin = false, bool isLHAin = false) {
0658     nameSave = nameIn; codeSave = codeIn; nFinalSave = nFinalIn;
0659     isND = isNonDiffIn; isRes = isResolvedIn; isDiffA = isDiffractiveAin;
0660     isDiffB = isDiffractiveBin; isDiffC = isDiffractiveCin; isLH = isLHAin;
0661     nTotal = 2 + nFinalSave; bIsSet = false; hasSubSave[0] = false;
0662     nameSubSave[0] = " "; codeSubSave[0] = 0; nFinalSubSave[0] = 0;
0663     evolIsSet = false;}
0664   void setSubType( int iDS, string nameSubIn, int codeSubIn,
0665     int nFinalSubIn) { hasSubSave[iDS] = true; nameSubSave[iDS] = nameSubIn;
0666     codeSubSave[iDS] = codeSubIn; nFinalSubSave[iDS] = nFinalSubIn;}
0667   void setPDFalpha( int iDS, int id1pdfIn, int id2pdfIn, double x1pdfIn,
0668     double x2pdfIn, double pdf1In, double pdf2In, double Q2FacIn,
0669     double alphaEMIn, double alphaSIn, double Q2RenIn, double scalupIn) {
0670     id1pdfSave[iDS] = id1pdfIn; id2pdfSave[iDS] = id2pdfIn;
0671     x1pdfSave[iDS] = x1pdfIn; x2pdfSave[iDS] = x2pdfIn;
0672     pdf1Save[iDS] = pdf1In; pdf2Save[iDS] = pdf2In; Q2FacSave[iDS] = Q2FacIn;
0673     alphaEMSave[iDS] = alphaEMIn; alphaSSave[iDS] = alphaSIn;
0674     Q2RenSave[iDS] = Q2RenIn; scalupSave[iDS] = scalupIn;}
0675   void setScalup( int iDS, double scalupIn) {scalupSave[iDS] = scalupIn;}
0676   void setKin( int iDS, int id1In, int id2In, double x1In, double x2In,
0677     double sHatIn, double tHatIn, double uHatIn, double pTHatIn,
0678     double m3HatIn, double m4HatIn, double thetaHatIn, double phiHatIn) {
0679     id1Save[iDS] = id1In; id2Save[iDS] = id2In; x1Save[iDS] = x1In;
0680     x2Save[iDS] = x2In; sH[iDS] = sHatIn; tH[iDS] = tHatIn;
0681     uH[iDS] = uHatIn; pTH[iDS] = pTHatIn; m3H[iDS] = m3HatIn;
0682     m4H[iDS] = m4HatIn; thetaH[iDS] = thetaHatIn; phiH[iDS] = phiHatIn;}
0683   void setTypeMPI( int codeMPIIn, double pTMPIIn, int iAMPIIn = 0,
0684     int iBMPIIn = 0, double eMPIIn = 1.) {codeMPISave.push_back(codeMPIIn);
0685     pTMPISave.push_back(pTMPIIn); iAMPISave.push_back(iAMPIIn);
0686     iBMPISave.push_back(iBMPIIn); eMPISave.push_back(eMPIIn); }
0687 
0688   // Set info on cross section: from ProcessLevel (reset in Pythia).
0689   void sigmaReset() { nTry = nSel = nAcc = 0; sigGen = sigErr = wtAccSum = 0.;
0690     procNameM.clear(); nTryM.clear(); nSelM.clear(); nAccM.clear();
0691     sigGenM.clear(); sigErrM.clear();}
0692   void setSigma( int i, string procNameIn, long nTryIn, long nSelIn,
0693     long nAccIn, double sigGenIn, double sigErrIn, double wtAccSumIn) {
0694     if (i == 0) {nTry = nTryIn; nSel = nSelIn; nAcc = nAccIn;
0695       sigGen = sigGenIn; sigErr = sigErrIn; wtAccSum = wtAccSumIn;}
0696     else {procNameM[i] = procNameIn; nTryM[i] = nTryIn; nSelM[i] = nSelIn;
0697       nAccM[i] = nAccIn; sigGenM[i] = sigGenIn; sigErrM[i] = sigErrIn;} }
0698   void addSigma( int i, long nTryIn, long nSelIn, long nAccIn, double sigGenIn,
0699     double sigErrIn) { nTryM[i] += nTryIn; nSelM[i] += nSelIn;
0700     nAccM[i] += nAccIn; sigGenM[i] += sigGenIn;
0701     sigErrM[i] = sqrtpos(sigErrM[i]*sigErrM[i] + sigErrIn*sigErrIn); }
0702 
0703   // Set info on impact parameter: from PartonLevel.
0704   void setImpact( double bMPIIn, double enhanceMPIIn, double enhanceMPIavgIn,
0705     bool bIsSetIn = true, bool pushBack = false) {
0706     if (pushBack) {bMPIoldSave = bMPISave; enhanceMPIoldSave = enhanceMPISave;
0707       enhanceMPIoldavgSave = enhanceMPIavgSave;}
0708     bMPISave = bMPIIn; enhanceMPISave = eMPISave[0] = enhanceMPIIn,
0709     enhanceMPIavgSave = enhanceMPIavgIn, bIsSet = bIsSetIn;}
0710 
0711   // Set info on pTmax scales and number of evolution steps: from PartonLevel.
0712   void setPartEvolved( int nMPIIn, int nISRIn) {
0713     nMPISave = nMPIIn; nISRSave = nISRIn;}
0714   void setEvolution( double pTmaxMPIIn, double pTmaxISRIn, double pTmaxFSRIn,
0715     int nMPIIn, int nISRIn, int nFSRinProcIn, int nFSRinResIn) {
0716     pTmaxMPISave = pTmaxMPIIn; pTmaxISRSave = pTmaxISRIn;
0717     pTmaxFSRSave = pTmaxFSRIn; nMPISave = nMPIIn; nISRSave = nISRIn;
0718     nFSRinProcSave = nFSRinProcIn; nFSRinResSave = nFSRinResIn;
0719     evolIsSet = true;}
0720 
0721   // Set current pT evolution scale for MPI/ISR/FSR; from PartonLevel.
0722   void setPTnow( double pTnowIn) {pTnowSave = pTnowIn;}
0723 
0724   // Set a0 from MultipartonInteractions.
0725   void seta0MPI(double a0MPIin) {a0MPISave = a0MPIin;}
0726 
0727   // Set info whether reading of Les Houches Accord file at end.
0728   void setEndOfFile( bool atEOFin) {atEOF = atEOFin;}
0729 
0730   // Set event weight, either for LHEF3 or for uncertainty bands. If weight
0731   // has units (lhaStrategy 4 or -4), input in mb
0732   void setWeight( double weightIn, int lhaStrategyIn) {
0733     for (int i = 0; i < nWeights(); ++i)
0734       weightContainerPtr->weightsShowerPtr->setValueByIndex(i,1.);
0735     // Nominal weight in weightContainer saved in pb for lhaStrategy +-4
0736     weightContainerPtr->setWeightNominal(
0737         abs(lhaStrategyIn) == 4 ? CONVERTMB2PB*weightIn : weightIn);
0738     lhaStrategySave = lhaStrategyIn;}
0739   //
0740   // Set info on resolved processes.
0741   void setIsResolved(bool isResIn) {isRes = isResIn;}
0742 
0743   // Set info on hard diffraction.
0744   void setHardDiff(bool hasUnresBeamsIn = false, bool hasPomPsysIn = false,
0745     bool isHardDiffAIn = false, bool isHardDiffBIn = false,
0746     double xPomAIn = 0., double xPomBIn = 0., double tPomAIn = 0.,
0747     double tPomBIn = 0.) { hasUnresBeams = hasUnresBeamsIn;
0748     hasPomPsys = hasPomPsysIn; isHardDiffA = isHardDiffAIn;
0749     isHardDiffB = isHardDiffBIn;
0750     xPomA = xPomAIn; xPomB = xPomBIn;
0751     tPomA = tPomAIn; tPomB = tPomBIn;}
0752 
0753   // Move process information to an another diffractive system.
0754   void reassignDiffSystem(int iDSold, int iDSnew);
0755 
0756   // Set information in hard diffractive events.
0757   void setHasUnresolvedBeams(bool hasUnresBeamsIn)
0758     {hasUnresBeams = hasUnresBeamsIn;}
0759   void setHasPomPsystem(bool hasPomPsysIn) {hasPomPsys = hasPomPsysIn;}
0760 
0761   // Variables for weak and onia shower setup.
0762   vector<int> weakModes, weak2to2lines;
0763   vector<Vec4> weakMomenta;
0764   vector<pair<int, int> > weakDipoles;
0765   bool oniumShower{false};
0766 
0767   int numberOfWeights() const {
0768     return weightContainerPtr->numberOfWeights(); }
0769   double weightValueByIndex(int key=0) const {
0770     return weightContainerPtr->weightValueByIndex(key); }
0771   string weightNameByIndex(int key=0) const {
0772     return weightContainerPtr->weightNameByIndex(key); }
0773   vector<double> weightValueVector() const {
0774     return weightContainerPtr->weightValueVector(); }
0775   vector<string> weightNameVector() const {
0776     return weightContainerPtr->weightNameVector(); }
0777 
0778   // The random state is saved before the generation of an event
0779   // starts. To be used for debugging purposes.
0780   RndmState currentEventRndmState{};
0781   void dumpRandomState(string fileName="EventRandomState.dat") const;
0782   void readRandomState(string fileName="EventRandomState.dat") const;
0783 
0784   // Special variables for (below-threshold) toponium production.
0785   double toponiumE, toponiumm3, toponiumm4;
0786 
0787 };
0788 
0789 //==========================================================================
0790 
0791 } // end namespace Pythia8
0792 
0793 #endif // Pythia8_Info_H