Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-15 09:06:43

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