Back to home page

EIC code displayed by LXR

 
 

    


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

0001 // PartonDistributions.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 files for parton densities.
0007 // PDF:           base class.
0008 // LHAGrid1:      internal read and use files in the LHAPDF6 lhagrid1 format.
0009 // GRV94L:        GRV 94L parton densities.
0010 // CTEQ5L:        CTEQ 5L parton densities.
0011 // MSTWpdf:       MRST LO*, LO**, MSTW 2008 LO, NLO.
0012 // CTEQ6pdf:      CTEQ 6L, 6L1, 66, CT09 MC1, MC2, MCS.
0013 // ProtonPoint:   unresolved proton with equivalent photon spectrum.
0014 // GRVpiL:        GRV LO pion parton densities.
0015 // PomFix:        Q2-independent Pomeron parton densities.
0016 // PomH1FitAB:    H1 2006 Fit A and Fit B Pomeron PDFs.
0017 // PomH1Jets:     H1 2007 Jets Pomeron PDFs.
0018 // PomHISASD:     a proton masked as a Pomeron for heavy ions applications.
0019 // Lepton:        parton densities inside a lepton.
0020 // LeptonPoint:   an unresolved lepton (mainly dummy).
0021 // NeutrinoPoint: an unresolved neutrino (mainly dummy).
0022 // CJKL:          CJKL parton densities for photons.
0023 // Lepton2gamma:  convolution of photon flux from leptons and photon PDFs.
0024 // PhotonPoint:   an unresolved photon.
0025 // Proton2gammaDZ: Photon flux from protons according to Drees-Zeppenfeld.
0026 // Nucleus2gamma: Photon flux from heavy nuclei.
0027 // EPAexternal:   approximated photon flux used for sampling of external flux.
0028 // nPDF:          a nuclear PDF, derived from a proton ditto.
0029 // Isospin:       isospin modification for nuclear pDF
0030 // EPS09, EPPS16: nuclear modification factors.
0031 
0032 
0033 #ifndef Pythia8_PartonDistributions_H
0034 #define Pythia8_PartonDistributions_H
0035 
0036 #include "Pythia8/Basics.h"
0037 #include "Pythia8/Info.h"
0038 #include "Pythia8/MathTools.h"
0039 #include "Pythia8/ParticleData.h"
0040 #include "Pythia8/PythiaStdlib.h"
0041 #include "Pythia8/SharedPointers.h"
0042 
0043 namespace Pythia8 {
0044 
0045 //==========================================================================
0046 
0047 // Base class for parton distribution functions.
0048 
0049 class PDF {
0050 
0051 public:
0052 
0053   // Constructor.
0054   PDF(int idBeamIn = 2212) : idBeam(idBeamIn), idBeamAbs(abs(idBeam)),
0055   idSav(9), xSav(-1), Q2Sav(-1.), isSet(true), isInit(false),
0056   hasGammaInLepton(false) { resetValenceContent(); }
0057 
0058   // Virtual destructor.
0059   virtual ~PDF() {}
0060 
0061   // Perform additional initialization (LHPADF only).
0062   // Arguments are:
0063   // int idBeamIn, string setName, int member, Logger* loggerPtr
0064   virtual bool init(int, string, int, Logger*) {return true;}
0065 
0066   // Confirm that PDF has been set up (important for LHAPDF and H1 Pomeron).
0067   bool isSetup() {return isSet;}
0068 
0069   // Switch to new beam particle identities; for similar hadrons only.
0070   virtual void setBeamID(int idBeamIn) { idBeam = idBeamIn;
0071     idBeamAbs = abs(idBeam); idSav = 9; xSav = -1.; Q2Sav = -1.;
0072     resetValenceContent();}
0073 
0074   // Set valence content.
0075   void resetValenceContent();
0076   void setValenceContent(int idVal1In, int idVal2In, int idVal3In) {
0077     idVal1 = idVal1In; idVal2 = idVal2In; idVal3 = idVal3In;}
0078 
0079   // Allow extrapolation beyond boundaries. This is optional.
0080   virtual void setExtrapolate(bool) {}
0081 
0082   // Read out parton density.
0083   double xf(int id, double x, double Q2);
0084 
0085   // Read out valence and sea part of parton densities.
0086   double xfVal(int id, double x, double Q2);
0087   double xfSea(int id, double x, double Q2);
0088 
0089   // Check whether x and Q2 values fall inside the fit bounds (LHAPDF6 only).
0090   virtual bool insideBounds(double, double) {return true;}
0091 
0092   // Access the running alpha_s of a PDF set (LHAPDF6 only).
0093   virtual double alphaS(double) { return 1.;}
0094 
0095   // Return quark masses used in the PDF fit (LHAPDF6 only).
0096   virtual double mQuarkPDF(int) { return -1.;}
0097 
0098   // Return number of members of this PDF family (LHAPDF6 only).
0099   virtual int nMembers() { return 1;}
0100 
0101   // Error envelope from PDF uncertainty.
0102   struct PDFEnvelope {
0103     double centralPDF, errplusPDF, errminusPDF, errsymmPDF, scalePDF;
0104     vector<double> pdfMemberVars;
0105     PDFEnvelope() : centralPDF(-1.0), errplusPDF(0.0), errminusPDF(0.0),
0106       errsymmPDF(0.0), scalePDF(-1.0), pdfMemberVars(0.0) {};
0107   };
0108 
0109   // Calculate PDF envelope.
0110   virtual void calcPDFEnvelope(int, double, double, int) {}
0111   virtual void calcPDFEnvelope(pair<int,int>, pair<double,double>, double,
0112     int) {}
0113   virtual PDFEnvelope getPDFEnvelope() { return PDFEnvelope(); }
0114 
0115   // Approximate photon PDFs by decoupling the scale and x-dependence.
0116   virtual double gammaPDFxDependence(int, double) { return 0.; }
0117 
0118   // Provide the reference scale for logarithmic Q^2 evolution for photons.
0119   virtual double gammaPDFRefScale(int) { return 0.; }
0120 
0121   // Sample the valence content for photons.
0122   virtual int sampleGammaValFlavor(double) { return 0.; }
0123 
0124   // The total x-integrated PDFs. Relevant for MPIs with photon beams.
0125   virtual double xfIntegratedTotal(double) { return 0.; }
0126 
0127   // Return the sampled value for x_gamma.
0128   virtual double xGamma() { return 1.; }
0129 
0130   // Keep track of pomeron momentum fraction.
0131   virtual void xPom(double = -1.0) {}
0132 
0133   // Return accurate and approximated photon fluxes and PDFs.
0134   virtual double xfFlux(int , double , double )   { return 0.; }
0135   virtual double xfApprox(int , double , double ) { return 0.; }
0136   virtual double xfGamma(int , double , double )  { return 0.; }
0137   virtual double intFluxApprox()                  { return 0.; }
0138   virtual bool hasApproxGammaFlux()               { return false; }
0139 
0140   // Return the kinematical limits and sample Q2 and x.
0141   virtual double getXmin()                 { return 0.; }
0142   virtual double getXhadr()                { return 0.; }
0143   virtual double sampleXgamma(double )     { return 0.; }
0144   virtual double sampleQ2gamma(double )    { return 0.; }
0145   virtual double fluxQ2dependence(double ) { return 0.; }
0146 
0147   // Normal PDFs unless gamma inside lepton -> an overestimate for sampling.
0148   virtual double xfMax(int id, double x, double Q2) { return xf( id, x, Q2); }
0149 
0150   // Normal PDFs unless gamma inside lepton -> Do not sample x_gamma.
0151   virtual double xfSame(int id, double x, double Q2) { return xf( id, x, Q2); }
0152 
0153   // Allow for new scaling factor for VMD PDFs.
0154   virtual void setVMDscale(double = 1.) {}
0155 
0156   // Return if s/sbar, c/cbar, and b/bbar PDFs are symmetric.
0157   bool sSymmetric() const { return sSymmetricSave; }
0158   bool cSymmetric() const { return cSymmetricSave; }
0159   bool bSymmetric() const { return bSymmetricSave; }
0160 
0161   // Set s/sbar, c/cbar, and b/bbar PDFs symmetric.
0162   void sSymmetric(bool sSymmetricIn) { sSymmetricSave = sSymmetricIn; }
0163   void cSymmetric(bool cSymmetricIn) { cSymmetricSave = cSymmetricIn; }
0164   void bSymmetric(bool bSymmetricIn) { bSymmetricSave = bSymmetricIn; }
0165 
0166 protected:
0167 
0168   // Store relevant quantities.
0169   int idBeam, idBeamAbs, idSav, idVal1, idVal2, idVal3;
0170   double xSav, Q2Sav;
0171   // Stored quantities.
0172   double xu, xd, xs, xubar, xdbar, xsbar, xc, xb, xcbar, xbbar,
0173          xg, xlepton, xgamma;
0174   bool   isSet, isInit;
0175 
0176   // For hadrons, beamType defines special cases and determines how
0177   // to handle isospin symmetries.
0178   //  1: no rearrangement (e.g. p, Sigma+, Omega-, pi+)
0179   // -1: switch u <-> d (e.g. n, Sigma-, Xi-, K0)
0180   //  0: take average of u and d (e.g. Sigma0, Lambda0)
0181   //  2/-2: Delta++/Delta-
0182   //  111: pi0-like special case (pi0, rho0, omega, etc.)
0183   //  221: Other diagonal meson cases (eta, eta', phi, J/psi, Upsilon, etc.)
0184   //  130: K_S,L special cases
0185   int beamType;
0186 
0187   // True if a photon beam inside a lepton beam, otherwise set false.
0188   bool hasGammaInLepton;
0189 
0190   // Whether to treat flavoured PDFs as symmetric, for efficiency.
0191   bool sSymmetricSave = false;
0192   bool cSymmetricSave = true, bSymmetricSave = true;
0193 
0194   // Update parton densities.
0195   virtual void xfUpdate(int id, double x, double Q2) = 0;
0196 
0197   // Small routine for error printout, depending on loggerPtr existing or not.
0198   void printErr(string loc, string errMsg, Logger* loggerPtr = nullptr) {
0199     if (loggerPtr) loggerPtr->errorMsg(loc, errMsg);
0200     else cout << "Error in " + loc + ": "  + errMsg << endl;
0201   }
0202 
0203   // Get the raw stored value for the quark variable corresponding to the id.
0204   double xfRaw(int id) const;
0205 
0206   // Check whether the specified id is a valence quark.
0207   bool isValence(int id) const {
0208     return id != 0 && (id == idVal1 || id == idVal2 || id == idVal3); }
0209 
0210 };
0211 
0212 //==========================================================================
0213 
0214 // The LHAGrid1 can be used to read files in the LHAPDF6 lhagrid1 format,
0215 // assuming that the same x grid is used for all Q subgrids.
0216 // Results are not identical with LHAPDF6, owing to different interpolation.
0217 
0218 class LHAGrid1 : public PDF {
0219 
0220 public:
0221 
0222   // Constructor.
0223   LHAGrid1(int idBeamIn = 2212, string pdfWord = "void",
0224     string xmlPath = "../share/Pythia8/xmldoc/", Logger* loggerPtr = 0)
0225     : PDF(idBeamIn), doExtraPol(false), nx(), nq(), nqSub(), xMin(), xMax(),
0226     qMin(), qMax(), pdfVal(), pdfGrid(), pdfSlope(nullptr) {
0227     init( pdfWord, xmlPath, loggerPtr); };
0228 
0229   // Constructor with a stream.
0230   LHAGrid1(int idBeamIn, istream& is, Logger* loggerPtr = 0)
0231     : PDF(idBeamIn), doExtraPol(false), nx(), nq(), nqSub(), xMin(), xMax(),
0232     qMin(), qMax(), pdfVal(), pdfGrid(), pdfSlope(nullptr) {
0233     init( is, loggerPtr); };
0234 
0235   // Destructor.
0236   ~LHAGrid1() { for (int iid = 0; iid < 12; ++iid) {
0237     for (int iq = 0; iq < nq; ++iq) delete[] pdfGrid[iid][iq];
0238     delete[] pdfGrid[iid]; }
0239     if (pdfSlope) { for (int iid = 0; iid < 12; ++iid) delete[] pdfSlope[iid];
0240     delete[] pdfSlope;} };
0241 
0242   // Allow extrapolation beyond boundaries. This is optional.
0243   void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
0244 
0245 private:
0246 
0247   // Variables to be set during code initialization.
0248   bool   doExtraPol;
0249   int    nx, nq, nqSub;
0250   vector<int> nqSum;
0251   double xMin, xMax, qMin, qMax, pdfVal[12];
0252   vector<double> xGrid, lnxGrid, qGrid, lnqGrid, qDiv;
0253   double** pdfGrid[12];
0254   double** pdfSlope;
0255 
0256   // These inits do not overwrite PDF init (prevents Clang warnings).
0257   using PDF::init;
0258 
0259   // Initialization of data array.
0260   void init( string pdfSet, string pdfdataPath, Logger* loggerPtr);
0261 
0262   // Initialization through a stream.
0263   void init( istream& is, Logger* loggerPtr);
0264 
0265   // Update PDF values.
0266   void xfUpdate(int id, double x, double Q2) override;
0267 
0268   // Interpolation in the grid for a given PDF flavour.
0269   void xfxevolve(double x, double Q2);
0270 
0271 };
0272 
0273 //==========================================================================
0274 
0275 // Gives the GRV 94L (leading order) parton distribution function set
0276 // in parametrized form. Authors: M. Glueck, E. Reya and A. Vogt.
0277 
0278 class GRV94L : public PDF {
0279 
0280 public:
0281 
0282   // Constructor.
0283   GRV94L(int idBeamIn = 2212) : PDF(idBeamIn) {}
0284 
0285 private:
0286 
0287   // Update PDF values.
0288   void xfUpdate(int , double x, double Q2) override;
0289 
0290   // Auxiliary routines used during the updating.
0291   double grvv (double x, double n, double ak, double bk, double a,
0292     double b, double c, double d);
0293   double grvw (double x, double s, double al, double be, double ak,
0294     double bk, double a, double b, double c, double d, double e, double es);
0295   double grvs (double x, double s, double sth, double al, double be,
0296     double ak, double ag, double b, double d, double e, double es);
0297 
0298 };
0299 
0300 //==========================================================================
0301 
0302 // Gives the CTEQ 5L (leading order) parton distribution function set
0303 // in parametrized form. Parametrization by J. Pumplin. Authors: CTEQ.
0304 
0305 class CTEQ5L : public PDF {
0306 
0307 public:
0308 
0309   // Constructor.
0310   CTEQ5L(int idBeamIn = 2212) : PDF(idBeamIn) {}
0311 
0312 private:
0313 
0314   // Update PDF values.
0315   void xfUpdate(int , double x, double Q2) override;
0316 
0317 };
0318 
0319 //==========================================================================
0320 
0321 // The MSTWpdf class.
0322 // MRST LO*(*) and MSTW 2008 PDF's, specifically the LO one.
0323 // Original C++ version by Jeppe Andersen.
0324 // Modified by Graeme Watt <watt(at)hep.ucl.ac.uk>.
0325 // Sets available:
0326 // iFit = 1 : MRST LO*  (2007).
0327 // iFit = 2 : MRST LO** (2008).
0328 // iFit = 3 : MSTW 2008 LO, central member.
0329 // iFit = 4 : MSTW 2008 NLO, central member. (Warning!)
0330 
0331 class MSTWpdf : public PDF {
0332 
0333 public:
0334 
0335   // Constructor.
0336   MSTWpdf(int idBeamIn = 2212, int iFitIn = 1,
0337     string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
0338     : PDF(idBeamIn), iFit(), alphaSorder(), alphaSnfmax(), mCharm(), mBottom(),
0339     alphaSQ0(), alphaSMZ(), distance(), tolerance(), xx(), qq(),
0340     c() {init( iFitIn,  pdfdataPath, loggerPtr);}
0341 
0342   // Constructor with a stream.
0343   MSTWpdf(int idBeamIn, istream& is, Logger* loggerPtr = 0)
0344     : PDF(idBeamIn), iFit(), alphaSorder(), alphaSnfmax(), mCharm(), mBottom(),
0345     alphaSQ0(), alphaSMZ(), distance(), tolerance(), xx(), qq(),
0346     c() {init( is, loggerPtr);}
0347 
0348 private:
0349 
0350   // Constants: could only be changed in the code itself.
0351   static const int    np, nx, nq, nqc0, nqb0;
0352   static const double xmin, xmax, qsqmin, qsqmax, xxInit[65], qqInit[49];
0353 
0354   // Data read in from grid file or set at initialization.
0355   int    iFit, alphaSorder, alphaSnfmax;
0356   double mCharm, mBottom, alphaSQ0, alphaSMZ, distance, tolerance,
0357          xx[65], qq[49], c[13][64][48][5][5];
0358 
0359   // These inits do not overwrite PDF init (prevents Clang warnings).
0360   using PDF::init;
0361 
0362   // Initialization of data array.
0363   void init( int iFitIn, string pdfdataPath, Logger* loggerPtr);
0364 
0365   // Initialization through a stream.
0366   void init( istream& is, Logger* loggerPtr);
0367 
0368   // Update PDF values.
0369   void xfUpdate(int , double x, double Q2) override;
0370 
0371   // Evaluate PDF of one flavour species.
0372   double parton(int flavour,double x,double q);
0373   double parton_interpolate(int flavour,double xxx,double qqq);
0374   double parton_extrapolate(int flavour,double xxx,double qqq);
0375 
0376   // Auxiliary routines for evaluation.
0377   int locate(double xx[],int n,double x);
0378   double polderivative1(double x1, double x2, double x3, double y1,
0379     double y2, double y3);
0380   double polderivative2(double x1, double x2, double x3, double y1,
0381     double y2, double y3);
0382   double polderivative3(double x1, double x2, double x3, double y1,
0383     double y2, double y3);
0384 
0385 };
0386 
0387 //==========================================================================
0388 
0389 // The CTEQ6pdf class.
0390 // Proton sets available:
0391 // iFit = 1 : CTEQ6L
0392 // iFit = 2 : CTEQ6L1
0393 // iFit = 3 : CTEQ66.00 (NLO, central member)
0394 // iFit = 4 : CT09MC1
0395 // iFit = 5 : CT09MC2
0396 // iFit = 6 : CT09MCS
0397 // Pomeron sets available (uses same .pds file format as CTEQ6pdf) :
0398 // iFit = 11: ACTWB14
0399 // iFit = 12: ACTWD14
0400 // iFit = 13: ACTWSG14
0401 // iFit = 14: ACTWD19
0402 
0403 class CTEQ6pdf : public PDF {
0404 
0405 public:
0406 
0407   // Constructor.
0408   CTEQ6pdf(int idBeamIn = 2212, int iFitIn = 1, double rescaleIn = 1.,
0409     string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
0410     : PDF(idBeamIn), doExtraPol(false), iFit(), order(), nQuark(), nfMx(),
0411     mxVal(), nX(), nT(), nG(), iGridX(), iGridQ(), iGridLX(), iGridLQ(),
0412     rescale(rescaleIn), lambda(), mQ(), qIni(), qMax(), tv(), xMin(), xv(),
0413     upd(), xvpow(), xMinEps(), xMaxEps(), qMinEps(), qMaxEps(), fVec(),
0414     tConst(), xConst(), dlx(), xLast(),
0415     qLast() {init( iFitIn, pdfdataPath, loggerPtr);}
0416 
0417   // Constructor with a stream.
0418   CTEQ6pdf(int idBeamIn, istream& is, bool isPdsGrid = false,
0419     Logger* loggerPtr = 0) : PDF(idBeamIn), doExtraPol(false), iFit(),
0420     order(), nQuark(), nfMx(), mxVal(), nX(), nT(), nG(), iGridX(), iGridQ(),
0421     iGridLX(), iGridLQ(), rescale(), lambda(), mQ(), qIni(), qMax(), tv(),
0422     xMin(), xv(), upd(), xvpow(), xMinEps(), xMaxEps(), qMinEps(), qMaxEps(),
0423     fVec(), tConst(), xConst(), dlx(), xLast(),
0424     qLast() { init( is, isPdsGrid, loggerPtr); }
0425 
0426   // Allow extrapolation beyond boundaries. This is optional.
0427   void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
0428 
0429 private:
0430 
0431   // Constants: could only be changed in the code itself.
0432   static const double EPSILON, XPOWER;
0433 
0434   // Data read in from grid file or set at initialization.
0435   bool   doExtraPol;
0436   int    iFit, order, nQuark, nfMx, mxVal, nX, nT, nG,
0437          iGridX, iGridQ, iGridLX, iGridLQ;
0438   double rescale, lambda, mQ[7], qIni, qMax, tv[26], xMin, xv[202], upd[57773],
0439          xvpow[202], xMinEps, xMaxEps, qMinEps, qMaxEps, fVec[5],
0440          tConst[9], xConst[9], dlx, xLast, qLast;
0441 
0442   // These inits do not overwrite PDF init (prevents Clang warnings).
0443   using PDF::init;
0444 
0445   // Initialization of data array.
0446   void init( int iFitIn, string pdfdataPath, Logger* loggerPtrIn);
0447 
0448   // Initialization through a stream.
0449   void init( istream& is, bool isPdsGrid, Logger* loggerPtrIn);
0450 
0451   // Update PDF values.
0452   void xfUpdate(int id, double x, double Q2) override;
0453 
0454   // Evaluate PDF of one flavour species.
0455   double parton6(int iParton, double x, double q);
0456 
0457   // Interpolation in grid.
0458   double polint4F(double xgrid[], double fgrid[], double xin);
0459 
0460 };
0461 
0462 //==========================================================================
0463 
0464 // SA Unresolved proton: equivalent photon spectrum from
0465 // V.M. Budnev, I.F. Ginzburg, G.V. Meledin and V.G. Serbo,
0466 // Phys. Rept. 15 (1974/1975) 181.
0467 
0468 class ProtonPoint : public PDF {
0469 
0470 public:
0471 
0472   // Constructor.
0473   ProtonPoint(int idBeamIn = 2212, Logger* loggerPtrIn = 0)
0474     : PDF(idBeamIn), loggerPtr(loggerPtrIn) {}
0475 
0476 private:
0477 
0478   // Stored value for PDF choice.
0479   static const double ALPHAEM, Q2MAX, Q20, A, B, C;
0480 
0481   // Update PDF values.
0482   void xfUpdate(int , double x, double Q2) override;
0483 
0484   // phi function from Q2 integration.
0485   double phiFunc(double x, double Q);
0486 
0487   // Pointer to logger.
0488   Logger* loggerPtr;
0489 
0490 };
0491 
0492 //==========================================================================
0493 
0494 // Gives the GRV 1992 pi+ (leading order) parton distribution function set
0495 // in parametrized form. Authors: Glueck, Reya and Vogt.
0496 
0497 class GRVpiL : public PDF {
0498 
0499 public:
0500 
0501   // Constructor.
0502   GRVpiL(int idBeamIn = 211, double vmdScaleIn = 1.) :
0503     PDF(idBeamIn) {vmdScale = vmdScaleIn;}
0504 
0505   // Allow for new rescaling factor of the PDF for VMD beams.
0506   void setVMDscale(double vmdScaleIn = 1.) override {vmdScale = vmdScaleIn;}
0507 
0508 private:
0509 
0510   // Rescaling of pion PDF for VMDs.
0511   double vmdScale;
0512 
0513   // Update PDF values.
0514   void xfUpdate(int , double x, double Q2) override;
0515 
0516 };
0517 
0518 //==========================================================================
0519 
0520 // Gives the GRS 1999 pi+ (leading order) parton distribution function set
0521 // in parametrized form. Authors: Glueck, Reya and Schienbein.
0522 
0523 class GRSpiL : public PDF {
0524 
0525 public:
0526 
0527   // Constructor.
0528   GRSpiL(int idBeamIn = 211, double vmdScaleIn = 1.) :
0529     PDF(idBeamIn) {vmdScale = vmdScaleIn;}
0530 
0531   // Allow for new rescaling factor of the PDF for VMD beams.
0532   void setVMDscale(double vmdScaleIn = 1.) override {vmdScale = vmdScaleIn;}
0533 
0534 private:
0535 
0536   // Rescaling of pion PDF for VMDs.
0537   double vmdScale;
0538 
0539   // Update PDF values.
0540   void xfUpdate(int , double x, double Q2) override;
0541 
0542 };
0543 
0544 //==========================================================================
0545 
0546 // Gives generic Q2-independent Pomeron PDF.
0547 
0548 class PomFix : public PDF {
0549 
0550 public:
0551 
0552   // Constructor.
0553   PomFix(int idBeamIn = 990, double PomGluonAIn = 0.,
0554     double PomGluonBIn = 0., double PomQuarkAIn = 0.,
0555     double PomQuarkBIn = 0., double PomQuarkFracIn = 0.,
0556     double PomStrangeSuppIn = 0.) : PDF(idBeamIn),
0557     PomGluonA(PomGluonAIn), PomGluonB(PomGluonBIn),
0558     PomQuarkA(PomQuarkAIn), PomQuarkB(PomQuarkBIn),
0559     PomQuarkFrac(PomQuarkFracIn), PomStrangeSupp(PomStrangeSuppIn),
0560     normGluon(), normQuark() { init(); }
0561 
0562 private:
0563 
0564   // Stored value for PDF choice.
0565   double PomGluonA, PomGluonB, PomQuarkA, PomQuarkB, PomQuarkFrac,
0566          PomStrangeSupp, normGluon, normQuark;
0567 
0568   // This init does not overwrite PDF init (prevents Clang warnings).
0569   using PDF::init;
0570 
0571   // Initialization of some constants.
0572   void init();
0573 
0574   // Update PDF values.
0575   void xfUpdate(int , double x, double) override;
0576 
0577 };
0578 
0579 //==========================================================================
0580 
0581 // The H1 2006 Fit A and Fit B Pomeron parametrization.
0582 // H1 Collaboration, A. Aktas et al., "Measurement and QCD Analysis of
0583 // the Diffractive Deep-Inelastic Scattering Cross Section at HERA",
0584 // DESY-06-049, Eur. Phys. J. C48 (2006) 715. e-Print: hep-ex/0606004.
0585 
0586 class PomH1FitAB : public PDF {
0587 
0588 public:
0589 
0590   // Constructor.
0591  PomH1FitAB(int idBeamIn = 990, int iFit = 1, double rescaleIn = 1.,
0592    string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
0593    : PDF(idBeamIn), doExtraPol(false), nx(), nQ2(), rescale(rescaleIn), xlow(),
0594     xupp(), dx(), Q2low(), Q2upp(), dQ2(), gluonGrid(), quarkGrid()
0595     { init( iFit, pdfdataPath, loggerPtr); }
0596 
0597   // Constructor with a stream.
0598  PomH1FitAB(int idBeamIn, double rescaleIn, istream& is,
0599    Logger* loggerPtr = 0) : PDF(idBeamIn), doExtraPol(false), nx(), nQ2(),
0600     rescale(rescaleIn), xlow(),xupp(), dx(), Q2low(), Q2upp(), dQ2(),
0601     gluonGrid(), quarkGrid() { init( is, loggerPtr); }
0602 
0603   // Allow extrapolation beyond boundaries. This is optional.
0604   void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
0605 
0606 private:
0607 
0608   // Limits for grid in x, in Q2, and data in (x, Q2).
0609   bool   doExtraPol;
0610   int    nx, nQ2;
0611   double rescale, xlow, xupp, dx, Q2low, Q2upp, dQ2;
0612   double gluonGrid[100][30];
0613   double quarkGrid[100][30];
0614 
0615   // These inits do not overwrite PDF init (prevents Clang warnings).
0616   using PDF::init;
0617 
0618   // Initialization of data array.
0619   void init( int iFit, string pdfdataPath, Logger* loggerPtr);
0620 
0621   // Initialization through a stream.
0622   void init( istream& is, Logger* loggerPtr);
0623 
0624   // Update PDF values.
0625   void xfUpdate(int , double x, double ) override;
0626 
0627 };
0628 
0629 //==========================================================================
0630 
0631 // The H1 2007 Jets Pomeron parametrization..
0632 // H1 Collaboration, A. Aktas et al., "Dijet Cross Sections and Parton
0633 // Densities in Diffractive DIS at HERA", DESY-07-115, Aug 2007. 33pp.
0634 // Published in JHEP 0710:042,2007. e-Print: arXiv:0708.3217 [hep-ex]
0635 
0636 class PomH1Jets : public PDF {
0637 
0638 public:
0639 
0640   // Constructor.
0641   PomH1Jets(int idBeamIn = 990, int iFit = 1, double rescaleIn = 1.,
0642     string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
0643     : PDF(idBeamIn), doExtraPol(false), rescale(rescaleIn), xGrid(), Q2Grid(),
0644     gluonGrid(), singletGrid(), charmGrid()
0645     {init( iFit, pdfdataPath, loggerPtr);}
0646 
0647   // Constructor with a stream.
0648   PomH1Jets(int idBeamIn, double rescaleIn, istream& is,
0649     Logger* loggerPtr = 0) : PDF(idBeamIn), doExtraPol(false),
0650     rescale(rescaleIn), xGrid(), Q2Grid(), gluonGrid(), singletGrid(),
0651     charmGrid() { init( is, loggerPtr); }
0652 
0653   // Allow extrapolation beyond boundaries. This is optional.
0654   void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
0655 
0656 private:
0657 
0658   // Arrays for grid in x, in Q2, and data in (x, Q2).
0659   bool   doExtraPol;
0660   double rescale;
0661   double xGrid[100];
0662   double Q2Grid[88];
0663   double gluonGrid[100][88];
0664   double singletGrid[100][88];
0665   double charmGrid[100][88];
0666 
0667   // These inits do not overwrite PDF init (prevents Clang warnings).
0668   using PDF::init;
0669 
0670   // Initialization of data array.
0671   void init( int iFit, string pdfdataPath, Logger* loggerPtr);
0672 
0673   // Initialization through a stream.
0674   void init( istream& is, Logger* loggerPtr);
0675 
0676   // Update PDF values.
0677   void xfUpdate(int id, double x, double ) override;
0678 
0679 };
0680 
0681 //==========================================================================
0682 
0683 // A proton masked as a Pomeron for use within the Heavy Ion machinery
0684 
0685 class PomHISASD : public PDF {
0686 
0687 public:
0688 
0689   // Basic constructor
0690   PomHISASD(int idBeamIn, PDFPtr ppdf, Settings & settings,
0691     Logger* loggerPtrIn = 0) : PDF(idBeamIn), pPDFPtr(ppdf),
0692     xPomNow(-1.0), hixpow(4.0), newfac(1.0) {
0693     loggerPtr = loggerPtrIn;
0694     hixpow = settings.parm("PDF:PomHixSupp");
0695     if ( settings.mode("Angantyr:SASDmode") == 3 ) newfac =
0696       log(settings.parm("Beams:eCM")/settings.parm("Diffraction:mMinPert"));
0697     if ( settings.mode("Angantyr:SASDmode") == 4 ) newfac = 0.0;
0698   }
0699 
0700   // Delete also the proton PDF
0701   ~PomHISASD() { }
0702 
0703   // (re-)Set the x_pomeron value.
0704   void xPom(double xpom = -1.0) override { xPomNow = xpom; }
0705 
0706 private:
0707 
0708   // The proton PDF.
0709   PDFPtr pPDFPtr;
0710 
0711   // The momentum fraction if the corresponding pomeron.
0712   double xPomNow;
0713 
0714   // The high-x suppression power.
0715   double hixpow;
0716 
0717   // Special options.
0718   double newfac;
0719 
0720   // Report possible errors.
0721   Logger* loggerPtr;
0722 
0723   // Update PDF values.
0724   void xfUpdate(int , double x, double Q2) override;
0725 
0726 };
0727 
0728 //==========================================================================
0729 
0730 // Gives electron (or muon, or tau) parton distribution.
0731 
0732 class Lepton : public PDF {
0733 
0734 public:
0735 
0736   // Constructor.
0737   Lepton(int idBeamIn = 11) : PDF(idBeamIn), m2Lep(), Q2maxGamma(),
0738     infoPtr(), rndmPtr() {}
0739 
0740   // Constructor with further info.
0741   Lepton(int idBeamIn, double Q2maxGammaIn, Info* infoPtrIn)
0742     : PDF(idBeamIn), m2Lep() { Q2maxGamma = Q2maxGammaIn;
0743     infoPtr = infoPtrIn; rndmPtr = infoPtrIn->rndmPtr; }
0744 
0745   // Sample the Q2 value.
0746   double sampleQ2gamma(double Q2min) override
0747     { return Q2min * pow(Q2maxGamma / Q2min, rndmPtr->flat()); }
0748 
0749 private:
0750 
0751   // Constants: could only be changed in the code itself.
0752   static const double ALPHAEM, ME, MMU, MTAU;
0753 
0754   // Update PDF values.
0755   void xfUpdate(int id, double x, double Q2) override;
0756 
0757   // The squared lepton mass, set at initialization.
0758   double m2Lep, Q2maxGamma;
0759 
0760   // Pointer to info, needed to get sqrt(s) to fix x_gamma limits.
0761   Info* infoPtr;
0762 
0763   // Pointer to random number generator, needed to sample Q2.
0764   Rndm* rndmPtr;
0765 
0766 };
0767 
0768 //==========================================================================
0769 
0770 // Gives electron (or other lepton) parton distribution when unresolved.
0771 
0772 class LeptonPoint : public PDF {
0773 
0774 public:
0775 
0776   // Constructor.
0777   LeptonPoint(int idBeamIn = 11) : PDF(idBeamIn) {}
0778 
0779 private:
0780 
0781   // Update PDF values in trivial way.
0782   void xfUpdate(int , double , double ) override {xlepton = 1; xgamma = 0.;}
0783 
0784 };
0785 
0786 //==========================================================================
0787 
0788 // Gives neutrino parton distribution when unresolved (only choice for now).
0789 // Note that the extra factor of 2 - wrt. charged leptons as there is no need
0790 // for spin averaging since neutrinos always lefthanded -  is taken care in
0791 // cross sections and not in the PDFs.
0792 
0793 class NeutrinoPoint : public PDF {
0794 
0795 public:
0796 
0797   // Constructor.
0798   NeutrinoPoint(int idBeamIn = 12) : PDF(idBeamIn) {}
0799 
0800 private:
0801 
0802   // Update PDF values in trivial way.
0803   void xfUpdate(int , double , double ) override {xlepton = 1; xgamma = 0.;}
0804 
0805 };
0806 
0807 //==========================================================================
0808 
0809 // Gives the CJKL leading order parton distribution function set
0810 // in parametrized form for the real photons. Authors: F.Cornet, P.Jankowski,
0811 // M.Krawczyk and A.Lorca, Phys. Rev. D68: 014010, 2003.
0812 
0813 class CJKL : public PDF {
0814 
0815 public:
0816 
0817   // Constructor. Needs the randon number generator to sample valence content.
0818   CJKL(int idBeamIn = 22, Rndm* rndmPtrIn = 0 ) : PDF(idBeamIn) {
0819     rndmPtr = rndmPtrIn; }
0820 
0821   // Functions to approximate pdfs for ISR.
0822   double gammaPDFxDependence(int id, double) override;
0823   double gammaPDFRefScale(int) override;
0824 
0825   // Set the valence content for photons.
0826   int sampleGammaValFlavor(double Q2) override;
0827 
0828   // The total x-integrated PDFs. Relevant for MPIs with photon beams.
0829   double xfIntegratedTotal(double Q2) override;
0830 
0831 private:
0832 
0833   // Parameters related to the fit.
0834   static const double ALPHAEM, Q02, Q2MIN, Q2REF, LAMBDA, MC, MB;
0835 
0836   // Pointer to random number generator used for valence sampling.
0837   Rndm *rndmPtr;
0838 
0839   // Update PDF values.
0840   void xfUpdate(int , double x, double Q2) override;
0841 
0842   // Functions for updating the point-like part.
0843   double pointlikeG(double x, double s);
0844   double pointlikeU(double x, double s);
0845   double pointlikeD(double x, double s);
0846   double pointlikeC(double x, double s, double Q2);
0847   double pointlikeB(double x, double s, double Q2);
0848 
0849   // Functions for updating the hadron-like part.
0850   double hadronlikeG(double x, double s);
0851   double hadronlikeSea(double x, double s);
0852   double hadronlikeVal(double x, double s);
0853   double hadronlikeC(double x, double s, double Q2);
0854   double hadronlikeB(double x, double s, double Q2);
0855 
0856 };
0857 
0858 //==========================================================================
0859 
0860 // Convolution with photon flux from leptons and photon PDFs.
0861 // Photon flux from equivalent photon approximation (EPA).
0862 // Contains a pointer to a photon PDF set and samples the
0863 // convolution integral event-by-event basis.
0864 // Includes also a overestimate for the PDF set in order to set up
0865 // the phase-space sampling correctly.
0866 
0867 class Lepton2gamma : public PDF {
0868 
0869 public:
0870 
0871   // Constructor.
0872   Lepton2gamma(int idBeamIn, double m2leptonIn, double Q2maxGamma,
0873     PDFPtr gammaPDFPtrIn, Info* infoPtrIn)
0874     : PDF(idBeamIn), m2lepton(m2leptonIn), Q2max(Q2maxGamma), xGm(),
0875     sampleXgamma(true), gammaPDFPtr(gammaPDFPtrIn),rndmPtr(infoPtrIn->rndmPtr),
0876     infoPtr(infoPtrIn) { hasGammaInLepton = true; }
0877 
0878   // Overload the member function definitions where relevant.
0879   void xfUpdate(int id, double x, double Q2) override;
0880   double xGamma() override { return xGm; }
0881   double xfMax(int id, double x, double Q2) override;
0882   double xfSame(int id, double x, double Q2) override;
0883 
0884   // Sample the Q2 value.
0885   double sampleQ2gamma(double Q2min) override
0886     { return Q2min * pow(Q2max / Q2min, rndmPtr->flat()); }
0887 
0888 private:
0889 
0890   // Parameters for convolution.
0891   static const double ALPHAEM, Q2MIN;
0892   double m2lepton, Q2max, xGm;
0893 
0894   // Sample new value for x_gamma.
0895   bool sampleXgamma;
0896 
0897   // Photon PDFs which the photon flux is convoluted with.
0898   PDFPtr gammaPDFPtr;
0899 
0900   // Pointer to random number generator used for sampling x_gamma.
0901   Rndm* rndmPtr;
0902 
0903   // Pointer to info, needed to get sqrt(s) to fix x_gamma limits.
0904   Info* infoPtr;
0905 
0906 };
0907 
0908 //==========================================================================
0909 
0910 // Gives photon parton distribution when unresolved.
0911 
0912 class GammaPoint : public PDF {
0913 
0914 public:
0915 
0916   // Constructor.
0917   GammaPoint(int idBeamIn = 22) : PDF(idBeamIn) {}
0918 
0919 private:
0920 
0921   // Update PDF values in trivial way.
0922   void xfUpdate(int , double , double ) override { xgamma = 1.;}
0923 
0924 };
0925 
0926 //==========================================================================
0927 
0928 // Unresolved proton: equivalent photon spectrum according
0929 // to the approximation by Drees and Zeppenfeld,
0930 // Phys.Rev. D39 (1989) 2536.
0931 
0932 class Proton2gammaDZ : public PDF {
0933 
0934 public:
0935 
0936   // Constructor.
0937   Proton2gammaDZ(int idBeamIn = 2212) : PDF(idBeamIn) {}
0938 
0939 private:
0940 
0941   // Stored parameters.
0942   static const double ALPHAEM, Q20;
0943 
0944   // Update PDF values.
0945   void xfUpdate(int , double x, double Q2) override;
0946   double fluxQ2dependence(double Q2) override;
0947 
0948 };
0949 
0950 //==========================================================================
0951 
0952 // Unresolved nucleus: equivalent photon approximation
0953 // for impact parameter integrated flux according to standard
0954 // form introduced in J.D. Jackson, Classical Electrodynamics,
0955 // 2nd edition, John Wiley & Sons (1975).
0956 
0957 class Nucleus2gamma : public PDF {
0958 
0959 public:
0960 
0961   // Constructor.
0962   Nucleus2gamma(int idBeamIn, double bMinIn, double mNucleonIn) :
0963     PDF(idBeamIn), a(), z(), bMin(bMinIn), mNucleon(mNucleonIn)
0964     { initNucleus(idBeamIn); }
0965 
0966 private:
0967 
0968   // Stored constant parameters.
0969   static const double ALPHAEM;
0970 
0971   // Initialize flux parameters.
0972   void initNucleus(int idBeamIn);
0973 
0974   // Update PDF values.
0975   void xfUpdate(int , double x, double Q2) override;
0976 
0977   // Mass number and electric charge.
0978   int a, z;
0979 
0980   // Minimum impact parameter for integration and per-nucleon mass.
0981   double bMin, mNucleon;
0982 
0983 };
0984 
0985 //==========================================================================
0986 
0987 // Equivalent photon approximation for sampling with external photon flux.
0988 
0989 class EPAexternal : public PDF {
0990 
0991 public:
0992 
0993   // Constructor.
0994   EPAexternal(int idBeamIn, double m2In, PDFPtr gammaFluxPtrIn,
0995     PDFPtr gammaPDFPtrIn, Info* infoPtrIn, Logger* loggerPtrIn = 0)
0996     : PDF(idBeamIn), m2(m2In), Q2max(), Q2min(), xMax(), xMin(), xHadr(),
0997     norm(), xPow(), xCut(), norm1(), norm2(), integral1(), integral2(),
0998     bmhbarc(), approxMode(0), isLHA(false), gammaFluxPtr(gammaFluxPtrIn),
0999     gammaPDFPtr(gammaPDFPtrIn), infoPtr(infoPtrIn),
1000     rndmPtr(infoPtrIn->rndmPtr), settingsPtr(infoPtrIn->settingsPtr),
1001     loggerPtr(loggerPtrIn) { hasGammaInLepton = true; init(); }
1002 
1003   // Update PDFs.
1004   void xfUpdate(int , double x, double Q2) override;
1005 
1006   // External flux and photon PDFs, and approximated flux for sampling.
1007   double xfFlux(int id, double x, double Q2 = 1.) override;
1008   double xfGamma(int id, double x, double Q2) override;
1009   double xfApprox(int id, double x, double Q2) override;
1010   double intFluxApprox() override;
1011 
1012   // This derived class use approximated flux for sampling.
1013   bool hasApproxGammaFlux() override { return true; }
1014 
1015   // Kinematics.
1016   double getXmin()  override { return xMin; }
1017   double getXhadr() override { return xHadr; }
1018 
1019   // Sampling of the x and Q2 according to differential flux.
1020   double sampleXgamma(double xMinIn) override;
1021   double sampleQ2gamma(double Q2minIn) override;
1022 
1023 private:
1024 
1025   // This init does not overwrite PDF init (prevents Clang warnings).
1026   using PDF::init;
1027 
1028   // Initialization.
1029   void init();
1030 
1031   // Kinematics.
1032   static const double ALPHAEM;
1033   double m2, Q2max, Q2min, xMax, xMin, xHadr, norm, xPow, xCut,
1034          norm1, norm2, integral1, integral2, bmhbarc;
1035   int    approxMode;
1036   bool   isLHA;
1037 
1038   // Photon Flux and PDF.
1039   PDFPtr gammaFluxPtr;
1040   PDFPtr gammaPDFPtr;
1041 
1042   // Pointer to info, needed to get sqrt(s) to fix x_gamma limits.
1043   Info* infoPtr;
1044 
1045   // Pointer to random number generator used for sampling x_gamma.
1046   Rndm* rndmPtr;
1047 
1048   // Pointer to settings to get Q2max.
1049   Settings* settingsPtr;
1050 
1051   // Pointer to logger.
1052   Logger* loggerPtr;
1053 
1054 };
1055 
1056 //==========================================================================
1057 
1058 // A derived class for nuclear PDFs. Needs a pointer for (free) proton PDFs.
1059 
1060 class nPDF : public PDF {
1061 
1062 public:
1063 
1064   // Constructor.
1065   nPDF(int idBeamIn = 2212, PDFPtr protonPDFPtrIn = 0) : PDF(idBeamIn), ruv(),
1066     rdv(), ru(), rd(), rs(), rc(), rb(), rg(), a(), z(), za(), na(),
1067     protonPDFPtr() { initNPDF(idBeamIn, protonPDFPtrIn); }
1068 
1069   // Update parton densities.
1070   void xfUpdate(int id, double x, double Q2) override;
1071 
1072   // Update nuclear modifications.
1073   virtual void rUpdate(int, double, double) = 0;
1074 
1075   // Initialize the nPDF-related members.
1076   void initNPDF(int idBeamIn, PDFPtr protonPDFPtrIn = 0);
1077 
1078   // Return the number of protons and nucleons.
1079   int getA() {return a;}
1080   int getZ() {return z;}
1081 
1082   // Set (and reset) the ratio of protons to nucleons to study nuclear
1083   // modifications of protons (= 1.0) and neutrons (= 0.0). By default Z/A.
1084   void setMode(double zaIn) { za = zaIn; na = 1. - za; }
1085   void resetMode() { za = double(z)/double(a); na = double(a-z)/double(a); }
1086 
1087 protected:
1088 
1089   // The nuclear modifications for each flavour, modified by derived nPDF
1090   // classes.
1091   double ruv, rdv, ru, rd, rs, rc, rb, rg;
1092 
1093 private:
1094 
1095   // The nuclear mass number and number of protons (charge) and normalized
1096   // number of protons and neutrons.
1097   int a, z;
1098   double za, na;
1099 
1100   // Pointer to (free) proton PDF.
1101   PDFPtr protonPDFPtr;
1102 
1103 };
1104 
1105 //==========================================================================
1106 
1107 // Isospin modification with nuclear beam, i.e. no other modifications
1108 // but correct number of protons and neutrons.
1109 
1110 class Isospin : public nPDF {
1111 
1112 public:
1113 
1114   // Constructor.
1115   Isospin(int idBeamIn = 2212, PDFPtr protonPDFPtrIn = 0)
1116     : nPDF(idBeamIn, protonPDFPtrIn) {}
1117 
1118   // Only the Isospin effect so no need to do anything here.
1119   void rUpdate(int , double , double ) override {}
1120 };
1121 
1122 //==========================================================================
1123 
1124 // Nuclear modifications from EPS09 fit.
1125 
1126 class EPS09 : public nPDF {
1127 
1128 public:
1129 
1130   // Constructor.
1131   EPS09(int idBeamIn = 2212, int iOrderIn = 1, int iSetIn = 1,
1132     string pdfdataPath = "../share/Pythia8/pdfdata/",
1133     PDFPtr protonPDFPtrIn = 0, Logger* loggerPtrIn = 0)
1134     : nPDF(idBeamIn, protonPDFPtrIn), iSet(), iOrder(), grid(),
1135     loggerPtr(loggerPtrIn) { init(iOrderIn, iSetIn, pdfdataPath);}
1136 
1137   // Update parton densities.
1138   void rUpdate(int id, double x, double Q2) override;
1139 
1140   // Use other than central set to study uncertainties.
1141   void setErrorSet(int iSetIn) {iSet = iSetIn;}
1142 
1143 private:
1144 
1145   // Parameters related to the fit.
1146   static const double Q2MIN, Q2MAX, XMIN, XMAX, XCUT;
1147   static const int Q2STEPS, XSTEPS;
1148 
1149   // Set parameters and the grid.
1150   int iSet, iOrder;
1151   double grid[31][51][51][8];
1152 
1153   // Pointer to logger for possible error messages.
1154   Logger* loggerPtr;
1155 
1156   // This init does not overwrite PDF init (prevents Clang warnings).
1157   using PDF::init;
1158 
1159   // Initialize with given inputs.
1160   void init(int iOrderIn, int iSetIn, string pdfdataPath);
1161 
1162   // Interpolation algorithm.
1163   double polInt(double* fi, double* xi, int n, double x);
1164 };
1165 
1166 //==========================================================================
1167 
1168 // Nuclear modifications from EPPS16 fit.
1169 
1170 class EPPS16 : public nPDF {
1171 
1172 public:
1173 
1174   // Constructor.
1175   EPPS16(int idBeamIn = 2212, int iSetIn = 1,
1176     string pdfdataPath = "../share/Pythia8/pdfdata/",
1177     PDFPtr protonPDFPtrIn = 0, Logger* loggerPtrIn = 0)
1178     : nPDF(idBeamIn, protonPDFPtrIn), iSet(), grid(), logQ2min(),
1179     loglogQ2maxmin(), logX2min(), loggerPtr(loggerPtrIn)
1180     { init(iSetIn, pdfdataPath); }
1181 
1182   // Update parton densities.
1183   void rUpdate(int id, double x, double Q2) override;
1184 
1185   // Use other than central set to study uncertainties.
1186   void setErrorSet(int iSetIn) {iSet = iSetIn;}
1187 
1188 private:
1189 
1190   // Parameters related to the fit.
1191   static const double Q2MIN, Q2MAX, XMIN, XMAX, XCUT;
1192   static const int Q2STEPS, XSTEPS, NINTQ2, NINTX, NSETS;
1193 
1194   // Set parameters and the grid.
1195   int iSet;
1196   double grid[41][31][80][8];
1197   double logQ2min, loglogQ2maxmin, logX2min;
1198 
1199   // Pointer to logger.
1200   Logger* loggerPtr;
1201 
1202   // This init does not overwrite PDF init (prevents Clang warnings).
1203   using PDF::init;
1204 
1205   // Initialize with given inputs.
1206   void init(int iSetIn, string pdfdataPath);
1207 
1208   // Interpolation algorithm.
1209   double polInt(double* fi, double* xi, int n, double x);
1210 };
1211 
1212 //==========================================================================
1213 
1214 } // end namespace Pythia8
1215 
1216 #endif // Pythia8_PartonDistributions_H