StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
LHAPDF6.h
1 // LHAPDF6.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 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 // This file contains the LHAPDF6 PDF plugin class.
7 
8 #ifndef Pythia8_LHAPDF6_H
9 #define Pythia8_LHAPDF6_H
10 
11 #include "Pythia8/PartonDistributions.h"
12 #include "LHAPDF/LHAPDF.h"
13 
14 namespace Pythia8 {
15 
16 //==========================================================================
17 
18 // Containers for PDF sets.
19 
20 //--------------------------------------------------------------------------
21 
22 // Class to hold a PDF set, its information, and its uncertainty sets.
23 
24 class PdfSets {
25 
26 public:
27 
28  // Constructors.
29  PdfSets() {;}
30  PdfSets(string setName) : info(::LHAPDF::PDFSet(setName)),
31  pdfs(vector< ::LHAPDF::PDF* >(info.size(), 0)) {;}
32 
33  // Access a PDF set.
34  ::LHAPDF::PDF *operator[](unsigned int member) {
35  if (!pdfs[member]) pdfs[member] = info.mkPDF(member);
36  return pdfs[member];
37  }
38 
39  // Get number of PDF sets.
40  int size() {return pdfs.size();}
41 
42  // PDF sets and info.
43  ::LHAPDF::PDFSet info;
44  vector< ::LHAPDF::PDF* > pdfs;
45 
46 };
47 
48 
49 //==========================================================================
50 
51 // Provide interface to the LHAPDF6 library of parton densities.
52 
53 class LHAPDF6 : public PDF {
54 
55 public:
56 
57  // Constructor.
58  LHAPDF6(int idBeamIn, string setName, int member, int)
59  : PDF(idBeamIn), pdf(nullptr), extrapol(false)
60  { init(setName, member); }
61 
62  // Allow extrapolation beyond boundaries (not implemented).
63  void setExtrapolate(bool extrapolIn) {extrapol = extrapolIn;}
64 
65 private:
66 
67  // The LHAPDF objects.
68  PdfSets pdfs;
69  ::LHAPDF::PDF *pdf;
70  ::LHAPDF::Extrapolator *ext;
71  bool extrapol;
72 
73  // Initialization of PDF set.
74  void init(string setName, int member);
75 
76  // Update parton densities.
77  void xfUpdate(int id, double x, double Q2);
78 
79  // Check whether x and Q2 values fall inside the fit bounds.
80  bool insideBounds(double x, double Q2) {
81  return (x > xMin && x < xMax && Q2 > q2Min && Q2 < q2Max);}
82 
83  // Return the running alpha_s shipped with the LHAPDF set.
84  double alphaS(double Q2) { return pdf->alphasQ2(Q2); }
85 
86  // Return quark masses used in the PDF fit.
87  double muPDFSave, mdPDFSave, mcPDFSave, msPDFSave, mbPDFSave,
88  xMin, xMax, q2Min, q2Max;
89  double mQuarkPDF(int id) {
90  switch(abs(id)){
91  case 1: return mdPDFSave;
92  case 2: return muPDFSave;
93  case 3: return msPDFSave;
94  case 4: return mcPDFSave;
95  case 5: return mbPDFSave;
96  }
97  return -1.;
98  }
99 
100  // Calculate uncertainties using the LHAPDF prescription.
101  void calcPDFEnvelope(int, double, double, int);
102  void calcPDFEnvelope(pair<int,int>, pair<double,double>, double, int);
103  PDFEnvelope pdfEnvelope;
104  PDFEnvelope getPDFEnvelope() {return pdfEnvelope;}
105  static const double PDFMINVALUE;
106 
107  int nMembersSave;
108  int nMembers() { return nMembersSave; }
109 
110 };
111 
112 //--------------------------------------------------------------------------
113 
114 // Constants.
115 
116 const double LHAPDF6::PDFMINVALUE = 1e-10;
117 
118 //--------------------------------------------------------------------------
119 
120 // Initialize a parton density function from LHAPDF6.
121 
122 void LHAPDF6::init(string setName, int member) {
123  isSet = false;
124 
125 
126  // Find the PDF set.
127  int id = ::LHAPDF::lookupLHAPDFID(setName, 0);
128  if (id < 0) {
129  cout << "Error in LHAPDF6::init: unknown PDF "
130  << setName << endl;
131  return;
132  }
133  pdfs = PdfSets(setName);
134  if (pdfs.size() == 0) {
135  cout << "Error in LHAPDF6::init: could not initialize PDF "
136  << setName << endl;
137  return;
138  } else if (member >= pdfs.size()) {
139  cout << "Error in LHAPDF6::init: " << setName
140  << " does not contain requested member" << endl;
141  return;
142  }
143  pdf = pdfs[member];
144  isSet = true;
145 
146  // Save x and Q2 limits.
147  xMax = pdf->xMax();
148  xMin = pdf->xMin();
149  q2Max = pdf->q2Max();
150  q2Min = pdf->q2Min();
151 
152  // Store quark masses used in PDF fit.
153  muPDFSave = pdf->info().get_entry_as<double>("MUp");
154  mdPDFSave = pdf->info().get_entry_as<double>("MDown");
155  mcPDFSave = pdf->info().get_entry_as<double>("MCharm");
156  msPDFSave = pdf->info().get_entry_as<double>("MStrange");
157  mbPDFSave = pdf->info().get_entry_as<double>("MBottom");
158 
159  nMembersSave = pdf->info().get_entry_as<int>("NumMembers");
160 
161 }
162 
163 //--------------------------------------------------------------------------
164 
165 // Give the parton distribution function set from LHAPDF6.
166 
167 void LHAPDF6::xfUpdate(int, double x, double Q2) {
168 
169  // Freeze at boundary value if PDF is evaluated outside the fit region.
170  if (x < xMin && !extrapol) x = xMin;
171  if (x > xMax) x = xMax;
172  if (Q2 < q2Min) Q2 = q2Min;
173  if (Q2 > q2Max) Q2 = q2Max;
174 
175  // Update values.
176  xg = pdf->xfxQ2(21, x, Q2);
177  xu = pdf->xfxQ2(2, x, Q2);
178  xd = pdf->xfxQ2(1, x, Q2);
179  xs = pdf->xfxQ2(3, x, Q2);
180  xubar = pdf->xfxQ2(-2, x, Q2);
181  xdbar = pdf->xfxQ2(-1, x, Q2);
182  xsbar = pdf->xfxQ2(-3, x, Q2);
183  xc = pdf->xfxQ2(4, x, Q2);
184  xb = pdf->xfxQ2(5, x, Q2);
185  xgamma = pdf->xfxQ2(22, x, Q2);
186 
187  // Subdivision of valence and sea.
188  xuVal = xu - xubar;
189  xuSea = xubar;
190  xdVal = xd - xdbar;
191  xdSea = xdbar;
192 
193  // idSav = 9 to indicate that all flavours reset.
194  idSav = 9;
195 
196 }
197 
198 //--------------------------------------------------------------------------
199 
200 // Calculate uncertainties using the LHAPDF prescription.
201 
202 void LHAPDF6::calcPDFEnvelope(int idNow, double xNow, double Q2NowIn,
203  int valSea) {
204 
205  // Freeze at boundary value if PDF is evaluated outside the fit region.
206  double x1 = (xNow < xMin && !extrapol) ? xMin : xNow;
207  if (x1 > xMax) x1 = xMax;
208  double Q2Now = (Q2NowIn < q2Min) ? q2Min : Q2NowIn;
209  if (Q2Now > q2Max) Q2Now = q2Max;
210 
211  // Loop over the members.
212  vector<double> xfCalc(pdfs.size());
213  for(int iMem = 0; iMem < pdfs.size(); ++iMem) {
214  if (valSea==0 || (idNow != 1 && idNow != 2)) {
215  xfCalc[iMem] = pdfs[iMem]->xfxQ2(idNow, x1, Q2Now);
216  } else if (valSea==1 && (idNow == 1 || idNow == 2 )) {
217  xfCalc[iMem] = pdfs[iMem]->xfxQ2(idNow, x1, Q2Now) -
218  pdfs[iMem]->xfxQ2(-idNow, x1, Q2Now);
219  } else if (valSea==2 && (idNow == 1 || idNow == 2 )) {
220  xfCalc[iMem] = pdfs[iMem]->xfxQ2(-idNow, x1, Q2Now);
221  }
222  }
223 
224  // Calculate the uncertainty.
225  ::LHAPDF::PDFUncertainty xfErr = pdfs.info.uncertainty(xfCalc);
226  pdfEnvelope.centralPDF = xfErr.central;
227  pdfEnvelope.errplusPDF = xfErr.errplus;
228  pdfEnvelope.errminusPDF = xfErr.errminus;
229  pdfEnvelope.errsymmPDF = xfErr.errsymm;
230  pdfEnvelope.scalePDF = xfErr.scale;
231 }
232 
233 //--------------------------------------------------------------------------
234 
235 // Calculate uncertainties using the LHAPDF prescription.
236 
237 void LHAPDF6::calcPDFEnvelope(pair<int,int> idNows, pair<double,double> xNows,
238  double Q2NowIn, int valSea) {
239 
240  // Freeze at boundary value if PDF is evaluated outside the fit region.
241  double x1 = (xNows.first < xMin && !extrapol) ? xMin : xNows.first;
242  if (x1 > xMax) x1 = xMax;
243  double x2 = (xNows.second < xMin && !extrapol) ? xMin : xNows.second;
244  if (x2 > xMax) x2 = xMax;
245  double Q2Now = (Q2NowIn < q2Min) ? q2Min : Q2NowIn;
246  if (Q2Now > q2Max) Q2Now = q2Max;
247 
248  // Loop over the members.
249  vector<double> xfCalc(pdfs.size());
250  pdfEnvelope.pdfMemberVars.resize(pdfs.size());
251  for(int iMem = 0; iMem < pdfs.size(); ++iMem) {
252  if (valSea == 0 || (idNows.first != 1 && idNows.first != 2 ) ) {
253  xfCalc[iMem] = pdfs[iMem]->xfxQ2(idNows.first, x1, Q2Now);
254  } else if (valSea == 1 && (idNows.first == 1 || idNows.first == 2)) {
255  xfCalc[iMem] = pdfs[iMem]->xfxQ2(idNows.first, x1, Q2Now)
256  - pdfs[iMem]->xfxQ2(-idNows.first, x1, Q2Now);
257  } else if (valSea == 2 && (idNows.first == 1 || idNows.first == 2 )) {
258  xfCalc[iMem] = pdfs[iMem]->xfxQ2(-idNows.first, x1, Q2Now);
259  }
260  xfCalc[iMem] = max(0.0, xfCalc[iMem]);
261  if (valSea == 0 || (idNows.second != 1 && idNows.second != 2)) {
262  xfCalc[iMem] /= max
263  (PDFMINVALUE, pdfs[iMem]->xfxQ2(idNows.second, x2, Q2Now));
264  } else if (valSea == 1 && (idNows.second == 1 || idNows.second == 2 )) {
265  xfCalc[iMem] /= max
266  (pdfs[iMem]->xfxQ2(idNows.second, x2, Q2Now) - pdfs[iMem]->xfxQ2
267  (-idNows.second, x2, Q2Now), PDFMINVALUE);
268  } else if (valSea == 2 && (idNows.second == 1 || idNows.second == 2 )) {
269  xfCalc[iMem] /= max
270  (pdfs[iMem]->xfxQ2(-idNows.second, x2, Q2Now), PDFMINVALUE);
271  }
272  pdfEnvelope.pdfMemberVars[iMem] = xfCalc[iMem];
273  }
274 
275  // Calculate the uncertainty.
276  ::LHAPDF::PDFUncertainty xfErr = pdfs.info.uncertainty(xfCalc);
277  pdfEnvelope.centralPDF = xfErr.central;
278  pdfEnvelope.errplusPDF = xfErr.errplus;
279  pdfEnvelope.errminusPDF = xfErr.errminus;
280  pdfEnvelope.errsymmPDF = xfErr.errsymm;
281  pdfEnvelope.scalePDF = xfErr.scale;
282 
283 }
284 
285 //--------------------------------------------------------------------------
286 
287 // Define external handles to the plugin for dynamic loading.
288 
289 extern "C" {
290 
291  LHAPDF6* newPDF(int idBeamIn, string setName, int member) {
292  return new LHAPDF6(idBeamIn, setName, member, 1);}
293 
294  void deletePDF(LHAPDF6* pdf) {delete pdf;}
295 
296 }
297 
298 //==========================================================================
299 
300 } // end namespace Pythia8
301 
302 #endif // end Pythia8_LHAPDF6_H