source: HiSusy/trunk/Pythia8/pythia8170/include/PartonDistributions.h @ 1

Last change on this file since 1 was 1, checked in by zerwas, 11 years ago

first import of structure, PYTHIA8 and DELPHES

File size: 13.2 KB
Line 
1// PartonDistributions.h is a part of the PYTHIA event generator.
2// Copyright (C) 2012 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.
5
6// Header file for parton densities.
7// PDF: base class.
8// LHAPDF: derived class for interface to the LHAPDF library.
9// GRV94L: derived class for the GRV 94L parton densities.
10// CTEQ5L: derived class for the CTEQ 5L parton densities.
11// MSTWpdf: derived class for MRST LO*, LO**, MSTW 2008 LO, NLO.
12// CTEQ6pdf: derived class for CTEQ 6L, 6L1, 66, CT09 MC1, MC2, (MCS?).
13// ProtonPoint: unresolved proton with equivalent photon spectrum. 
14// GRVpiL: derived class for the GRV LO pion parton densities.
15// PomFix: derived class for Q2-independent Pomeron parton densities.
16// PomH1FitAB: derived class for the H1 2006 Fit A and Fit B Pomeron PDFs.
17// PomH1Jets: derived class for the H1 2007 Jets Pomeron PDFs.
18// Lepton: derived class for parton densities inside a lepton.
19// LeptonPoint: derived class for unresolved lepton (mainly dummy).
20
21#ifndef Pythia8_PartonDistributions_H
22#define Pythia8_PartonDistributions_H
23
24#include "Basics.h"
25#include "Info.h"
26#include "ParticleData.h"
27#include "PythiaStdlib.h"
28
29namespace Pythia8 {
30
31//==========================================================================
32
33// Base class for parton distribution functions.
34
35class PDF {
36
37public:
38
39  // Constructor.
40  PDF(int idBeamIn = 2212) {idBeam = idBeamIn; idBeamAbs = abs(idBeam);
41    setValenceContent(); idSav = 9; xSav = -1.; Q2Sav = -1.;
42    xu = 0.; xd = 0.; xs = 0.; xubar = 0.; xdbar = 0.; xsbar = 0.; xc = 0.; 
43    xb = 0.; xg = 0.; xlepton = 0.; xgamma = 0.; xuVal = 0.; xuSea = 0.;
44    xdVal = 0.; xdSea = 0.; isSet = true; isInit = false;}
45
46  // Destructor.
47  virtual ~PDF() {}
48
49  // Confirm that PDF has been set up (important for LHAPDF and H1 Pomeron).
50  bool isSetup() {return isSet;}
51
52  // Dynamic choice of meson valence flavours for pi0, K0S, K0L, Pomeron.
53  void newValenceContent(int idVal1In, int idVal2In) {
54    idVal1 = idVal1In; idVal2 = idVal2In;}
55
56  // Allow extrapolation beyond boundaries. This is optional.
57  virtual void setExtrapolate(bool) {}
58
59  // Read out parton density
60  double xf(int id, double x, double Q2);
61
62  // Read out valence and sea part of parton densities.
63  double xfVal(int id, double x, double Q2);
64  double xfSea(int id, double x, double Q2);
65 
66protected:
67
68  // Store relevant quantities.
69  int    idBeam, idBeamAbs, idSav, idVal1, idVal2;
70  double xSav, Q2Sav;
71  double xu, xd, xs, xubar, xdbar, xsbar, xc, xb, xg, xlepton, xgamma,
72         xuVal, xuSea, xdVal, xdSea;
73  bool   isSet, isInit; 
74
75  // Resolve valence content for assumed meson. Possibly modified later.
76  void setValenceContent();
77
78  // Update parton densities.
79  virtual void xfUpdate(int id, double x, double Q2) = 0; 
80
81};
82 
83//==========================================================================
84
85// Provide interface to the LHAPDF library of parton densities.
86
87class LHAPDF : public PDF {
88
89public:
90
91  // Constructor.
92  LHAPDF(int idBeamIn, string setName, int member,  int nSetIn = 1, 
93    Info* infoPtr = 0) : PDF(idBeamIn), nSet(nSetIn) 
94    {init( setName, member, infoPtr);} 
95
96  // Allow extrapolation beyond boundaries. This is optional.
97  void setExtrapolate(bool extrapol); 
98
99private:
100
101  // Initialization of PDF set.
102  void init(string setName, int member, Info* infoPtr);
103
104  // Update all PDF values.
105  void xfUpdate(int , double x, double Q2);
106
107  // Current set and pdf values.
108  int    nSet;
109  double xfArray[13];
110  double xPhoton;
111
112  // Keep track of latest initialized PDF, so does not have to repeat.
113  static string latestSetName;
114  static int    latestMember, latestNSet;   
115
116};
117 
118//==========================================================================
119
120// Gives the GRV 94L (leading order) parton distribution function set
121// in parametrized form. Authors: M. Glueck, E. Reya and A. Vogt.
122
123class GRV94L : public PDF {
124
125public:
126
127  // Constructor.
128  GRV94L(int idBeamIn = 2212) : PDF(idBeamIn) {}
129
130private:
131
132  // Update PDF values.
133  void xfUpdate(int , double x, double Q2);
134
135  // Auxiliary routines used during the updating.
136  double grvv (double x, double n, double ak, double bk, double a, 
137    double b, double c, double d);
138  double grvw (double x, double s, double al, double be, double ak, 
139    double bk, double a, double b, double c, double d, double e, double es);
140  double grvs (double x, double s, double sth, double al, double be, 
141    double ak, double ag, double b, double d, double e, double es);
142
143};
144 
145//==========================================================================
146
147// Gives the CTEQ 5L (leading order) parton distribution function set
148// in parametrized form. Parametrization by J. Pumplin. Authors: CTEQ.
149
150class CTEQ5L : public PDF {
151
152public:
153
154  // Constructor.
155  CTEQ5L(int idBeamIn = 2212) : PDF(idBeamIn) {}
156
157private:
158
159  // Update PDF values.
160  void xfUpdate(int , double x, double Q2);
161
162};
163 
164//==========================================================================
165
166// The MSTWpdf class.
167// MRST LO*(*) and MSTW 2008 PDF's, specifically the LO one.
168// Original C++ version by Jeppe Andersen.
169// Modified by Graeme Watt <watt(at)hep.ucl.ac.uk>.
170// Sets available:
171// iFit = 1 : MRST LO*  (2007).
172// iFit = 2 : MRST LO** (2008).
173// iFit = 3 : MSTW 2008 LO, central member.
174// iFit = 4 : MSTW 2008 NLO, central member. (Warning!)
175
176class MSTWpdf : public PDF {
177
178public:
179
180  // Constructor.
181  MSTWpdf(int idBeamIn = 2212, int iFitIn = 1, string xmlPath = "../xmldoc/", 
182    Info* infoPtr = 0) : PDF(idBeamIn) {init( iFitIn,  xmlPath, infoPtr);}
183
184private:
185
186  // Constants: could only be changed in the code itself.
187  static const int    np, nx, nq, nqc0, nqb0;
188  static const double xmin, xmax, qsqmin, qsqmax, xxInit[65], qqInit[49];
189
190  // Data read in from grid file or set at initialization.
191  int    iFit, alphaSorder, alphaSnfmax;
192  double mCharm, mBottom, alphaSQ0, alphaSMZ, distance, tolerance, 
193         xx[65], qq[49], c[13][64][48][5][5];
194
195  // Initialization of data array.
196  void init( int iFitIn, string xmlPath, Info* infoPtr);
197
198  // Update PDF values.
199  void xfUpdate(int , double x, double Q2);
200
201  // Evaluate PDF of one flavour species.
202  double parton(int flavour,double x,double q);
203  double parton_interpolate(int flavour,double xxx,double qqq);
204  double parton_extrapolate(int flavour,double xxx,double qqq);
205
206  // Auxiliary routines for evaluation.
207  int locate(double xx[],int n,double x);
208  double polderivative1(double x1, double x2, double x3, double y1,
209    double y2, double y3);
210  double polderivative2(double x1, double x2, double x3, double y1, 
211    double y2, double y3);
212  double polderivative3(double x1, double x2, double x3, double y1, 
213    double y2, double y3);
214
215};
216 
217//==========================================================================
218
219// The CTEQ6pdf class.
220// Sets available:
221// iFit = 1 : CTEQ6L
222// iFit = 2 : CTEQ6L1
223// iFit = 3 : CTEQ66.00 (NLO, central member)
224// iFit = 4 : CT09MC1
225// iFit = 5 : CT09MC2
226// iFit = 6 : CT09MCS (not yet implemented)
227
228class CTEQ6pdf : public PDF {
229
230public:
231
232  // Constructor.
233  CTEQ6pdf(int idBeamIn = 2212, int iFitIn = 1, string xmlPath = "../xmldoc/", 
234    Info* infoPtr = 0) : PDF(idBeamIn) {init( iFitIn, xmlPath, infoPtr);}
235
236private:
237
238  // Constants: could only be changed in the code itself.
239  static const double EPSILON, XPOWER;
240
241  // Data read in from grid file or set at initialization.
242  int    iFit, order, nQuark, nfMx, mxVal, nX, nT, nG,
243         iGridX, iGridQ, iGridLX, iGridLQ;
244  double lambda, mQ[7], qIni, qMax, tv[26], xMin, xv[202], upd[57773],
245         xvpow[202], xMinEps, xMaxEps, qMinEps, qMaxEps, fVec[5],
246         tConst[9], xConst[9], xLast, qLast;
247
248  // Initialization of data array.
249  void init( int iFitIn, string xmlPath, Info* infoPtr);
250
251  // Update PDF values.
252  void xfUpdate(int id, double x, double Q2);
253
254  // Evaluate PDF of one flavour species.
255  double parton6(int iParton, double x, double q);
256
257  // Interpolation in grid.
258  double polint4F(double xgrid[], double fgrid[], double xin);
259
260};
261
262//==========================================================================
263
264// SA Unresolved proton: equivalent photon spectrum from
265// V.M. Budnev, I.F. Ginzburg, G.V. Meledin and V.G. Serbo,
266// Phys. Rept. 15 (1974/1975) 181.
267
268class ProtonPoint : public PDF {
269
270public:
271
272  // Constructor.
273  ProtonPoint(int idBeamIn = 2212, Info* infoPtrIn = 0) : 
274              PDF(idBeamIn), m_infoPtr(infoPtrIn) {}
275
276private:
277
278  // Stored value for PDF choice.
279  static const double ALPHAEM, Q2MAX, Q20, A, B, C;
280
281  // Update PDF values.
282  void xfUpdate(int , double x, double Q2);
283
284  // phi function from Q2 integration.
285  double phiFunc(double x, double Q);
286
287  // Info and errors
288  Info* m_infoPtr;
289
290};
291 
292//==========================================================================
293
294// Gives the GRV 1992 pi+ (leading order) parton distribution function set
295// in parametrized form. Authors: Glueck, Reya and Vogt.
296
297class GRVpiL : public PDF {
298
299public:
300
301  // Constructor.
302  GRVpiL(int idBeamIn = 221) : PDF(idBeamIn) {}
303
304private:
305
306  // Update PDF values.
307  void xfUpdate(int , double x, double Q2);
308
309};
310
311//==========================================================================
312
313// Gives generic Q2-independent Pomeron PDF.
314
315class PomFix : public PDF {
316
317public:
318
319  // Constructor.
320  PomFix(int idBeamIn = 990, double PomGluonAIn = 0., 
321    double PomGluonBIn = 0., double PomQuarkAIn = 0., 
322    double PomQuarkBIn = 0., double PomQuarkFracIn = 0., 
323    double PomStrangeSuppIn = 0.) : PDF(idBeamIn), 
324    PomGluonA(PomGluonAIn), PomGluonB(PomGluonBIn), 
325    PomQuarkA(PomQuarkAIn), PomQuarkB(PomQuarkBIn), 
326    PomQuarkFrac(PomQuarkFracIn), PomStrangeSupp(PomStrangeSuppIn) 
327    {init();}
328
329private:
330
331  // Stored value for PDF choice.
332  double PomGluonA, PomGluonB, PomQuarkA, PomQuarkB, PomQuarkFrac, 
333         PomStrangeSupp, normGluon, normQuark;
334
335  // Initialization of some constants.
336  void init();
337
338  // Update PDF values.
339  void xfUpdate(int , double x, double);
340
341};
342 
343//==========================================================================
344
345// The H1 2006 Fit A and Fit B Pomeron parametrization.
346// H1 Collaboration, A. Aktas et al., "Measurement and QCD Analysis of
347// the Diffractive Deep-Inelastic Scattering Cross Section at HERA",
348// DESY-06-049, Eur. Phys. J. C48 (2006) 715. e-Print: hep-ex/0606004.
349
350class PomH1FitAB : public PDF {
351
352public:
353
354  // Constructor.
355 PomH1FitAB(int idBeamIn = 990, int iFit = 1, double rescaleIn = 1.,
356   string xmlPath = "../xmldoc/", Info* infoPtr = 0) : PDF(idBeamIn) 
357   {rescale = rescaleIn; init( iFit, xmlPath, infoPtr);}
358
359private:
360
361  // Limits for grid in x, in Q2, and data in (x, Q2).
362  int    nx, nQ2;
363  double rescale, xlow, xupp, dx, Q2low, Q2upp, dQ2;
364  double gluonGrid[100][30];
365  double quarkGrid[100][30];
366
367  // Initialization of data array.
368  void init( int iFit, string xmlPath, Info* infoPtr);
369
370  // Update PDF values.
371  void xfUpdate(int , double x, double );
372
373};
374 
375//==========================================================================
376
377// The H1 2007 Jets Pomeron parametrization..
378// H1 Collaboration, A. Aktas et al., "Dijet Cross Sections and Parton     
379// Densities in Diffractive DIS at HERA", DESY-07-115, Aug 2007. 33pp.   
380// Published in JHEP 0710:042,2007. e-Print: arXiv:0708.3217 [hep-ex] 
381
382class PomH1Jets : public PDF {
383
384public:
385
386  // Constructor.
387  PomH1Jets(int idBeamIn = 990,  double rescaleIn = 1., 
388   string xmlPath = "../xmldoc/", Info* infoPtr = 0) : PDF(idBeamIn) 
389   {rescale = rescaleIn; init( xmlPath, infoPtr);}
390
391private:
392
393  // Arrays for grid in x, in Q2, and data in (x, Q2).
394  double rescale;
395  double xGrid[100];
396  double Q2Grid[88];
397  double gluonGrid[100][88];
398  double singletGrid[100][88];
399  double charmGrid[100][88];
400
401  // Initialization of data array.
402  void init( string xmlPath, Info* infoPtr);
403
404  // Update PDF values.
405  void xfUpdate(int id, double x, double );
406
407};
408 
409//==========================================================================
410
411// Gives electron (or muon, or tau) parton distribution.
412 
413class Lepton : public PDF {
414
415public:
416
417  // Constructor.
418  Lepton(int idBeamIn = 11) : PDF(idBeamIn) {}
419
420private:
421
422  // Constants: could only be changed in the code itself.
423  static const double ALPHAEM, ME, MMU, MTAU;
424
425  // Update PDF values.
426  void xfUpdate(int id, double x, double Q2);
427
428  // The squared lepton mass, set at initialization.
429  double m2Lep;
430
431};
432 
433//==========================================================================
434
435// Gives electron (or other lepton) parton distribution when unresolved.
436 
437class LeptonPoint : public PDF {
438
439public:
440
441  // Constructor.
442  LeptonPoint(int idBeamIn = 11) : PDF(idBeamIn) {}
443
444private:
445
446  // Update PDF values in trivial way.
447  void xfUpdate(int , double , double ) {xlepton = 1; xgamma = 0.;}
448
449};
450 
451//==========================================================================
452
453// Gives neutrino parton distribution when unresolved (only choice for now).
454// Note factor of 2 since only lefthanded implies no spin averaging.
455 
456class NeutrinoPoint : public PDF {
457
458public:
459
460  // Constructor.
461  NeutrinoPoint(int idBeamIn = 12) : PDF(idBeamIn) {}
462
463private:
464
465  // Update PDF values, with spin factor of 2.
466  void xfUpdate(int , double , double ) {xlepton = 2; xgamma = 0.;}
467
468};
469 
470//==========================================================================
471
472} // end namespace Pythia8
473
474#endif // Pythia8_PartonDistributions_H
Note: See TracBrowser for help on using the repository browser.