8 #ifndef Pythia8_LHAPDF6_H
9 #define Pythia8_LHAPDF6_H
11 #include "Pythia8/PartonDistributions.h"
12 #include "LHAPDF/LHAPDF.h"
30 PdfSets(
string setName) : info(::LHAPDF::PDFSet(setName)),
31 pdfs(vector< ::LHAPDF::PDF* >(info.size(), 0)) {;}
34 ::LHAPDF::PDF *operator[](
unsigned int member) {
35 if (!pdfs[member]) pdfs[member] = info.mkPDF(member);
40 int size() {
return pdfs.size();}
43 ::LHAPDF::PDFSet info;
44 vector< ::LHAPDF::PDF* > pdfs;
58 LHAPDF6(
int idBeamIn,
string setName,
int member,
int)
59 :
PDF(idBeamIn), pdf(nullptr), extrapol(false)
60 { init(setName, member); }
63 void setExtrapolate(
bool extrapolIn) {extrapol = extrapolIn;}
70 ::LHAPDF::Extrapolator *ext;
74 void init(
string setName,
int member);
77 void xfUpdate(
int id,
double x,
double Q2);
80 bool insideBounds(
double x,
double Q2) {
81 return (x > xMin && x < xMax && Q2 > q2Min && Q2 < q2Max);}
84 double alphaS(
double Q2) {
return pdf->alphasQ2(Q2); }
87 double muPDFSave, mdPDFSave, mcPDFSave, msPDFSave, mbPDFSave,
88 xMin, xMax, q2Min, q2Max;
89 double mQuarkPDF(
int 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;
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;
108 int nMembers() {
return nMembersSave; }
116 const double LHAPDF6::PDFMINVALUE = 1e-10;
122 void LHAPDF6::init(
string setName,
int member) {
127 int id = ::LHAPDF::lookupLHAPDFID(setName, 0);
129 cout <<
"Error in LHAPDF6::init: unknown PDF "
133 pdfs = PdfSets(setName);
134 if (pdfs.size() == 0) {
135 cout <<
"Error in LHAPDF6::init: could not initialize PDF "
138 }
else if (member >= pdfs.size()) {
139 cout <<
"Error in LHAPDF6::init: " << setName
140 <<
" does not contain requested member" << endl;
149 q2Max = pdf->q2Max();
150 q2Min = pdf->q2Min();
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");
159 nMembersSave = pdf->info().get_entry_as<
int>(
"NumMembers");
167 void LHAPDF6::xfUpdate(
int,
double x,
double Q2) {
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;
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);
202 void LHAPDF6::calcPDFEnvelope(
int idNow,
double xNow,
double Q2NowIn,
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;
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);
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;
237 void LHAPDF6::calcPDFEnvelope(pair<int,int> idNows, pair<double,double> xNows,
238 double Q2NowIn,
int valSea) {
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;
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);
260 xfCalc[iMem] = max(0.0, xfCalc[iMem]);
261 if (valSea == 0 || (idNows.second != 1 && idNows.second != 2)) {
263 (PDFMINVALUE, pdfs[iMem]->xfxQ2(idNows.second, x2, Q2Now));
264 }
else if (valSea == 1 && (idNows.second == 1 || idNows.second == 2 )) {
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 )) {
270 (pdfs[iMem]->xfxQ2(-idNows.second, x2, Q2Now), PDFMINVALUE);
272 pdfEnvelope.pdfMemberVars[iMem] = xfCalc[iMem];
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;
291 LHAPDF6* newPDF(
int idBeamIn,
string setName,
int member) {
292 return new LHAPDF6(idBeamIn, setName, member, 1);}
294 void deletePDF(LHAPDF6* pdf) {
delete pdf;}
302 #endif // end Pythia8_LHAPDF6_H