1 // PartonDistributions.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2008 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
6 // Header file for parton densities.
8 // GRV94L: derived class for the GRV 94L parton densities.
9 // CTEQ5L: derived class for the CTEQ 5L parton densities.
10 // LHAPDFinterface: derived class for interface to the LHAPDF library.
11 // Lepton: derived class for parton densities inside a lepton.
12 // LeptonPoint: derived class for unresolved lepton (mainly dummy).
14 #ifndef Pythia8_PartonDistributions_H
15 #define Pythia8_PartonDistributions_H
19 #include "ParticleData.h"
20 #include "PythiaStdlib.h"
24 //**************************************************************************
26 // Base class for parton distribution functions.
33 PDF(int idBeamIn = 2212) {idBeam = idBeamIn; idSav = 9; xSav = -1.;
34 Q2Sav = -1.; isSet = true; isInit = false;}
39 // Confirm that PDF has been set up (important for LHAPDF).
40 bool isSetup() {return isSet;}
42 // Allow extrapolation beyond boundaries. This is optional.
43 virtual void setExtrapolate(bool) {}
45 // Read out parton density
46 double xf(int id, double x, double Q2);
48 // Read out valence and sea part of parton densities.
49 double xfVal(int id, double x, double Q2);
50 double xfSea(int id, double x, double Q2);
54 // Store relevant quantities.
57 double xu, xd, xubar, xdbar, xs, xc, xb, xg, xlepton, xgamma,
58 xuVal, xuSea, xdVal, xdSea;
59 bool isSet, isInit, hasLimits;
61 // Update parton densities.
62 virtual void xfUpdate(int id, double x, double Q2) = 0;
66 //**************************************************************************
68 // Gives the GRV 94 L (leading order) parton distribution function set
69 // in parametrized form. Authors: M. Glueck, E. Reya and A. Vogt.
71 class GRV94L : public PDF {
76 GRV94L(int idBeamIn = 2212) : PDF(idBeamIn) {}
81 void xfUpdate(int id, double x, double Q2);
83 // Auxiliary routines used during the updating.
84 double grvv (double x, double n, double ak, double bk, double a,
85 double b, double c, double d);
86 double grvw (double x, double s, double al, double be, double ak,
87 double bk, double a, double b, double c, double d, double e, double es);
88 double grvs (double x, double s, double sth, double al, double be,
89 double ak, double ag, double b, double d, double e, double es);
93 //**************************************************************************
95 // Gives the GRV 94 L (leading order) parton distribution function set
96 // in parametrized form. Authors: M. Glueck, E. Reya and A. Vogt.
98 class CTEQ5L : public PDF {
103 CTEQ5L(int idBeamIn = 2212) : PDF(idBeamIn) {}
107 // Update PDF values.
108 void xfUpdate(int id, double x, double Q2);
112 //**************************************************************************
114 // Provide interface to the LHAPDF library of parton densities.
116 class LHAPDF : public PDF {
121 LHAPDF(int idBeamIn, string setName, int member, int nSetIn = 1,
122 Info* infoPtr = 0) : PDF(idBeamIn), nSet(nSetIn)
123 {init( setName, member, infoPtr);}
125 // Allow extrapolation beyond boundaries. This is optional.
126 void setExtrapolate(bool extrapol);
130 // Initialization of PDF set.
131 void init(string setName, int member, Info* infoPtr);
133 // Update all PDF values.
134 void xfUpdate(int , double x, double Q2);
136 // Current set and pdf values.
140 // Keep track of latest initialized PDF, so does not have to repeat.
141 static string latestSetName;
142 static int latestMember, latestNSet;
146 //**************************************************************************
148 // Gives electron (or muon, or tau) parton distribution.
150 class Lepton : public PDF {
155 Lepton(int idBeamIn = 11) : PDF(idBeamIn) {}
159 // Constants: could only be changed in the code itself.
160 static const double ALPHAEM;
162 // Update PDF values.
163 void xfUpdate(int id, double x, double Q2);
165 // The lepton mass, set at initialization.
170 //**************************************************************************
172 // Gives electron (or other lepton) parton distribution when unresolved.
174 class LeptonPoint : public PDF {
179 LeptonPoint(int idBeamIn = 11) : PDF(idBeamIn) {}
183 // Update PDF values in trivial way.
184 void xfUpdate(int , double , double ) {xlepton = 1; xgamma = 0.;}
188 } // end namespace Pythia8
190 //**************************************************************************
192 #endif // Pythia8_PartonDistributions_H