Back to home page

EIC code displayed by LXR

 
 

    


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

0001 // SimpleSpaceShower.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 // Header file for the original simple spacelike initial-state showers.
0007 // SpaceDipoleEnd: data on a radiating dipole end in ISR.
0008 // SimpleSpaceShower: handles the showering description.
0009 
0010 #ifndef Pythia8_SimpleSpaceShower_H
0011 #define Pythia8_SimpleSpaceShower_H
0012 
0013 #include "Pythia8/SpaceShower.h"
0014 #include "Pythia8/SimpleWeakShowerMEs.h"
0015 
0016 namespace Pythia8 {
0017 
0018 //==========================================================================
0019 
0020 // Data on radiating dipole ends, only used inside SimpleSpaceShower.
0021 
0022 class SpaceDipoleEnd {
0023 
0024 public:
0025 
0026   // Constructor.
0027   SpaceDipoleEnd( int systemIn = 0, int sideIn = 0, int iRadiatorIn = 0,
0028     int iRecoilerIn = 0, double pTmaxIn = 0., int colTypeIn = 0,
0029     int chgTypeIn = 0, int weakTypeIn = 0,  int MEtypeIn = 0,
0030     bool normalRecoilIn = true, int weakPolIn = 0,
0031     int iColPartnerIn = 0, int idColPartnerIn = 0) :
0032     system(systemIn), side(sideIn), iRadiator(iRadiatorIn),
0033     iRecoiler(iRecoilerIn), pTmax(pTmaxIn), colType(colTypeIn),
0034     chgType(chgTypeIn), weakType(weakTypeIn), MEtype(MEtypeIn),
0035     normalRecoil(normalRecoilIn), weakPol(weakPolIn),
0036     iColPartner(iColPartnerIn), idColPartner(idColPartnerIn),
0037     nBranch(0), idDaughter(), idMother(), idSister(), iFinPol(), x1(), x2(),
0038     m2Dip(), pT2(), z(), xMo(), Q2(), mSister(), m2Sister(), pT2corr(),
0039     pT2Old(0.), zOld(0.5), asymPol(), m2IF(), mColPartner(),
0040     pAccept() { }
0041 
0042   // Store values for trial emission.
0043   void store( int idDaughterIn, int idMotherIn, int idSisterIn,
0044     double x1In, double x2In, double m2DipIn, double pT2In, double zIn,
0045     double xMoIn, double Q2In, double mSisterIn, double m2SisterIn,
0046     double pT2corrIn, int iColPartnerIn, double m2IFIn, double mColPartnerIn)
0047     {idDaughter = idDaughterIn; idMother = idMotherIn;
0048     idSister = idSisterIn; x1 = x1In; x2 = x2In; m2Dip = m2DipIn;
0049     pT2 = pT2In; z = zIn; xMo = xMoIn; Q2 = Q2In; mSister = mSisterIn;
0050     m2Sister = m2SisterIn; pT2corr = pT2corrIn; iColPartner = iColPartnerIn;
0051     m2IF = m2IFIn; mColPartner = mColPartnerIn;}
0052 
0053   // Basic properties related to evolution and matrix element corrections.
0054   int    system, side, iRadiator, iRecoiler;
0055   double pTmax;
0056   int    colType, chgType, weakType, MEtype;
0057   bool   normalRecoil;
0058   int    weakPol, iColPartner, idColPartner;
0059 
0060   // Properties specific to current trial emission.
0061   int    nBranch, idDaughter, idMother, idSister, iFinPol;
0062   double x1, x2, m2Dip, pT2, z, xMo, Q2, mSister, m2Sister, pT2corr,
0063          pT2Old, zOld, asymPol, m2IF, mColPartner;
0064 
0065   // Properties needed for the evaluation of parameter variations
0066   double pAccept;
0067 
0068 } ;
0069 
0070 //==========================================================================
0071 
0072 // The SimpleSpaceShower class does spacelike showers.
0073 
0074 class SimpleSpaceShower : public SpaceShower {
0075 
0076 public:
0077 
0078   // Constructor.
0079   SimpleSpaceShower() : rescatterFail(), gamma2qqbar(), hasWeaklyRadiated(),
0080     iSysSel(), pTmaxFudge(), doQCDshower(), doQEDshowerByQ(), doQEDshowerByL(),
0081     useSamePTasMPI(), doWeakShower(), doMEcorrections(), doMEafterFirst(),
0082     doPhiPolAsym(), doPhiPolAsymHard(), doPhiIntAsym(), doRapidityOrder(),
0083     useFixedFacScale(), doSecondHard(), canVetoEmission(), hasUserHooks(),
0084     alphaSuseCMW(), singleWeakEmission(), vetoWeakJets(), weakExternal(),
0085     doRapidityOrderMPI(), doMPI(), doDipoleRecoil(), doPartonVertex(),
0086     pTmaxMatch(), pTdampMatch(), alphaSorder(), alphaSnfmax(), alphaEMorder(),
0087     nQuarkIn(), enhanceScreening(), weakMode(), pT0paramMode(), pTdampFudge(),
0088     mc(), mb(), m2c(), m2b(), renormMultFac(), factorMultFac(),
0089     fixedFacScale2(), alphaSvalue(), alphaS2pi(), Lambda3flav(), Lambda4flav(),
0090     Lambda5flav(), Lambda3flav2(), Lambda4flav2(), Lambda5flav2(), pT0Ref(),
0091     ecmRef(), ecmPow(), pTmin(), sCM(), eCM(), pT0(), pTminChgQ(), pTminChgL(),
0092     pT20(), pT2min(), pT2minChgQ(), pT2minChgL(), pTweakCut(), pT2weakCut(),
0093     pTmaxFudgeMPI(), strengthIntAsym(), weakEnhancement(), mZ(), gammaZ(),
0094     thetaWRat(), mW(), gammaW(), weakMaxWt(), vetoWeakDeltaR2(), sideA(),
0095     twoHard(), dopTlimit1(), dopTlimit2(), dopTdamp(),  tChannel(),
0096     doUncertaintiesNow(), iNow(), iRec(), idDaughter(), nRad(), idResFirst(),
0097     idResSecond(), xDaughter(), x1Now(), x2Now(),
0098     m2ColPair(), mColPartner(), m2ColPartner(), m2Dip(), m2Rec(), pT2damp(),
0099     pTbegRef(), pdfScale2(), doTrialNow(), canEnhanceEmission(),
0100     canEnhanceTrial(), canEnhanceET(), iDipNow(), iSysNow(), dipEndNow(),
0101     iDipSel(), dipEndSel() { beamOffset = 0; pdfMode = 0; }
0102 
0103   // Destructor.
0104   virtual ~SimpleSpaceShower() override {}
0105 
0106   // Initialize generation. Possibility to force re-initialization by hand.
0107   virtual void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn)
0108     override;
0109 
0110   // Find whether to limit maximum scale of emissions, and whether to dampen.
0111   virtual bool limitPTmax( Event& event, double Q2Fac = 0.,
0112     double Q2Ren = 0.) override;
0113 
0114   // Prepare system for evolution; identify ME.
0115   virtual void prepare( int iSys, Event& event, bool limitPTmaxIn = true)
0116     override;
0117 
0118   // Update dipole list after each FSR emission.
0119   virtual void update( int iSys, Event& event, bool hasWeakRad = false)
0120     override;
0121 
0122   // Select next pT in downwards evolution.
0123   virtual double pTnext( Event& event, double pTbegAll, double pTendAll,
0124     int nRadIn = -1, bool doTrialIn = false) override;
0125 
0126   // ME corrections and kinematics that may give failure.
0127   virtual bool branch( Event& event) override;
0128 
0129   // Print dipole list; for debug mainly.
0130   virtual void list() const override;
0131 
0132   // Initialize data members for calculation of uncertainty bands.
0133   virtual bool initUncertainties() override;
0134 
0135    // Initialize data members for application of enhancements.
0136   virtual bool initEnhancements() override;
0137 
0138   // Flag for failure in branch(...) that will force a retry of parton level.
0139   virtual bool doRestart() const override {return rescatterFail;}
0140 
0141   // Tell if latest scattering was a gamma->qqbar.
0142   virtual bool wasGamma2qqbar() override { return gamma2qqbar; }
0143 
0144   // Tell whether ISR has done a weak emission.
0145   virtual bool getHasWeaklyRadiated() override {return hasWeaklyRadiated;}
0146 
0147   // Tell which system was the last processed one.
0148   virtual int system() const override {return iSysSel;}
0149 
0150   // Potential enhancement factor of pTmax scale for hardest emission.
0151   virtual double enhancePTmax() const override {return pTmaxFudge;}
0152 
0153   // Functions to directly extract the probability of no emission between two
0154   // scales. These functions are not used in the Pythia core code, but can be
0155   // used by external programs to interface with the shower directly.
0156   double noEmissionProbability( double pTbegAll, double pTendAll, double m2dip,
0157     int id, int type, double s = -1., double x = -1.) override;
0158   double pTnext( vector<SpaceDipoleEnd> dipEnds, Event event, double pTbegAll,
0159     double pTendAll, double m2dip, int id, int type, double s = -1.,
0160     double x = -1.);
0161   int pdfMode;
0162 
0163 private:
0164 
0165   // Constants: could only be changed in the code itself.
0166   static const int    MAXLOOPTINYPDF;
0167   static const double MCMIN, MBMIN, CTHRESHOLD, BTHRESHOLD, EVALPDFSTEP,
0168          TINYPDF, TINYKERNELPDF, TINYPT2, HEAVYPT2EVOL, HEAVYXEVOL,
0169          EXTRASPACEQ, LAMBDA3MARGIN, PT2MINWARN, LEPTONXMIN, LEPTONXMAX,
0170          LEPTONPT2MIN, LEPTONFUDGE, WEAKPSWEIGHT, HEADROOMQ2Q, HEADROOMQ2G,
0171          HEADROOMG2G, HEADROOMG2Q, HEADROOMHQG, REJECTFACTOR, PROBLIMIT;
0172 
0173   // Store properties to be returned by methods.
0174   bool   rescatterFail, gamma2qqbar, hasWeaklyRadiated;
0175   int    iSysSel;
0176   double pTmaxFudge;
0177 
0178   // Initialization data, normally only set once.
0179   bool   doQCDshower, doQEDshowerByQ, doQEDshowerByL, useSamePTasMPI,
0180          doWeakShower, doMEcorrections, doMEafterFirst, doPhiPolAsym,
0181          doPhiPolAsymHard, doPhiIntAsym, doRapidityOrder, useFixedFacScale,
0182          doSecondHard, canVetoEmission, hasUserHooks, alphaSuseCMW,
0183          singleWeakEmission, vetoWeakJets, weakExternal, doRapidityOrderMPI,
0184          doMPI, doDipoleRecoil, doPartonVertex;
0185   int    pdfModeSave;
0186   int    pTmaxMatch, pTdampMatch, alphaSorder, alphaSnfmax, alphaEMorder,
0187          nQuarkIn, enhanceScreening, weakMode, pT0paramMode;
0188   double pTdampFudge, mc, mb, m2c, m2b, renormMultFac, factorMultFac,
0189          fixedFacScale2, alphaSvalue, alphaS2pi, Lambda3flav, Lambda4flav,
0190          Lambda5flav, Lambda3flav2, Lambda4flav2, Lambda5flav2, pT0Ref,
0191          ecmRef, ecmPow, pTmin, sCM, eCM, pT0, pTminChgQ, pTminChgL, pT20,
0192          pT2min, pT2minChgQ, pT2minChgL, pTweakCut, pT2weakCut, pTmaxFudgeMPI,
0193          strengthIntAsym, weakEnhancement, mZ, gammaZ, thetaWRat, mW, gammaW,
0194          weakMaxWt, vetoWeakDeltaR2;
0195 
0196   // alphaStrong and alphaEM calculations.
0197   AlphaStrong alphaS;
0198   AlphaEM alphaEM;
0199 
0200   // Weak matrix elements used for corrections both of ISR and FSR.
0201   SimpleWeakShowerMEs  simpleWeakShowerMEs;
0202 
0203   // Some current values.
0204   bool   sideA, twoHard, dopTlimit1, dopTlimit2, dopTdamp, tChannel,
0205          doUncertaintiesNow;
0206   int    iNow, iRec, idDaughter, nRad, idResFirst, idResSecond;
0207   double xDaughter, x1Now, x2Now, m2ColPair, mColPartner, m2ColPartner,
0208          m2Dip, m2Rec, pT2damp, pTbegRef, pdfScale2;
0209 
0210   // Bookkeeping of enhanced  actual or trial emissions (see EPJC (2013) 73).
0211   bool doTrialNow, canEnhanceEmission, canEnhanceTrial, canEnhanceET;
0212   string splittingNameNow, splittingNameSel;
0213   map< double, pair<string,double> > enhanceFactors;
0214   void storeEnhanceFactor(double pT2, string name, double enhanceFactorIn)
0215     { enhanceFactors.insert(make_pair(pT2,make_pair(name,enhanceFactorIn)));}
0216 
0217   // List of emissions in different sides in different systems:
0218   vector<int> nRadA,nRadB;
0219 
0220   // All dipole ends
0221   vector<SpaceDipoleEnd> dipEnd;
0222 
0223   // List of 2 -> 2 momenta for external weak setup.
0224   vector<Vec4> weakMomenta;
0225 
0226   // Pointers to the current and hardest (so far) dipole ends.
0227   int iDipNow, iSysNow;
0228   SpaceDipoleEnd* dipEndNow;
0229   int iDipSel;
0230   SpaceDipoleEnd* dipEndSel;
0231 
0232   // Evolve a QCD dipole end.
0233   void pT2nextQCD( double pT2begDip, double pT2endDip);
0234 
0235   // Evolve a QCD and QED dipole end near heavy quark threshold region.
0236   void pT2nearThreshold( BeamParticle& beam, double m2Massive,
0237     double m2Threshold, double xMaxAbs, double zMinAbs,
0238     double zMaxMassive, int iColPartner);
0239 
0240   // Evolve a QED dipole end.
0241   void pT2nextQED( double pT2begDip, double pT2endDip);
0242 
0243   // Evolve a Weak dipole end.
0244   void pT2nextWeak( double pT2begDip, double pT2endDip);
0245 
0246   // Find class of ME correction.
0247   int findMEtype( int iSys, Event& event, bool weakRadiation = false);
0248 
0249   // Provide maximum of expected ME weight; for preweighting of evolution.
0250   double calcMEmax( int MEtype, int idMother, int idDaughterIn);
0251 
0252   // Provide actual ME weight for current branching.
0253   double calcMEcorr(int MEtype, int idMother, int idDaughterIn, double M2,
0254     double z, double Q2,double m2Sister);
0255 
0256   // Provide actual ME weight for t-channel weak emissions.
0257   double calcMEcorrWeak(int MEtype, double m2, double z,
0258     double pT2, Vec4 pMother, Vec4 pB, Vec4 pDaughter,
0259     Vec4 pB0, Vec4 p1, Vec4 p2, Vec4 pSister);
0260 
0261   // Find coefficient of azimuthal asymmetry from gluon polarization.
0262   void findAsymPol( Event& event, SpaceDipoleEnd* dip);
0263 
0264   // Find a possible colour partner in the case of dipole recoil.
0265   int findColPartner(Event& event, int iSideA, int iSideB, int iSystem);
0266 
0267   // Calculate uncertainty-band weights for accepted/rejected trial branching.
0268   void calcUncertainties(bool accept, double pAcceptIn, double pT20in,
0269     double enhance, SpaceDipoleEnd* dip, Particle* motherPtr,
0270     Particle* sisterPtr);
0271 
0272 };
0273 
0274 //==========================================================================
0275 
0276 } // end namespace Pythia8
0277 
0278 #endif // Pythia8_SimpleSpaceShower_H