Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:06:25

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