PYTHIA  8.311
PartonDistributions.h
1 // PartonDistributions.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2024 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Header files for parton densities.
7 // PDF: base class.
8 // LHAGrid1: internal read and use files in the LHAPDF6 lhagrid1 format.
9 // GRV94L: GRV 94L parton densities.
10 // CTEQ5L: CTEQ 5L parton densities.
11 // MSTWpdf: MRST LO*, LO**, MSTW 2008 LO, NLO.
12 // CTEQ6pdf: CTEQ 6L, 6L1, 66, CT09 MC1, MC2, MCS.
13 // ProtonPoint: unresolved proton with equivalent photon spectrum.
14 // GRVpiL: GRV LO pion parton densities.
15 // PomFix: Q2-independent Pomeron parton densities.
16 // PomH1FitAB: H1 2006 Fit A and Fit B Pomeron PDFs.
17 // PomH1Jets: H1 2007 Jets Pomeron PDFs.
18 // PomHISASD: a proton masked as a Pomeron for heavy ions applications.
19 // Lepton: parton densities inside a lepton.
20 // LeptonPoint: an unresolved lepton (mainly dummy).
21 // NeutrinoPoint: an unresolved neutrino (mainly dummy).
22 // CJKL: CJKL parton densities for photons.
23 // Lepton2gamma: convolution of photon flux from leptons and photon PDFs.
24 // PhotonPoint: an unresolved photon.
25 // Proton2gammaDZ: Photon flux from protons according to Drees-Zeppenfeld.
26 // Nucleus2gamma: Photon flux from heavy nuclei.
27 // EPAexternal: approximated photon flux used for sampling of external flux.
28 // nPDF: a nuclear PDF, derived from a proton ditto.
29 // Isospin: isospin modification for nuclear pDF
30 // EPS09, EPPS16: nuclear modification factors.
31 
32 
33 #ifndef Pythia8_PartonDistributions_H
34 #define Pythia8_PartonDistributions_H
35 
36 #include "Pythia8/Basics.h"
37 #include "Pythia8/Info.h"
38 #include "Pythia8/MathTools.h"
39 #include "Pythia8/ParticleData.h"
40 #include "Pythia8/PythiaStdlib.h"
41 #include "Pythia8/SharedPointers.h"
42 
43 namespace Pythia8 {
44 
45 //==========================================================================
46 
47 // Base class for parton distribution functions.
48 
49 class PDF {
50 
51 public:
52 
53  // Constructor.
54  PDF(int idBeamIn = 2212) : idBeam(idBeamIn), idBeamAbs(abs(idBeam)),
55  idSav(9), xSav(-1), Q2Sav(-1.), isSet(true), isInit(false),
57 
58  // Virtual destructor.
59  virtual ~PDF() {}
60 
61  // Perform additional initialization (LHPADF only).
62  // Arguments are:
63  // int idBeamIn, string setName, int member, Logger* loggerPtr
64  virtual bool init(int, string, int, Logger*) {return true;}
65 
66  // Confirm that PDF has been set up (important for LHAPDF and H1 Pomeron).
67  bool isSetup() {return isSet;}
68 
69  // Switch to new beam particle identities; for similar hadrons only.
70  virtual void setBeamID(int idBeamIn) { idBeam = idBeamIn;
71  idBeamAbs = abs(idBeam); idSav = 9; xSav = -1.; Q2Sav = -1.;
73 
74  // Set valence content.
75  void resetValenceContent();
76  void setValenceContent(int idVal1In, int idVal2In, int idVal3In) {
77  idVal1 = idVal1In; idVal2 = idVal2In; idVal3 = idVal3In;}
78 
79  // Allow extrapolation beyond boundaries. This is optional.
80  virtual void setExtrapolate(bool) {}
81 
82  // Read out parton density.
83  double xf(int id, double x, double Q2);
84 
85  // Read out valence and sea part of parton densities.
86  double xfVal(int id, double x, double Q2);
87  double xfSea(int id, double x, double Q2);
88 
89  // Check whether x and Q2 values fall inside the fit bounds (LHAPDF6 only).
90  virtual bool insideBounds(double, double) {return true;}
91 
92  // Access the running alpha_s of a PDF set (LHAPDF6 only).
93  virtual double alphaS(double) { return 1.;}
94 
95  // Return quark masses used in the PDF fit (LHAPDF6 only).
96  virtual double mQuarkPDF(int) { return -1.;}
97 
98  // Return number of members of this PDF family (LHAPDF6 only).
99  virtual int nMembers() { return 1;}
100 
101  // Error envelope from PDF uncertainty.
102  struct PDFEnvelope {
103  double centralPDF, errplusPDF, errminusPDF, errsymmPDF, scalePDF;
104  vector<double> pdfMemberVars;
105  PDFEnvelope() : centralPDF(-1.0), errplusPDF(0.0), errminusPDF(0.0),
106  errsymmPDF(0.0), scalePDF(-1.0), pdfMemberVars(0.0) {};
107  };
108 
109  // Calculate PDF envelope.
110  virtual void calcPDFEnvelope(int, double, double, int) {}
111  virtual void calcPDFEnvelope(pair<int,int>, pair<double,double>, double,
112  int) {}
113  virtual PDFEnvelope getPDFEnvelope() { return PDFEnvelope(); }
114 
115  // Approximate photon PDFs by decoupling the scale and x-dependence.
116  virtual double gammaPDFxDependence(int, double) { return 0.; }
117 
118  // Provide the reference scale for logarithmic Q^2 evolution for photons.
119  virtual double gammaPDFRefScale(int) { return 0.; }
120 
121  // Sample the valence content for photons.
122  virtual int sampleGammaValFlavor(double) { return 0.; }
123 
124  // The total x-integrated PDFs. Relevant for MPIs with photon beams.
125  virtual double xfIntegratedTotal(double) { return 0.; }
126 
127  // Return the sampled value for x_gamma.
128  virtual double xGamma() { return 1.; }
129 
130  // Keep track of pomeron momentum fraction.
131  virtual void xPom(double = -1.0) {}
132 
133  // Return accurate and approximated photon fluxes and PDFs.
134  virtual double xfFlux(int , double , double ) { return 0.; }
135  virtual double xfApprox(int , double , double ) { return 0.; }
136  virtual double xfGamma(int , double , double ) { return 0.; }
137  virtual double intFluxApprox() { return 0.; }
138  virtual bool hasApproxGammaFlux() { return false; }
139 
140  // Return the kinematical limits and sample Q2 and x.
141  virtual double getXmin() { return 0.; }
142  virtual double getXhadr() { return 0.; }
143  virtual double sampleXgamma(double ) { return 0.; }
144  virtual double sampleQ2gamma(double ) { return 0.; }
145  virtual double fluxQ2dependence(double ) { return 0.; }
146 
147  // Normal PDFs unless gamma inside lepton -> an overestimate for sampling.
148  virtual double xfMax(int id, double x, double Q2) { return xf( id, x, Q2); }
149 
150  // Normal PDFs unless gamma inside lepton -> Do not sample x_gamma.
151  virtual double xfSame(int id, double x, double Q2) { return xf( id, x, Q2); }
152 
153  // Allow for new scaling factor for VMD PDFs.
154  virtual void setVMDscale(double = 1.) {}
155 
156  // Return if s/sbar, c/cbar, and b/bbar PDFs are symmetric.
157  bool sSymmetric() const { return sSymmetricSave; }
158  bool cSymmetric() const { return cSymmetricSave; }
159  bool bSymmetric() const { return bSymmetricSave; }
160 
161  // Set s/sbar, c/cbar, and b/bbar PDFs symmetric.
162  void sSymmetric(bool sSymmetricIn) { sSymmetricSave = sSymmetricIn; }
163  void cSymmetric(bool cSymmetricIn) { cSymmetricSave = cSymmetricIn; }
164  void bSymmetric(bool bSymmetricIn) { bSymmetricSave = bSymmetricIn; }
165 
166 protected:
167 
168  // Store relevant quantities.
169  int idBeam, idBeamAbs, idSav, idVal1, idVal2, idVal3;
170  double xSav, Q2Sav;
171  // Stored quantities.
172  double xu, xd, xs, xubar, xdbar, xsbar, xc, xb, xcbar, xbbar,
173  xg, xlepton, xgamma;
174  bool isSet, isInit;
175 
176  // For hadrons, beamType defines special cases and determines how
177  // to handle isospin symmetries.
178  // 1: no rearrangement (e.g. p, Sigma+, Omega-, pi+)
179  // -1: switch u <-> d (e.g. n, Sigma-, Xi-, K0)
180  // 0: take average of u and d (e.g. Sigma0, Lambda0)
181  // 2/-2: Delta++/Delta-
182  // 111: pi0-like special case (pi0, rho0, omega, etc.)
183  // 221: Other diagonal meson cases (eta, eta', phi, J/psi, Upsilon, etc.)
184  // 130: K_S,L special cases
185  int beamType;
186 
187  // True if a photon beam inside a lepton beam, otherwise set false.
189 
190  // Whether to treat flavoured PDFs as symmetric, for efficiency.
191  bool sSymmetricSave = false;
192  bool cSymmetricSave = true, bSymmetricSave = true;
193 
194  // Update parton densities.
195  virtual void xfUpdate(int id, double x, double Q2) = 0;
196 
197  // Small routine for error printout, depending on loggerPtr existing or not.
198  void printErr(string loc, string errMsg, Logger* loggerPtr = nullptr) {
199  if (loggerPtr) loggerPtr->errorMsg(loc, errMsg);
200  else cout << "Error in " + loc + ": " + errMsg << endl;
201  }
202 
203  // Get the raw stored value for the quark variable corresponding to the id.
204  double xfRaw(int id) const;
205 
206  // Check whether the specified id is a valence quark.
207  bool isValence(int id) const {
208  return id != 0 && (id == idVal1 || id == idVal2 || id == idVal3); }
209 
210 };
211 
212 //==========================================================================
213 
214 // The LHAGrid1 can be used to read files in the LHAPDF6 lhagrid1 format,
215 // assuming that the same x grid is used for all Q subgrids.
216 // Results are not identical with LHAPDF6, owing to different interpolation.
217 
218 class LHAGrid1 : public PDF {
219 
220 public:
221 
222  // Constructor.
223  LHAGrid1(int idBeamIn = 2212, string pdfWord = "void",
224  string xmlPath = "../share/Pythia8/xmldoc/", Logger* loggerPtr = 0)
225  : PDF(idBeamIn), doExtraPol(false), nx(), nq(), nqSub(), xMin(), xMax(),
226  qMin(), qMax(), pdfVal(), pdfGrid(), pdfSlope(nullptr) {
227  init( pdfWord, xmlPath, loggerPtr); };
228 
229  // Constructor with a stream.
230  LHAGrid1(int idBeamIn, istream& is, Logger* loggerPtr = 0)
231  : PDF(idBeamIn), doExtraPol(false), nx(), nq(), nqSub(), xMin(), xMax(),
232  qMin(), qMax(), pdfVal(), pdfGrid(), pdfSlope(nullptr) {
233  init( is, loggerPtr); };
234 
235  // Destructor.
236  ~LHAGrid1() { for (int iid = 0; iid < 12; ++iid) {
237  for (int iq = 0; iq < nq; ++iq) delete[] pdfGrid[iid][iq];
238  delete[] pdfGrid[iid]; }
239  if (pdfSlope) { for (int iid = 0; iid < 12; ++iid) delete[] pdfSlope[iid];
240  delete[] pdfSlope;} };
241 
242  // Allow extrapolation beyond boundaries. This is optional.
243  void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
244 
245 private:
246 
247  // Variables to be set during code initialization.
248  bool doExtraPol;
249  int nx, nq, nqSub;
250  vector<int> nqSum;
251  double xMin, xMax, qMin, qMax, pdfVal[12];
252  vector<double> xGrid, lnxGrid, qGrid, lnqGrid, qDiv;
253  double** pdfGrid[12];
254  double** pdfSlope;
255 
256  // These inits do not overwrite PDF init (prevents Clang warnings).
257  using PDF::init;
258 
259  // Initialization of data array.
260  void init( string pdfSet, string pdfdataPath, Logger* loggerPtr);
261 
262  // Initialization through a stream.
263  void init( istream& is, Logger* loggerPtr);
264 
265  // Update PDF values.
266  void xfUpdate(int id, double x, double Q2) override;
267 
268  // Interpolation in the grid for a given PDF flavour.
269  void xfxevolve(double x, double Q2);
270 
271 };
272 
273 //==========================================================================
274 
275 // Gives the GRV 94L (leading order) parton distribution function set
276 // in parametrized form. Authors: M. Glueck, E. Reya and A. Vogt.
277 
278 class GRV94L : public PDF {
279 
280 public:
281 
282  // Constructor.
283  GRV94L(int idBeamIn = 2212) : PDF(idBeamIn) {}
284 
285 private:
286 
287  // Update PDF values.
288  void xfUpdate(int , double x, double Q2) override;
289 
290  // Auxiliary routines used during the updating.
291  double grvv (double x, double n, double ak, double bk, double a,
292  double b, double c, double d);
293  double grvw (double x, double s, double al, double be, double ak,
294  double bk, double a, double b, double c, double d, double e, double es);
295  double grvs (double x, double s, double sth, double al, double be,
296  double ak, double ag, double b, double d, double e, double es);
297 
298 };
299 
300 //==========================================================================
301 
302 // Gives the CTEQ 5L (leading order) parton distribution function set
303 // in parametrized form. Parametrization by J. Pumplin. Authors: CTEQ.
304 
305 class CTEQ5L : public PDF {
306 
307 public:
308 
309  // Constructor.
310  CTEQ5L(int idBeamIn = 2212) : PDF(idBeamIn) {}
311 
312 private:
313 
314  // Update PDF values.
315  void xfUpdate(int , double x, double Q2) override;
316 
317 };
318 
319 //==========================================================================
320 
321 // The MSTWpdf class.
322 // MRST LO*(*) and MSTW 2008 PDF's, specifically the LO one.
323 // Original C++ version by Jeppe Andersen.
324 // Modified by Graeme Watt <watt(at)hep.ucl.ac.uk>.
325 // Sets available:
326 // iFit = 1 : MRST LO* (2007).
327 // iFit = 2 : MRST LO** (2008).
328 // iFit = 3 : MSTW 2008 LO, central member.
329 // iFit = 4 : MSTW 2008 NLO, central member. (Warning!)
330 
331 class MSTWpdf : public PDF {
332 
333 public:
334 
335  // Constructor.
336  MSTWpdf(int idBeamIn = 2212, int iFitIn = 1,
337  string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
338  : PDF(idBeamIn), iFit(), alphaSorder(), alphaSnfmax(), mCharm(), mBottom(),
339  alphaSQ0(), alphaSMZ(), distance(), tolerance(), xx(), qq(),
340  c() {init( iFitIn, pdfdataPath, loggerPtr);}
341 
342  // Constructor with a stream.
343  MSTWpdf(int idBeamIn, istream& is, Logger* loggerPtr = 0)
344  : PDF(idBeamIn), iFit(), alphaSorder(), alphaSnfmax(), mCharm(), mBottom(),
345  alphaSQ0(), alphaSMZ(), distance(), tolerance(), xx(), qq(),
346  c() {init( is, loggerPtr);}
347 
348 private:
349 
350  // Constants: could only be changed in the code itself.
351  static const int np, nx, nq, nqc0, nqb0;
352  static const double xmin, xmax, qsqmin, qsqmax, xxInit[65], qqInit[49];
353 
354  // Data read in from grid file or set at initialization.
355  int iFit, alphaSorder, alphaSnfmax;
356  double mCharm, mBottom, alphaSQ0, alphaSMZ, distance, tolerance,
357  xx[65], qq[49], c[13][64][48][5][5];
358 
359  // These inits do not overwrite PDF init (prevents Clang warnings).
360  using PDF::init;
361 
362  // Initialization of data array.
363  void init( int iFitIn, string pdfdataPath, Logger* loggerPtr);
364 
365  // Initialization through a stream.
366  void init( istream& is, Logger* loggerPtr);
367 
368  // Update PDF values.
369  void xfUpdate(int , double x, double Q2) override;
370 
371  // Evaluate PDF of one flavour species.
372  double parton(int flavour,double x,double q);
373  double parton_interpolate(int flavour,double xxx,double qqq);
374  double parton_extrapolate(int flavour,double xxx,double qqq);
375 
376  // Auxiliary routines for evaluation.
377  int locate(double xx[],int n,double x);
378  double polderivative1(double x1, double x2, double x3, double y1,
379  double y2, double y3);
380  double polderivative2(double x1, double x2, double x3, double y1,
381  double y2, double y3);
382  double polderivative3(double x1, double x2, double x3, double y1,
383  double y2, double y3);
384 
385 };
386 
387 //==========================================================================
388 
389 // The CTEQ6pdf class.
390 // Proton sets available:
391 // iFit = 1 : CTEQ6L
392 // iFit = 2 : CTEQ6L1
393 // iFit = 3 : CTEQ66.00 (NLO, central member)
394 // iFit = 4 : CT09MC1
395 // iFit = 5 : CT09MC2
396 // iFit = 6 : CT09MCS
397 // Pomeron sets available (uses same .pds file format as CTEQ6pdf) :
398 // iFit = 11: ACTWB14
399 // iFit = 12: ACTWD14
400 // iFit = 13: ACTWSG14
401 // iFit = 14: ACTWD19
402 
403 class CTEQ6pdf : public PDF {
404 
405 public:
406 
407  // Constructor.
408  CTEQ6pdf(int idBeamIn = 2212, int iFitIn = 1, double rescaleIn = 1.,
409  string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
410  : PDF(idBeamIn), doExtraPol(false), iFit(), order(), nQuark(), nfMx(),
411  mxVal(), nX(), nT(), nG(), iGridX(), iGridQ(), iGridLX(), iGridLQ(),
412  rescale(rescaleIn), lambda(), mQ(), qIni(), qMax(), tv(), xMin(), xv(),
413  upd(), xvpow(), xMinEps(), xMaxEps(), qMinEps(), qMaxEps(), fVec(),
414  tConst(), xConst(), dlx(), xLast(),
415  qLast() {init( iFitIn, pdfdataPath, loggerPtr);}
416 
417  // Constructor with a stream.
418  CTEQ6pdf(int idBeamIn, istream& is, bool isPdsGrid = false,
419  Logger* loggerPtr = 0) : PDF(idBeamIn), doExtraPol(false), iFit(),
420  order(), nQuark(), nfMx(), mxVal(), nX(), nT(), nG(), iGridX(), iGridQ(),
421  iGridLX(), iGridLQ(), rescale(), lambda(), mQ(), qIni(), qMax(), tv(),
422  xMin(), xv(), upd(), xvpow(), xMinEps(), xMaxEps(), qMinEps(), qMaxEps(),
423  fVec(), tConst(), xConst(), dlx(), xLast(),
424  qLast() { init( is, isPdsGrid, loggerPtr); }
425 
426  // Allow extrapolation beyond boundaries. This is optional.
427  void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
428 
429 private:
430 
431  // Constants: could only be changed in the code itself.
432  static const double EPSILON, XPOWER;
433 
434  // Data read in from grid file or set at initialization.
435  bool doExtraPol;
436  int iFit, order, nQuark, nfMx, mxVal, nX, nT, nG,
437  iGridX, iGridQ, iGridLX, iGridLQ;
438  double rescale, lambda, mQ[7], qIni, qMax, tv[26], xMin, xv[202], upd[57773],
439  xvpow[202], xMinEps, xMaxEps, qMinEps, qMaxEps, fVec[5],
440  tConst[9], xConst[9], dlx, xLast, qLast;
441 
442  // These inits do not overwrite PDF init (prevents Clang warnings).
443  using PDF::init;
444 
445  // Initialization of data array.
446  void init( int iFitIn, string pdfdataPath, Logger* loggerPtrIn);
447 
448  // Initialization through a stream.
449  void init( istream& is, bool isPdsGrid, Logger* loggerPtrIn);
450 
451  // Update PDF values.
452  void xfUpdate(int id, double x, double Q2) override;
453 
454  // Evaluate PDF of one flavour species.
455  double parton6(int iParton, double x, double q);
456 
457  // Interpolation in grid.
458  double polint4F(double xgrid[], double fgrid[], double xin);
459 
460 };
461 
462 //==========================================================================
463 
464 // SA Unresolved proton: equivalent photon spectrum from
465 // V.M. Budnev, I.F. Ginzburg, G.V. Meledin and V.G. Serbo,
466 // Phys. Rept. 15 (1974/1975) 181.
467 
468 class ProtonPoint : public PDF {
469 
470 public:
471 
472  // Constructor.
473  ProtonPoint(int idBeamIn = 2212, Logger* loggerPtrIn = 0)
474  : PDF(idBeamIn), loggerPtr(loggerPtrIn) {}
475 
476 private:
477 
478  // Stored value for PDF choice.
479  static const double ALPHAEM, Q2MAX, Q20, A, B, C;
480 
481  // Update PDF values.
482  void xfUpdate(int , double x, double Q2) override;
483 
484  // phi function from Q2 integration.
485  double phiFunc(double x, double Q);
486 
487  // Pointer to logger.
488  Logger* loggerPtr;
489 
490 };
491 
492 //==========================================================================
493 
494 // Gives the GRV 1992 pi+ (leading order) parton distribution function set
495 // in parametrized form. Authors: Glueck, Reya and Vogt.
496 
497 class GRVpiL : public PDF {
498 
499 public:
500 
501  // Constructor.
502  GRVpiL(int idBeamIn = 211, double vmdScaleIn = 1.) :
503  PDF(idBeamIn) {vmdScale = vmdScaleIn;}
504 
505  // Allow for new rescaling factor of the PDF for VMD beams.
506  void setVMDscale(double vmdScaleIn = 1.) override {vmdScale = vmdScaleIn;}
507 
508 private:
509 
510  // Rescaling of pion PDF for VMDs.
511  double vmdScale;
512 
513  // Update PDF values.
514  void xfUpdate(int , double x, double Q2) override;
515 
516 };
517 
518 //==========================================================================
519 
520 // Gives the GRS 1999 pi+ (leading order) parton distribution function set
521 // in parametrized form. Authors: Glueck, Reya and Schienbein.
522 
523 class GRSpiL : public PDF {
524 
525 public:
526 
527  // Constructor.
528  GRSpiL(int idBeamIn = 211, double vmdScaleIn = 1.) :
529  PDF(idBeamIn) {vmdScale = vmdScaleIn;}
530 
531  // Allow for new rescaling factor of the PDF for VMD beams.
532  void setVMDscale(double vmdScaleIn = 1.) override {vmdScale = vmdScaleIn;}
533 
534 private:
535 
536  // Rescaling of pion PDF for VMDs.
537  double vmdScale;
538 
539  // Update PDF values.
540  void xfUpdate(int , double x, double Q2) override;
541 
542 };
543 
544 //==========================================================================
545 
546 // Gives generic Q2-independent Pomeron PDF.
547 
548 class PomFix : public PDF {
549 
550 public:
551 
552  // Constructor.
553  PomFix(int idBeamIn = 990, double PomGluonAIn = 0.,
554  double PomGluonBIn = 0., double PomQuarkAIn = 0.,
555  double PomQuarkBIn = 0., double PomQuarkFracIn = 0.,
556  double PomStrangeSuppIn = 0.) : PDF(idBeamIn),
557  PomGluonA(PomGluonAIn), PomGluonB(PomGluonBIn),
558  PomQuarkA(PomQuarkAIn), PomQuarkB(PomQuarkBIn),
559  PomQuarkFrac(PomQuarkFracIn), PomStrangeSupp(PomStrangeSuppIn),
560  normGluon(), normQuark() { init(); }
561 
562 private:
563 
564  // Stored value for PDF choice.
565  double PomGluonA, PomGluonB, PomQuarkA, PomQuarkB, PomQuarkFrac,
566  PomStrangeSupp, normGluon, normQuark;
567 
568  // This init does not overwrite PDF init (prevents Clang warnings).
569  using PDF::init;
570 
571  // Initialization of some constants.
572  void init();
573 
574  // Update PDF values.
575  void xfUpdate(int , double x, double) override;
576 
577 };
578 
579 //==========================================================================
580 
581 // The H1 2006 Fit A and Fit B Pomeron parametrization.
582 // H1 Collaboration, A. Aktas et al., "Measurement and QCD Analysis of
583 // the Diffractive Deep-Inelastic Scattering Cross Section at HERA",
584 // DESY-06-049, Eur. Phys. J. C48 (2006) 715. e-Print: hep-ex/0606004.
585 
586 class PomH1FitAB : public PDF {
587 
588 public:
589 
590  // Constructor.
591  PomH1FitAB(int idBeamIn = 990, int iFit = 1, double rescaleIn = 1.,
592  string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
593  : PDF(idBeamIn), doExtraPol(false), nx(), nQ2(), rescale(rescaleIn), xlow(),
594  xupp(), dx(), Q2low(), Q2upp(), dQ2(), gluonGrid(), quarkGrid()
595  { init( iFit, pdfdataPath, loggerPtr); }
596 
597  // Constructor with a stream.
598  PomH1FitAB(int idBeamIn, double rescaleIn, istream& is,
599  Logger* loggerPtr = 0) : PDF(idBeamIn), doExtraPol(false), nx(), nQ2(),
600  rescale(rescaleIn), xlow(),xupp(), dx(), Q2low(), Q2upp(), dQ2(),
601  gluonGrid(), quarkGrid() { init( is, loggerPtr); }
602 
603  // Allow extrapolation beyond boundaries. This is optional.
604  void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
605 
606 private:
607 
608  // Limits for grid in x, in Q2, and data in (x, Q2).
609  bool doExtraPol;
610  int nx, nQ2;
611  double rescale, xlow, xupp, dx, Q2low, Q2upp, dQ2;
612  double gluonGrid[100][30];
613  double quarkGrid[100][30];
614 
615  // These inits do not overwrite PDF init (prevents Clang warnings).
616  using PDF::init;
617 
618  // Initialization of data array.
619  void init( int iFit, string pdfdataPath, Logger* loggerPtr);
620 
621  // Initialization through a stream.
622  void init( istream& is, Logger* loggerPtr);
623 
624  // Update PDF values.
625  void xfUpdate(int , double x, double ) override;
626 
627 };
628 
629 //==========================================================================
630 
631 // The H1 2007 Jets Pomeron parametrization..
632 // H1 Collaboration, A. Aktas et al., "Dijet Cross Sections and Parton
633 // Densities in Diffractive DIS at HERA", DESY-07-115, Aug 2007. 33pp.
634 // Published in JHEP 0710:042,2007. e-Print: arXiv:0708.3217 [hep-ex]
635 
636 class PomH1Jets : public PDF {
637 
638 public:
639 
640  // Constructor.
641  PomH1Jets(int idBeamIn = 990, int iFit = 1, double rescaleIn = 1.,
642  string pdfdataPath = "../share/Pythia8/pdfdata/", Logger* loggerPtr = 0)
643  : PDF(idBeamIn), doExtraPol(false), rescale(rescaleIn), xGrid(), Q2Grid(),
644  gluonGrid(), singletGrid(), charmGrid()
645  {init( iFit, pdfdataPath, loggerPtr);}
646 
647  // Constructor with a stream.
648  PomH1Jets(int idBeamIn, double rescaleIn, istream& is,
649  Logger* loggerPtr = 0) : PDF(idBeamIn), doExtraPol(false),
650  rescale(rescaleIn), xGrid(), Q2Grid(), gluonGrid(), singletGrid(),
651  charmGrid() { init( is, loggerPtr); }
652 
653  // Allow extrapolation beyond boundaries. This is optional.
654  void setExtrapolate(bool doExtraPolIn) override {doExtraPol = doExtraPolIn;}
655 
656 private:
657 
658  // Arrays for grid in x, in Q2, and data in (x, Q2).
659  bool doExtraPol;
660  double rescale;
661  double xGrid[100];
662  double Q2Grid[88];
663  double gluonGrid[100][88];
664  double singletGrid[100][88];
665  double charmGrid[100][88];
666 
667  // These inits do not overwrite PDF init (prevents Clang warnings).
668  using PDF::init;
669 
670  // Initialization of data array.
671  void init( int iFit, string pdfdataPath, Logger* loggerPtr);
672 
673  // Initialization through a stream.
674  void init( istream& is, Logger* loggerPtr);
675 
676  // Update PDF values.
677  void xfUpdate(int id, double x, double ) override;
678 
679 };
680 
681 //==========================================================================
682 
683 // A proton masked as a Pomeron for use within the Heavy Ion machinery
684 
685 class PomHISASD : public PDF {
686 
687 public:
688 
689  // Basic constructor
690  PomHISASD(int idBeamIn, PDFPtr ppdf, Settings & settings,
691  Logger* loggerPtrIn = 0) : PDF(idBeamIn), pPDFPtr(ppdf),
692  xPomNow(-1.0), hixpow(4.0), newfac(1.0) {
693  loggerPtr = loggerPtrIn;
694  hixpow = settings.parm("PDF:PomHixSupp");
695  if ( settings.mode("Angantyr:SASDmode") == 3 ) newfac =
696  log(settings.parm("Beams:eCM")/settings.parm("Diffraction:mMinPert"));
697  if ( settings.mode("Angantyr:SASDmode") == 4 ) newfac = 0.0;
698  }
699 
700  // Delete also the proton PDF
702 
703  // (re-)Set the x_pomeron value.
704  void xPom(double xpom = -1.0) override { xPomNow = xpom; }
705 
706 private:
707 
708  // The proton PDF.
709  PDFPtr pPDFPtr;
710 
711  // The momentum fraction if the corresponding pomeron.
712  double xPomNow;
713 
714  // The high-x suppression power.
715  double hixpow;
716 
717  // Special options.
718  double newfac;
719 
720  // Report possible errors.
721  Logger* loggerPtr;
722 
723  // Update PDF values.
724  void xfUpdate(int , double x, double Q2) override;
725 
726 };
727 
728 //==========================================================================
729 
730 // Gives electron (or muon, or tau) parton distribution.
731 
732 class Lepton : public PDF {
733 
734 public:
735 
736  // Constructor.
737  Lepton(int idBeamIn = 11) : PDF(idBeamIn), m2Lep(), Q2maxGamma(),
738  infoPtr(), rndmPtr() {}
739 
740  // Constructor with further info.
741  Lepton(int idBeamIn, double Q2maxGammaIn, Info* infoPtrIn)
742  : PDF(idBeamIn), m2Lep() { Q2maxGamma = Q2maxGammaIn;
743  infoPtr = infoPtrIn; rndmPtr = infoPtrIn->rndmPtr; }
744 
745  // Sample the Q2 value.
746  double sampleQ2gamma(double Q2min) override
747  { return Q2min * pow(Q2maxGamma / Q2min, rndmPtr->flat()); }
748 
749 private:
750 
751  // Constants: could only be changed in the code itself.
752  static const double ALPHAEM, ME, MMU, MTAU;
753 
754  // Update PDF values.
755  void xfUpdate(int id, double x, double Q2) override;
756 
757  // The squared lepton mass, set at initialization.
758  double m2Lep, Q2maxGamma;
759 
760  // Pointer to info, needed to get sqrt(s) to fix x_gamma limits.
761  Info* infoPtr;
762 
763  // Pointer to random number generator, needed to sample Q2.
764  Rndm* rndmPtr;
765 
766 };
767 
768 //==========================================================================
769 
770 // Gives electron (or other lepton) parton distribution when unresolved.
771 
772 class LeptonPoint : public PDF {
773 
774 public:
775 
776  // Constructor.
777  LeptonPoint(int idBeamIn = 11) : PDF(idBeamIn) {}
778 
779 private:
780 
781  // Update PDF values in trivial way.
782  void xfUpdate(int , double , double ) override {xlepton = 1; xgamma = 0.;}
783 
784 };
785 
786 //==========================================================================
787 
788 // Gives neutrino parton distribution when unresolved (only choice for now).
789 // Note that the extra factor of 2 - wrt. charged leptons as there is no need
790 // for spin averaging since neutrinos always lefthanded - is taken care in
791 // cross sections and not in the PDFs.
792 
793 class NeutrinoPoint : public PDF {
794 
795 public:
796 
797  // Constructor.
798  NeutrinoPoint(int idBeamIn = 12) : PDF(idBeamIn) {}
799 
800 private:
801 
802  // Update PDF values in trivial way.
803  void xfUpdate(int , double , double ) override {xlepton = 1; xgamma = 0.;}
804 
805 };
806 
807 //==========================================================================
808 
809 // Gives the CJKL leading order parton distribution function set
810 // in parametrized form for the real photons. Authors: F.Cornet, P.Jankowski,
811 // M.Krawczyk and A.Lorca, Phys. Rev. D68: 014010, 2003.
812 
813 class CJKL : public PDF {
814 
815 public:
816 
817  // Constructor. Needs the randon number generator to sample valence content.
818  CJKL(int idBeamIn = 22, Rndm* rndmPtrIn = 0 ) : PDF(idBeamIn) {
819  rndmPtr = rndmPtrIn; }
820 
821  // Functions to approximate pdfs for ISR.
822  double gammaPDFxDependence(int id, double) override;
823  double gammaPDFRefScale(int) override;
824 
825  // Set the valence content for photons.
826  int sampleGammaValFlavor(double Q2) override;
827 
828  // The total x-integrated PDFs. Relevant for MPIs with photon beams.
829  double xfIntegratedTotal(double Q2) override;
830 
831 private:
832 
833  // Parameters related to the fit.
834  static const double ALPHAEM, Q02, Q2MIN, Q2REF, LAMBDA, MC, MB;
835 
836  // Pointer to random number generator used for valence sampling.
837  Rndm *rndmPtr;
838 
839  // Update PDF values.
840  void xfUpdate(int , double x, double Q2) override;
841 
842  // Functions for updating the point-like part.
843  double pointlikeG(double x, double s);
844  double pointlikeU(double x, double s);
845  double pointlikeD(double x, double s);
846  double pointlikeC(double x, double s, double Q2);
847  double pointlikeB(double x, double s, double Q2);
848 
849  // Functions for updating the hadron-like part.
850  double hadronlikeG(double x, double s);
851  double hadronlikeSea(double x, double s);
852  double hadronlikeVal(double x, double s);
853  double hadronlikeC(double x, double s, double Q2);
854  double hadronlikeB(double x, double s, double Q2);
855 
856 };
857 
858 //==========================================================================
859 
860 // Convolution with photon flux from leptons and photon PDFs.
861 // Photon flux from equivalent photon approximation (EPA).
862 // Contains a pointer to a photon PDF set and samples the
863 // convolution integral event-by-event basis.
864 // Includes also a overestimate for the PDF set in order to set up
865 // the phase-space sampling correctly.
866 
867 class Lepton2gamma : public PDF {
868 
869 public:
870 
871  // Constructor.
872  Lepton2gamma(int idBeamIn, double m2leptonIn, double Q2maxGamma,
873  PDFPtr gammaPDFPtrIn, Info* infoPtrIn)
874  : PDF(idBeamIn), m2lepton(m2leptonIn), Q2max(Q2maxGamma), xGm(),
875  sampleXgamma(true), gammaPDFPtr(gammaPDFPtrIn),rndmPtr(infoPtrIn->rndmPtr),
876  infoPtr(infoPtrIn) { hasGammaInLepton = true; }
877 
878  // Overload the member function definitions where relevant.
879  void xfUpdate(int id, double x, double Q2) override;
880  double xGamma() override { return xGm; }
881  double xfMax(int id, double x, double Q2) override;
882  double xfSame(int id, double x, double Q2) override;
883 
884  // Sample the Q2 value.
885  double sampleQ2gamma(double Q2min) override
886  { return Q2min * pow(Q2max / Q2min, rndmPtr->flat()); }
887 
888 private:
889 
890  // Parameters for convolution.
891  static const double ALPHAEM, Q2MIN;
892  double m2lepton, Q2max, xGm;
893 
894  // Sample new value for x_gamma.
895  bool sampleXgamma;
896 
897  // Photon PDFs which the photon flux is convoluted with.
898  PDFPtr gammaPDFPtr;
899 
900  // Pointer to random number generator used for sampling x_gamma.
901  Rndm* rndmPtr;
902 
903  // Pointer to info, needed to get sqrt(s) to fix x_gamma limits.
904  Info* infoPtr;
905 
906 };
907 
908 //==========================================================================
909 
910 // Gives photon parton distribution when unresolved.
911 
912 class GammaPoint : public PDF {
913 
914 public:
915 
916  // Constructor.
917  GammaPoint(int idBeamIn = 22) : PDF(idBeamIn) {}
918 
919 private:
920 
921  // Update PDF values in trivial way.
922  void xfUpdate(int , double , double ) override { xgamma = 1.;}
923 
924 };
925 
926 //==========================================================================
927 
928 // Unresolved proton: equivalent photon spectrum according
929 // to the approximation by Drees and Zeppenfeld,
930 // Phys.Rev. D39 (1989) 2536.
931 
932 class Proton2gammaDZ : public PDF {
933 
934 public:
935 
936  // Constructor.
937  Proton2gammaDZ(int idBeamIn = 2212) : PDF(idBeamIn) {}
938 
939 private:
940 
941  // Stored parameters.
942  static const double ALPHAEM, Q20;
943 
944  // Update PDF values.
945  void xfUpdate(int , double x, double Q2) override;
946  double fluxQ2dependence(double Q2) override;
947 
948 };
949 
950 //==========================================================================
951 
952 // Unresolved nucleus: equivalent photon approximation
953 // for impact parameter integrated flux according to standard
954 // form introduced in J.D. Jackson, Classical Electrodynamics,
955 // 2nd edition, John Wiley & Sons (1975).
956 
957 class Nucleus2gamma : public PDF {
958 
959 public:
960 
961  // Constructor.
962  Nucleus2gamma(int idBeamIn, double bMinIn, double mNucleonIn) :
963  PDF(idBeamIn), a(), z(), bMin(bMinIn), mNucleon(mNucleonIn)
964  { initNucleus(idBeamIn); }
965 
966 private:
967 
968  // Stored constant parameters.
969  static const double ALPHAEM;
970 
971  // Initialize flux parameters.
972  void initNucleus(int idBeamIn);
973 
974  // Update PDF values.
975  void xfUpdate(int , double x, double Q2) override;
976 
977  // Mass number and electric charge.
978  int a, z;
979 
980  // Minimum impact parameter for integration and per-nucleon mass.
981  double bMin, mNucleon;
982 
983 };
984 
985 //==========================================================================
986 
987 // Equivalent photon approximation for sampling with external photon flux.
988 
989 class EPAexternal : public PDF {
990 
991 public:
992 
993  // Constructor.
994  EPAexternal(int idBeamIn, double m2In, PDFPtr gammaFluxPtrIn,
995  PDFPtr gammaPDFPtrIn, Info* infoPtrIn, Logger* loggerPtrIn = 0)
996  : PDF(idBeamIn), m2(m2In), Q2max(), Q2min(), xMax(), xMin(), xHadr(),
997  norm(), xPow(), xCut(), norm1(), norm2(), integral1(), integral2(),
998  bmhbarc(), approxMode(0), isLHA(false), gammaFluxPtr(gammaFluxPtrIn),
999  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
virtual int sampleGammaValFlavor(double)
Sample the valence content for photons.
Definition: PartonDistributions.h:122
virtual double xGamma()
Return the sampled value for x_gamma.
Definition: PartonDistributions.h:128
Lepton2gamma(int idBeamIn, double m2leptonIn, double Q2maxGamma, PDFPtr gammaPDFPtrIn, Info *infoPtrIn)
Constructor.
Definition: PartonDistributions.h:872
MSTWpdf(int idBeamIn, istream &is, Logger *loggerPtr=0)
Constructor with a stream.
Definition: PartonDistributions.h:343
void printErr(string loc, string errMsg, Logger *loggerPtr=nullptr)
Small routine for error printout, depending on loggerPtr existing or not.
Definition: PartonDistributions.h:198
PomHISASD(int idBeamIn, PDFPtr ppdf, Settings &settings, Logger *loggerPtrIn=0)
Basic constructor.
Definition: PartonDistributions.h:690
Definition: PartonDistributions.h:1110
void resetValenceContent()
Set valence content.
Definition: PartonDistributions.cc:24
Gives electron (or other lepton) parton distribution when unresolved.
Definition: PartonDistributions.h:772
CJKL(int idBeamIn=22, Rndm *rndmPtrIn=0)
Constructor. Needs the randon number generator to sample valence content.
Definition: PartonDistributions.h:818
PDF(int idBeamIn=2212)
Constructor.
Definition: PartonDistributions.h:54
Rndm * rndmPtr
Pointer to the random number generator.
Definition: Info.h:89
A proton masked as a Pomeron for use within the Heavy Ion machinery.
Definition: PartonDistributions.h:685
virtual bool insideBounds(double, double)
Check whether x and Q2 values fall inside the fit bounds (LHAPDF6 only).
Definition: PartonDistributions.h:90
CTEQ5L(int idBeamIn=2212)
Constructor.
Definition: PartonDistributions.h:310
Definition: PartonDistributions.h:957
CTEQ6pdf(int idBeamIn=2212, int iFitIn=1, double rescaleIn=1., string pdfdataPath="../share/Pythia8/pdfdata/", Logger *loggerPtr=0)
Constructor.
Definition: PartonDistributions.h:408
Definition: Info.h:45
bool sSymmetricSave
Whether to treat flavoured PDFs as symmetric, for efficiency.
Definition: PartonDistributions.h:191
void xPom(double xpom=-1.0) override
(re-)Set the x_pomeron value.
Definition: PartonDistributions.h:704
virtual void setExtrapolate(bool)
Allow extrapolation beyond boundaries. This is optional.
Definition: PartonDistributions.h:80
Definition: PartonDistributions.h:793
Nucleus2gamma(int idBeamIn, double bMinIn, double mNucleonIn)
Constructor.
Definition: PartonDistributions.h:962
Definition: PartonDistributions.h:331
LeptonPoint(int idBeamIn=11)
Constructor.
Definition: PartonDistributions.h:777
void setVMDscale(double vmdScaleIn=1.) override
Allow for new rescaling factor of the PDF for VMD beams.
Definition: PartonDistributions.h:506
virtual double xfMax(int id, double x, double Q2)
Normal PDFs unless gamma inside lepton -> an overestimate for sampling.
Definition: PartonDistributions.h:148
LHAGrid1(int idBeamIn=2212, string pdfWord="void", string xmlPath="../share/Pythia8/xmldoc/", Logger *loggerPtr=0)
Constructor.
Definition: PartonDistributions.h:223
~PomHISASD()
Delete also the proton PDF.
Definition: PartonDistributions.h:701
Gives photon parton distribution when unresolved.
Definition: PartonDistributions.h:912
double xf(int id, double x, double Q2)
Read out parton density.
Definition: PartonDistributions.cc:122
PomH1FitAB(int idBeamIn=990, int iFit=1, double rescaleIn=1., string pdfdataPath="../share/Pythia8/pdfdata/", Logger *loggerPtr=0)
Constructor.
Definition: PartonDistributions.h:591
double xfVal(int id, double x, double Q2)
Read out valence and sea part of parton densities.
Definition: PartonDistributions.cc:222
double sampleQ2gamma(double Q2min) override
Sample the Q2 value.
Definition: PartonDistributions.h:885
virtual double xfIntegratedTotal(double)
The total x-integrated PDFs. Relevant for MPIs with photon beams.
Definition: PartonDistributions.h:125
void setVMDscale(double vmdScaleIn=1.) override
Allow for new rescaling factor of the PDF for VMD beams.
Definition: PartonDistributions.h:532
virtual void xfUpdate(int id, double x, double Q2)=0
Update parton densities.
GRV94L(int idBeamIn=2212)
Constructor.
Definition: PartonDistributions.h:283
LHAGrid1(int idBeamIn, istream &is, Logger *loggerPtr=0)
Constructor with a stream.
Definition: PartonDistributions.h:230
GRSpiL(int idBeamIn=211, double vmdScaleIn=1.)
Constructor.
Definition: PartonDistributions.h:528
void setErrorSet(int iSetIn)
Use other than central set to study uncertainties.
Definition: PartonDistributions.h:1186
virtual double xfSame(int id, double x, double Q2)
Normal PDFs unless gamma inside lepton -> Do not sample x_gamma.
Definition: PartonDistributions.h:151
nPDF(int idBeamIn=2212, PDFPtr protonPDFPtrIn=0)
Constructor.
Definition: PartonDistributions.h:1065
virtual bool init(int, string, int, Logger *)
Definition: PartonDistributions.h:64
double xfRaw(int id) const
Get the raw stored value for the quark variable corresponding to the id.
Definition: PartonDistributions.cc:101
double xGamma() override
Return the sampled value for x_gamma.
Definition: PartonDistributions.h:880
virtual double xfFlux(int, double, double)
Return accurate and approximated photon fluxes and PDFs.
Definition: PartonDistributions.h:134
Gives electron (or muon, or tau) parton distribution.
Definition: PartonDistributions.h:732
bool sSymmetric() const
Return if s/sbar, c/cbar, and b/bbar PDFs are symmetric.
Definition: PartonDistributions.h:157
Definition: Logger.h:23
bool hasGammaInLepton
True if a photon beam inside a lepton beam, otherwise set false.
Definition: PartonDistributions.h:188
Definition: PartonDistributions.h:305
void setExtrapolate(bool doExtraPolIn) override
Allow extrapolation beyond boundaries. This is optional.
Definition: PartonDistributions.h:427
void setErrorSet(int iSetIn)
Use other than central set to study uncertainties.
Definition: PartonDistributions.h:1141
void sSymmetric(bool sSymmetricIn)
Set s/sbar, c/cbar, and b/bbar PDFs symmetric.
Definition: PartonDistributions.h:162
virtual double mQuarkPDF(int)
Return quark masses used in the PDF fit (LHAPDF6 only).
Definition: PartonDistributions.h:96
void rUpdate(int, double, double) override
Only the Isospin effect so no need to do anything here.
Definition: PartonDistributions.h:1119
Definition: Basics.h:385
Definition: PartonDistributions.h:218
Error envelope from PDF uncertainty.
Definition: PartonDistributions.h:102
int beamType
Definition: PartonDistributions.h:185
Equivalent photon approximation for sampling with external photon flux.
Definition: PartonDistributions.h:989
Definition: PartonDistributions.h:278
virtual void setBeamID(int idBeamIn)
Switch to new beam particle identities; for similar hadrons only.
Definition: PartonDistributions.h:70
MSTWpdf(int idBeamIn=2212, int iFitIn=1, string pdfdataPath="../share/Pythia8/pdfdata/", Logger *loggerPtr=0)
Constructor.
Definition: PartonDistributions.h:336
double sampleQ2gamma(double Q2min) override
Sample the Q2 value.
Definition: PartonDistributions.h:746
bool hasApproxGammaFlux() override
This derived class use approximated flux for sampling.
Definition: PartonDistributions.h:1013
~LHAGrid1()
Destructor.
Definition: PartonDistributions.h:236
void setExtrapolate(bool doExtraPolIn) override
Allow extrapolation beyond boundaries. This is optional.
Definition: PartonDistributions.h:243
double xfSea(int id, double x, double Q2)
Only sea part of parton densities.
Definition: PartonDistributions.cc:301
double xu
Stored quantities.
Definition: PartonDistributions.h:172
GRVpiL(int idBeamIn=211, double vmdScaleIn=1.)
Constructor.
Definition: PartonDistributions.h:502
A derived class for nuclear PDFs. Needs a pointer for (free) proton PDFs.
Definition: PartonDistributions.h:1060
EPAexternal(int idBeamIn, double m2In, PDFPtr gammaFluxPtrIn, PDFPtr gammaPDFPtrIn, Info *infoPtrIn, Logger *loggerPtrIn=0)
Constructor.
Definition: PartonDistributions.h:994
PomH1Jets(int idBeamIn, double rescaleIn, istream &is, Logger *loggerPtr=0)
Constructor with a stream.
Definition: PartonDistributions.h:648
bool isSetup()
Confirm that PDF has been set up (important for LHAPDF and H1 Pomeron).
Definition: PartonDistributions.h:67
Proton2gammaDZ(int idBeamIn=2212)
Constructor.
Definition: PartonDistributions.h:937
NeutrinoPoint(int idBeamIn=12)
Constructor.
Definition: PartonDistributions.h:798
void setExtrapolate(bool doExtraPolIn) override
Allow extrapolation beyond boundaries. This is optional.
Definition: PartonDistributions.h:654
virtual void xPom(double=-1.0)
Keep track of pomeron momentum fraction.
Definition: PartonDistributions.h:131
GammaPoint(int idBeamIn=22)
Constructor.
Definition: PartonDistributions.h:917
int idBeam
Store relevant quantities.
Definition: PartonDistributions.h:169
ProtonPoint(int idBeamIn=2212, Logger *loggerPtrIn=0)
Constructor.
Definition: PartonDistributions.h:473
Nuclear modifications from EPS09 fit.
Definition: PartonDistributions.h:1126
virtual double gammaPDFRefScale(int)
Provide the reference scale for logarithmic Q^2 evolution for photons.
Definition: PartonDistributions.h:119
void setMode(double zaIn)
Definition: PartonDistributions.h:1084
void setExtrapolate(bool doExtraPolIn) override
Allow extrapolation beyond boundaries. This is optional.
Definition: PartonDistributions.h:604
Definition: PartonDistributions.h:932
Gives generic Q2-independent Pomeron PDF.
Definition: PartonDistributions.h:548
virtual double gammaPDFxDependence(int, double)
Approximate photon PDFs by decoupling the scale and x-dependence.
Definition: PartonDistributions.h:116
int getA()
Return the number of protons and nucleons.
Definition: PartonDistributions.h:1079
double getXmin() override
Kinematics.
Definition: PartonDistributions.h:1016
double m2(const Vec4 &v1)
The squared invariant mass of one or more four-vectors.
Definition: Basics.cc:586
CTEQ6pdf(int idBeamIn, istream &is, bool isPdsGrid=false, Logger *loggerPtr=0)
Constructor with a stream.
Definition: PartonDistributions.h:418
virtual void calcPDFEnvelope(int, double, double, int)
Calculate PDF envelope.
Definition: PartonDistributions.h:110
Definition: PartonDistributions.h:523
Base class for parton distribution functions.
Definition: PartonDistributions.h:49
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
Definition: PartonDistributions.h:468
virtual void setVMDscale(double=1.)
Allow for new scaling factor for VMD PDFs.
Definition: PartonDistributions.h:154
bool isValence(int id) const
Check whether the specified id is a valence quark.
Definition: PartonDistributions.h:207
PomH1FitAB(int idBeamIn, double rescaleIn, istream &is, Logger *loggerPtr=0)
Constructor with a stream.
Definition: PartonDistributions.h:598
Nuclear modifications from EPPS16 fit.
Definition: PartonDistributions.h:1170
Definition: PartonDistributions.h:636
virtual int nMembers()
Return number of members of this PDF family (LHAPDF6 only).
Definition: PartonDistributions.h:99
PomH1Jets(int idBeamIn=990, int iFit=1, double rescaleIn=1., string pdfdataPath="../share/Pythia8/pdfdata/", Logger *loggerPtr=0)
Constructor.
Definition: PartonDistributions.h:641
Lepton(int idBeamIn=11)
Constructor.
Definition: PartonDistributions.h:737
virtual ~PDF()
Virtual destructor.
Definition: PartonDistributions.h:59
Definition: PartonDistributions.h:813
virtual double alphaS(double)
Access the running alpha_s of a PDF set (LHAPDF6 only).
Definition: PartonDistributions.h:93
Definition: PartonDistributions.h:586
Isospin(int idBeamIn=2212, PDFPtr protonPDFPtrIn=0)
Constructor.
Definition: PartonDistributions.h:1115
virtual double getXmin()
Return the kinematical limits and sample Q2 and x.
Definition: PartonDistributions.h:141
double ruv
Definition: PartonDistributions.h:1091
EPPS16(int idBeamIn=2212, int iSetIn=1, string pdfdataPath="../share/Pythia8/pdfdata/", PDFPtr protonPDFPtrIn=0, Logger *loggerPtrIn=0)
Constructor.
Definition: PartonDistributions.h:1175
Definition: PartonDistributions.h:497
Definition: PartonDistributions.h:403
PomFix(int idBeamIn=990, double PomGluonAIn=0., double PomGluonBIn=0., double PomQuarkAIn=0., double PomQuarkBIn=0., double PomQuarkFracIn=0., double PomStrangeSuppIn=0.)
Constructor.
Definition: PartonDistributions.h:553
Definition: Settings.h:195
Definition: PartonDistributions.h:867
EPS09(int idBeamIn=2212, int iOrderIn=1, int iSetIn=1, string pdfdataPath="../share/Pythia8/pdfdata/", PDFPtr protonPDFPtrIn=0, Logger *loggerPtrIn=0)
Constructor.
Definition: PartonDistributions.h:1131
Lepton(int idBeamIn, double Q2maxGammaIn, Info *infoPtrIn)
Constructor with further info.
Definition: PartonDistributions.h:741