StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
LHAPDF5.h
1 // LHAPDF5.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 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 LHAPDF5 PDF plugin class.
7 
8 #ifndef Pythia8_LHAPDF5_H
9 #define Pythia8_LHAPDF5_H
10 
11 #include "Pythia8/PartonDistributions.h"
12 
13 namespace Pythia8 {
14 
15 //==========================================================================
16 
17 // Interfaces to to make the C++ calls similar to f77.
18 
19 //--------------------------------------------------------------------------
20 
21 // Declare the LHAPDF5 f77 subroutines that are needed.
22 
23 extern "C" {
24 
25  extern void initpdfsetm_(int&, const char*, int);
26 
27  extern void initpdfsetbynamem_(int&, const char*, int);
28 
29  extern void initpdfm_(int&, int&);
30 
31  extern void evolvepdfm_(int&, double&, double&, double*);
32 
33  extern void evolvepdfpm_(int&, double&, double&, double&, double&, double*);
34 
35  extern void evolvepdfphotonm_(int&, double&, double&, double*, double&);
36 
37  extern void setlhaparm_(const char*, int);
38 
39  extern void getxminm_(int &, int &, double &);
40 
41  extern void getxmaxm_(int &, int &, double &);
42 
43  extern void getq2minm_(int &, int &, double &);
44 
45  extern void getq2maxm_(int &, int &, double &);
46 
47 }
48 
49 //--------------------------------------------------------------------------
50 
51 // Map the f77 routines to C++.
52 
53 namespace LHAPDF5Interface {
54 
55  // Initialize set with full pathname, allowing multiple sets.
56  void initPDFsetM( int& nSet, string name) {
57  const char* cName = name.c_str(); int lenName = name.length();
58  initpdfsetm_( nSet, cName, lenName);
59  }
60 
61  // Initialize set with simple name, allowing multiple sets.
62  void initPDFsetByNameM( int& nSet, string name) {
63  const char* cName = name.c_str(); int lenName = name.length();
64  initpdfsetbynamem_( nSet, cName, lenName);
65  }
66 
67  // Initialize member of set.
68  void initPDFM(int& nSet, int member) {
69  initpdfm_(nSet, member);
70  }
71 
72  // Evaluate x f_i(x, Q).
73  void evolvePDFM( int& nSet, double x, double Q, double* xfArray) {
74  evolvepdfm_( nSet, x, Q, xfArray);
75  }
76 
77  // Evaluate x f_i(x, Q) for photon beams.
78  void evolvePDFpM( int& nSet, double x, double Q, double P2, double IP2,
79  double* xfArray) {
80  evolvepdfpm_( nSet, x, Q, P2, IP2, xfArray);
81  }
82 
83  // Evaluate x f_i(x, Q) including photon
84  void evolvePDFPHOTONM(int& nSet, double x, double Q, double* xfArray,
85  double& xPhoton) {
86  evolvepdfphotonm_( nSet, x, Q, xfArray, xPhoton);
87  }
88 
89  // Extrapolate PDF set beyond boundaries, or freeze them there.
90  void setPDFparm(string name) {
91  const char* cName = name.c_str(); int lenName = name.length();
92  setlhaparm_( cName, lenName);
93  }
94 
95  // Simple structure to hold LHAPDF set information.
96  struct LHAPDFInfo {
97  string name;
98  int member;
99  bool photon;
100  };
101 
102  // Global tracking of opened PDF sets.
103  map<int, LHAPDFInfo> initializedSets;
104 
105  // Method to find the nSet number corresponding to a name and member.
106  // Returns -1 if no such LHAPDF5 set has been initialized.
107  int findNSet(string setName, int member) {
108  for (map<int, LHAPDFInfo>::const_iterator i = initializedSets.begin();
109  i != initializedSets.end(); ++i) {
110  int iSet = i->first;
111  string iName = i->second.name;
112  int iMember = i->second.member;
113  if (iName == setName && iMember == member) return iSet;
114  }
115  return -1;
116  }
117 
118  // Method to return the lowest non-occupied nSet number.
119  int freeNSet() {
120  for (int iSet = 1; iSet <= int(initializedSets.size()); ++iSet) {
121  if (initializedSets.find(iSet) == initializedSets.end())
122  return iSet;
123  }
124  return initializedSets.size() + 1;
125  }
126 
127 }
128 
129 //==========================================================================
130 
131 // Plugin interface to the LHAPDF5 library.
132 
133 //==========================================================================
134 
135 // Provide plugin interface to the LHAPDF5 library of parton densities.
136 
137 class LHAPDF5 : public PDF {
138 
139 public:
140 
141  // Constructor.
142  LHAPDF5(int idBeamIn, string setName, int member, int nSetIn = -1,
143  Info* infoPtr = 0) : PDF(idBeamIn), doExtraPol(false), nSet(nSetIn)
144  { init(setName, member, infoPtr);
145  isPhoton = (idBeamIn == 22) ? true : false; }
146 
147  // Allow extrapolation beyond boundaries. This is optional.
148  void setExtrapolate(bool extrapol);
149 
150 private:
151 
152  // Initialization of PDF set.
153  void init(string setName, int member, Info* infoPtr);
154 
155  // Update all PDF values.
156  void xfUpdate(int , double x, double Q2);
157 
158  // Current set and pdf values.
159  bool doExtraPol;
160  int nSet;
161  double xfArray[13];
162  bool hasPhoton, isPhoton;
163  double xPhoton;
164 
165 };
166 
167 //--------------------------------------------------------------------------
168 
169 // Initialize a parton density function from LHAPDF5.
170 
171 void LHAPDF5::init(string setName, int member, Info*) {
172 
173  // If already initialized then need not do anything further.
174  LHAPDF5Interface::LHAPDFInfo initializedInfo =
175  LHAPDF5Interface::initializedSets[nSet];
176  string initializedSetName = initializedInfo.name;
177  int initializedMember = initializedInfo.member;
178  hasPhoton = initializedInfo.photon;
179  if (setName == initializedSetName && member == initializedMember) return;
180 
181  // Initialize set. If first character is '/' then assume that name
182  // is given with path, else not.
183  if (setName[0] == '/') LHAPDF5Interface::initPDFsetM( nSet, setName);
184  else LHAPDF5Interface::initPDFsetByNameM( nSet, setName);
185  isSet = (nSet >= 0);
186 
187  // Initialize member.
188  LHAPDF5Interface::initPDFM(nSet, member);
189 
190  // Do not collect statistics on under/overflow to save time and space.
191  LHAPDF5Interface::setPDFparm( "NOSTAT" );
192  LHAPDF5Interface::setPDFparm( "LOWKEY" );
193 
194  // Check if photon PDF available (has_photon does not work properly).
195  xPhoton = 0;
196  LHAPDF5Interface::evolvePDFPHOTONM(nSet, 0.01, 1, xfArray, xPhoton);
197  hasPhoton = xPhoton != 0;
198 
199  // Save values to avoid unnecessary reinitializations.
200  initializedInfo.name = setName;
201  initializedInfo.member = member;
202  initializedInfo.photon = hasPhoton;
203  if (nSet > 0) LHAPDF5Interface::initializedSets[nSet] = initializedInfo;
204 
205 }
206 
207 //--------------------------------------------------------------------------
208 
209 // Allow optional extrapolation beyond boundaries.
210 
211 void LHAPDF5::setExtrapolate(bool extrapol) {
212 
213  doExtraPol = extrapol;
214  LHAPDF5Interface::setPDFparm( (extrapol) ? "EXTRAPOLATE" : "18" );
215 
216 }
217 
218 //--------------------------------------------------------------------------
219 
220 // Give the parton distribution function set from LHAPDF5.
221 
222 void LHAPDF5::xfUpdate(int, double x, double Q2) {
223 
224  // Freeze at boundary value if PDF is evaluated outside the fit region.
225  int member = LHAPDF5Interface::initializedSets[nSet].member;
226  double xMin, xMax, q2Min, q2Max;
227  getxminm_( nSet, member, xMin);
228  getxmaxm_( nSet, member, xMax);
229  getq2minm_(nSet, member, q2Min);
230  getq2maxm_(nSet, member, q2Max);
231  if (x < xMin && !doExtraPol) x = xMin;
232  if (x > xMax) x = xMax;
233  if (Q2 < q2Min) Q2 = q2Min;
234  if (Q2 > q2Max) Q2 = q2Max;
235  double Q = sqrt( max( 0., Q2));
236 
237  // Use special call if photon included in proton.
238  if (hasPhoton) {
239  LHAPDF5Interface::evolvePDFPHOTONM( nSet, x, Q, xfArray, xPhoton);
240  }
241 
242  // Use special call with photon beams. No virtualities implemented yet.
243  else if (isPhoton) {
244  LHAPDF5Interface::evolvePDFpM( nSet, x, Q, 0., 0., xfArray);
245  }
246 
247  // Else use default LHAPDF5 call.
248  else {
249  LHAPDF5Interface::evolvePDFM( nSet, x, Q, xfArray);
250  xPhoton=0.0;
251  }
252 
253  // Update values.
254  xg = xfArray[6];
255  xu = xfArray[8];
256  xd = xfArray[7];
257  xs = xfArray[9];
258  xubar = xfArray[4];
259  xdbar = xfArray[5];
260  xsbar = xfArray[3];
261  xc = xfArray[10];
262  xb = xfArray[11];
263  xgamma = xPhoton;
264 
265  // Subdivision of valence and sea.
266  xuVal = xu - xubar;
267  xuSea = xubar;
268  xdVal = xd - xdbar;
269  xdSea = xdbar;
270 
271  // idSav = 9 to indicate that all flavours reset.
272  idSav = 9;
273 
274 }
275 
276 //--------------------------------------------------------------------------
277 
278 // Define external handles to the plugin for dynamic loading.
279 
280 extern "C" {
281 
282  LHAPDF5* newLHAPDF(int idBeamIn, string setName, int member,
283  Info* infoPtr) {
284  int nSet = LHAPDF5Interface::findNSet(setName, member);
285  if (nSet == -1) nSet = LHAPDF5Interface::freeNSet();
286  return new LHAPDF5(idBeamIn, setName, member, nSet, infoPtr);
287  }
288 
289  void deleteLHAPDF(LHAPDF5* pdf) {
290  delete pdf;
291  }
292 
293 }
294 
295 //==========================================================================
296 
297 } // end namespace Pythia8
298 
299 #endif // end Pythia8_LHAPDF5_H