lhapdf is hosted by Hepforge, IPPP Durham
LHAPDF 6.5.5
Loading...
Searching...
No Matches
LHAGlue.h
Go to the documentation of this file.
1// -*- C++ -*-
2//
3// This file is part of LHAPDF
4// Copyright (C) 2012-2024 The LHAPDF collaboration (see AUTHORS for details)
5//
6#pragma once
7#ifndef LHAPDF_LHAGlue_H
8#define LHAPDF_LHAGlue_H
9
10/// @file LHAGlue.h
11/// A file that provides backwards compatibility for some C functions from LHAPDF 5.x
12
13#include "LHAPDF/Version.h"
14#if LHAPDF_LHA5CXX
15
16/// A special C++ function to return the PDF name + code currently being used via LHAGlue.
17std::string lhaglue_get_current_pdf(int nset=1);
18
19// Compatibility preprocessor-based aliasing of deprecated "M" function names
20#define initPDFSetM initPDFSet
21#define initPDFSetByNameM initPDFSetByName
22#define initPDFM initPDF
23#define initPDFByNameM initPDFByName
24#define getDescriptionM getDescription
25#define xfxM xfx
26#define xfxpM xfxp
27#define xfxaM xfxa
28#define xfxphotonM xfxphoton
29#define numberPDFM numberPDF
30#define alphasPDFM alphasPDF
31#define getOrderPDFM getOrderPDF
32#define getOrderAlphaSM getOrderAlphaS
33#define getQMassM getQMass
34#define getThresholdM getThreshold
35#define getNfM getNf
36#define getLam4M getLam4
37#define getLam5M getLam5
38#define getXminM getXmin
39#define getXmaxM getXmax
40#define getQ2minM getQ2min
41#define getQ2maxM getQ2max
42
43namespace LHAPDF {
44
45
46 /// @brief Only Provided for LHAPDF5 compatibility. Distinction between evolution
47 /// or interpolation PDF sets.
48 /// Enum to choose whether evolution (i.e. @c LHpdf data file) or
49 /// interpolation (i.e. @c LHgrid data file) is used.
50 /// This distinction ismeaningless in LHAPDF6 interpolation is always used.
51 enum SetType {
52 EVOLVE = 0, LHPDF = 0,
53 INTERPOLATE = 1, LHGRID = 1
54 };
55
56 /// Level of noisiness.
57 enum Verbosity { SILENT=0, LOWKEY=1, DEFAULT=2 };
58
59 /// Get LHAPDF version string (prefer LHAPDF::version())
60 inline std::string getVersion() {
61 return version();
62 }
63
64 /// Get max allowed number of concurrent sets (there is no limit anymore)
65 inline int getMaxNumSets() { return 1000; }
66
67 /// Global initialisation (there is none)
68 inline void initLHAPDF() {}
69
70 /// Extrapolate beyond grid edges (not an option at present)
71 /// @todo Use this to set the default extrapolator when there is a physical extrapolation option
72 // This form commented due to unused variable warnings, until this can actually have an effect
73 // inline void extrapolate(bool extrapolate=true) {}
74 inline void extrapolate(bool) {}
75 inline void extrapolate() {}
76
77 /// Choose level of noisiness.
78 void setVerbosity(Verbosity noiselevel);
79
80 /// Set a steering parameter (does nothing!)
81 inline void setParameter(const std::string&) {
82 std::cerr << "LHAPDF::setParameter() has no effect in LHAPDF6: "
83 << "please update your code to use the native C++ interface" << std::endl;
84 }
85
86 /// Check if the PDF includes a photon member
87 bool hasPhoton();
88
89 /// Prepends path to path list
90 void setPDFPath(const string& path);
91 std::string pdfsetsPath();
92
93
94 /// The PDF set by filename, see subdir @c PDFsets of LHAPDF for choices.
95 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
96 void initPDFSetByName(const std::string& filename);
97 void initPDFSetByName(const std::string& filename, SetType type);
98
99 /// The PDF set by filename, see subdir @c PDFsets of LHAPDF for choices.
100 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
101 void initPDFSetByName(int nset, const std::string& filename);
102 void initPDFSetByName(int nset, const std::string& filename, SetType type);
103
104
105 /// Number of members available in the current set.
106 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
107 int numberPDF();
108
109 /// Number of members available in the current set.
110 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
111 int numberPDF(int nset);
112
113 /// The choice of PDF member out of one distribution.
114 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
115 void initPDF(int memset);
116
117 /// The choice of PDF member out of one distribution.
118 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
119 void initPDF(int nset, int memset);
120
121
122 /// Nucleon PDF: returns \f$ x f(x, Q) \f$ for flavour @a fl - flavour encoding as in the LHAPDF manual.
123 /// @arg -6..-1 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
124 /// @arg 0 = \f$ g \f$
125 /// @arg 1..6 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
126 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
127 double xfx(double x, double Q, int fl);
128
129 /// Nucleon PDF: returns @c x f(x, Q) for flavour @a fl - flavour encoding as in the LHAPDF manual.
130 /// @arg -6..-1 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
131 /// @arg 0 = \f$ g \f$
132 /// @arg 1..6 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
133 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
134 double xfx(int nset, double x, double Q, int fl);
135
136 /// Nucleon PDF: fills primitive 13 element array pointed at by @a results with
137 /// \f$ x f(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
138 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
139 /// @arg 6 = \f$ g \f$;
140 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
141 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
142 void xfx(double x, double Q, double* results);
143
144 /// Nucleon PDF: fills primitive 13 element array pointed at by @a results with
145 /// \f$ x f(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
146 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
147 /// @arg 6 = \f$ g \f$;
148 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
149 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
150 void xfx(int nset, double x, double Q, double* results);
151
152 /// Nucleon PDF: returns a vector \f$ x f_i(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
153 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
154 /// @arg 6 = \f$ g \f$;
155 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
156 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
157 std::vector<double> xfx(double x, double Q);
158
159 /// Nucleon PDF: returns a vector @c x f_i(x, Q) with index \f$ 0 < i < 12 \f$.
160 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
161 /// @arg 6 = \f$ g \f$;
162 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
163 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
164 std::vector<double> xfx(int nset, double x, double Q);
165
166
167 /// MRST QED PDF: returns a vector \f$ x f_i(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
168 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
169 /// @arg 6 = \f$ g \f$;
170 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$;
171 /// @arg 13 = \f$ \gamma \f$.
172 ///
173 /// NB. Note extra element in this set for MRST photon.
174 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
175 std::vector<double> xfxphoton(double x, double Q);
176
177 /// MRST QED PDF: returns a vector \f$ x f_i(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
178 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
179 /// @arg 6 = \f$ g \f$;
180 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$;
181 /// @arg 13 = \f$ \gamma \f$.
182 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
183 std::vector<double> xfxphoton(int nset, double x, double Q);
184
185 /// MRST QED PDF: fills primitive 14 element array pointed at by @a results with
186 /// \f$ x f(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
187 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
188 /// @arg 6 = \f$ g \f$;
189 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
190 /// @arg 13 = \f$ \gamma \f$.
191 ///
192 /// NB. Note extra element in this set for MRST photon.
193 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
194 void xfxphoton(double x, double Q, double* results);
195
196 /// MRST QED PDF: fills primitive 14 element array pointed at by @a results with
197 /// \f$ x f(x, Q) \f$ with index \f$ 0 < i < 12 \f$.
198 /// @arg 0..5 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
199 /// @arg 6 = \f$ g \f$;
200 /// @arg 7..12 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$.
201 /// @arg 13 = \f$ \gamma \f$.
202 ///
203 /// NB. Note extra element in this set for MRST photon.
204 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
205 void xfxphoton(int nset, double x, double Q, double* results);
206
207 /// MRST QED PDF: returns \f$ x f(x, Q) \f$ for flavour @a fl - this time the flavour encoding
208 /// is as in the LHAPDF manual.
209 /// @arg -6..-1 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
210 /// @arg 0 = \f$ g \f$
211 /// @arg 1..6 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$;
212 /// @arg 7 = \f$ \gamma \f$.
213 ///
214 /// NB. Note extra element in this set for MRST photon.
215 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
216 double xfxphoton(double x, double Q, int fl);
217
218 /// MRST QED PDF: returns \f$ x f(x, Q) \f$ for flavour @a fl - this time the flavour encoding
219 /// is as in the LHAPDF manual.
220 /// @arg -6..-1 = \f$ \bar{t} \f$, ..., \f$ \bar{u} \f$, \f$ \bar{d} \f$;
221 /// @arg 0 = \f$ g \f$
222 /// @arg 1..6 = \f$ d \f$, \f$ u \f$, ..., \f$ t \f$;
223 /// @arg 7 = \f$ \gamma \f$.
224 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
225 double xfxphoton(int nset, double x, double Q, int fl);
226
227
228 /// Print PDF description to stdout
229 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
230 void getDescription();
231
232 /// Print PDF description to stdout
233 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
234 void getDescription(int nset);
235
236
237 /// Order of \f$ \alpha_\mathrm{s} \f$ used by the current PDF.
238 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
239 int getOrderAlphaS();
240
241 /// Order of \f$ \alpha_\mathrm{s} \f$ used by the current PDF.
242 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
243 int getOrderAlphaS(int nset);
244
245
246 /// Order of QCD used by the current PDF.
247 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
248 int getOrderPDF();
249
250 /// Order of QCD used by the current PDF.
251 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
252 int getOrderPDF(int nset);
253
254
255 /// Number of flavours used in current PDF
256 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
257 int getNf(int nset);
258
259 /// Number of flavours used in current PDF
260 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
261 int getNf();
262
263
264 /// 4-flavour LambdaQCD used in current PDF, if available else -1.0
265 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
266 double getLam4(int nset);
267
268 /// 4-flavour LambdaQCD used in current PDF, if available else -1.0
269 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
270 double getLam4(int nset, int nmem);
271
272
273 /// 5-flavour LambdaQCD used in current PDF, if available else -1.0
274 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
275 double getLam5(int nset);
276
277 /// 5-flavour LambdaQCD used in current PDF, if available else -1.0
278 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
279 double getLam5(int nset, int nmem);
280
281
282 /// Minimum X for current PDF
283 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
284 double getXmin(int nmem);
285
286 /// Minimum X for current PDF
287 /// @deprecated Use the proper C++ interface of LHAPDF6 -instead!
288 double getXmin(int nset, int nmem);
289
290 /// Maximum X for current PDF
291 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
292 double getXmax(int nset, int nmem);
293
294 /// Maximum X for current PDF
295 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
296 double getXmax(int nmem);
297
298
299 /// Minimum Q2 for current PDF
300 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
301 double getQ2min(int nset, int nmem);
302
303 /// Minimum Q2 for current PDF
304 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
305 double getQ2min(int nmem);
306
307 /// Maximum Q2 for current PDF
308 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
309 double getQ2max(int nset, int nmem);
310
311 /// Maximum Q2 for current PDF
312 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
313 double getQ2max(int nmem);
314
315
316 /// Mass of quarks for current PDF
317 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
318 double getQMass(int nset, int nf);
319
320 /// Mass of quarks for current PDF
321 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
322 double getQMass(int nf);
323
324
325 /// Mass of quarks for current PDF
326 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
327 double getThreshold(int nset, int nf);
328
329 /// Mass of quarks for current PDF
330 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
331 double getThreshold(int nf);
332
333
334 /// \f$ \alpha_\mathrm{s} \f$ used by the current PDF.
335 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
336 double alphasPDF(double Q);
337
338 /// \f$ \alpha_\mathrm{s} \f$ used by the current PDF.
339 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
340 double alphasPDF(int nset, double Q);
341
342 /// @brief Use @a member in current PDF set.
343 /// This operation is computationally cheap.
344 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
345 void usePDFMember(int member);
346
347 /// @brief Use @a member in PDF set @a nset (multi-set version).
348 /// This operation is computationally cheap.
349 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
350 void usePDFMember(int nset, int member);
351
352 /// Initialise @a member in PDF set @a setid.
353 /// @deprecated Use the proper C++ interface of LHAPDF6 instead!
354 void initPDFSet(const std::string& name, int member=0);
355 void initPDFSet(int nset ,const std::string& name, int nmem=0);
356 void initPDFSet(int setid, int member=0);
357 void initPDFSet(int nset , int setid, int nmem=0);
358 /// Initialise @a member in PDF set @a name, of type @a type. - LHAPDF5 compatibility
359 void initPDFSet(const std::string& name, SetType type, int member=0);
360 /// Initialise @a member in PDF set @a name, of type @a type (multi-set version) - LHAPDF5 compatibility
361 void initPDFSet(int nset, const std::string& name, SetType type, int member=0);
362
363}
364
365#endif
366#endif