StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
main54.cc
1 // main54.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 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 // This program compares the internal and LHAPDF implementations of the
7 // NNPDF 2.3 QCD+QED sets, for results and for timing.
8 // Author: Juan Rojo.
9 
10 #include "Pythia8/Pythia.h"
11 using namespace Pythia8;
12 
13 int main() {
14 
15  cout<<"\n NNPDF2.3 QED LO phenomenology \n "<<endl;
16  cout<<"\n Check access to NNPDF2.3 LO QED sets \n "<<endl;
17 
18  // Generator.
19  Pythia pythia;
20 
21  // Access the PDFs.
22  int idBeamIn = 2212;
23  string xmlPath = "../xmldoc/";
24  Info* infoPtr = 0;
25  int member=0;
26 
27  // Grid of studied points.
28  string xpdf[] = {"x*g","x*d","x*u","x*s"};
29  double xlha[] = {1e-5, 1e-1};
30  double Q2[] = { 2.0, 10000.0 };
31  string setName;
32  string setName_lha;
33 
34  // For timing checks.
35  int const nq = 200;
36  int const nx = 200;
37 
38  // Loop over all internal PDF sets in Pythia8
39  // and compare with their LHAPDF correspondents.
40  for (int iFitIn = 1; iFitIn < 5; iFitIn++) {
41 
42  // Constructor for internal PDFs.
43  NNPDF pdfs_nnpdf( idBeamIn, iFitIn, xmlPath, infoPtr);
44 
45  // Constructor for LHAPDF.
46  if (iFitIn == 1) setName = "NNPDF23_lo_as_0130_qed_mem0.LHgrid";
47  if (iFitIn == 2) setName = "NNPDF23_lo_as_0119_qed_mem0.LHgrid";
48  if (iFitIn == 3) setName = "NNPDF23_nlo_as_0119_qed_mc_mem0.LHgrid";
49  if (iFitIn == 4) setName = "NNPDF23_nnlo_as_0119_qed_mc_mem0.LHgrid";
50  LHAPDF pdfs_nnpdf_lha( idBeamIn, setName, member);
51  cout << "\n PDF set = " << setName << " \n" << endl;
52 
53  // Check quarks and gluons.
54  for (int f = 0; f < 4; f++) {
55  for (int iq = 0; iq < 2; iq++) {
56  cout << " " << xpdf[f] << ", Q2 = " << Q2[iq] << endl;
57  cout << " x \t Pythia8\t LHAPDF\t diff(%) " << endl;
58  for (int ix = 0; ix < 2; ix++) {
59  double a = pdfs_nnpdf.xf( f, xlha[ix], Q2[iq]);
60  double b = pdfs_nnpdf_lha.xf( f, xlha[ix], Q2[iq]);
61  double diff = 1e2 * fabs((a-b)/b);
62  cout << scientific << xlha[ix] << " " << a << " " << b
63  << " " << diff << endl;
64  if (diff > 1e-8) exit(-10);
65  }
66  }
67  }
68 
69  // Check photon.
70  cout << "\n Now checking the photon PDF \n" << endl;
71  for (int iq = 0; iq < 2; iq++) {
72  cout << " " << "x*gamma" << ", Q2 = " << Q2[iq] << endl;
73  cout << " x \t Pythia8\t LHAPDF\t diff(%)" << endl;
74  for (int ix = 0; ix < 2; ix++) {
75  double a = pdfs_nnpdf.xf( 22, xlha[ix], Q2[iq]);
76  double b = pdfs_nnpdf_lha.xf( 22, xlha[ix], Q2[iq]);
77  double diff = 1e2 * fabs((a-b)/b);
78  cout << scientific << xlha[ix] << " " << a << " " << b
79  << " " << diff << endl;
80  if(diff > 1e-8) exit(-10);
81  }
82  }
83 
84  // Now check the timings for evolution.
85  cout << "\n Checking timings " << endl;
86 
87  clock_t t1 = clock();
88  for (int f = -4; f < 4; f++) {
89  for (int iq = 0; iq < nq; iq++) {
90  double qq2 = 2.0 * pow( 1e6 / 2.0, double(iq)/nq);
91  for (int ix = 0; ix < nx; ix++) {
92  double xx = 1e-6 * pow( 9e-1 / 1e-6, double(ix)/nx);
93  pdfs_nnpdf.xf(f,xx,qq2);
94  }
95  }
96  }
97  clock_t t2 = clock();
98  cout << " NNPDF internal timing = " << (t2-t1)/(double)CLOCKS_PER_SEC
99  << endl;
100 
101  t1=clock();
102  for (int f = -4; f < 4; f++) {
103  for (int iq = 0; iq < nq; iq++) {
104  double qq2 = 2.0 * pow(1e6 / 2.0, double(iq)/nq);
105  for (int ix = 0; ix < nx; ix++) {
106  double xx = 1e-6 * pow( 9e-1 / 1e-6, double(ix)/nx);
107  pdfs_nnpdf_lha.xf(f,xx,qq2);
108  }
109  }
110  }
111  t2=clock();
112  cout << " NNPDF LHAPDF timing = " << (t2-t1)/(double)CLOCKS_PER_SEC
113  << endl;
114 
115  } // End loop over NNPDF internal sets
116 
117  // Done.
118  cout << "\n Checked that LHAPDF and Internal Pythia8 give identical"
119  << " results\n" << endl;
120 
121  return 0;
122 }