Python API
In the following example, we are going to showcase some of the main features of the NeoPDF
Python API. The outline of the tutorial is as follow:
1. Standard PDF Manipulation ¶
import numpy as np
from neopdf.pdf import PDF as NEOPDF
# To load NeoPDF sets, if available, append the name with `.neopdf.lz4`
pdf = NEOPDF.mkPDF("NNPDF40_nnlo_as_01180")
# Get the Metadata of the PDF set.
metadata = pdf.metadata()
metadata_dict = metadata.to_dict()
metadata_dict.keys()
dict_keys(['set_desc', 'set_index', 'num_members', 'x_min', 'x_max', 'q_min', 'q_max', 'flavors', 'format', 'alphas_q_values', 'alphas_vals', 'polarised', 'set_type', 'interpolator_type', 'error_type', 'hadron_pid', 'git_version', 'code_version', 'flavor_scheme', 'order_qcd', 'alphas_order_qcd', 'm_w', 'm_z', 'm_up', 'm_down', 'm_strange', 'm_charm', 'm_bottom', 'm_top'])
metadata_dict["set_desc"]
'NNPDF4.0 NNLO global fit, alphas(MZ)=0.1180. mem=0 => average on replicas; mem=1-100 => PDF replicas'
metadata_dict["alphas_vals"][:4]
[0.33074891, 0.3176246, 0.30507081, 0.29305875]
metadata_dict["set_type"]
'PDF'
We can also check the shape of the subgrids. In the following example, we are checking the shape of the first subgrid:
pdf_subgrids = pdf.subgrids()
pdf_subgrids[0].grid_shape()
(1, 1, 11, 1, 196, 12)
The way to compute the value of the interpolated PDF $x f(x, Q^2)$ and the strong coupling $\alpha_s$ is the same as in LHAPDF.
pdf.alphasQ2(q2=10)
0.2485925816007479
pdf.xfxQ2(id=21, x=1e-5, q2=1e3)
111.20984759980468
One can also compute $x f(x, Q^2)$ for values of PID, $x$, and $Q^2$:
xs = np.geomspace(1e-9, 1.0, 50, endpoint=False)
q2s = np.geomspace(1e2, 1e6, 50, endpoint=False)
pids = [nf for nf in range(-4, 5) if nf != 0]
result = pdf.xfxQ2s(pids=pids, xs=xs, q2s=q2s)
result = result.reshape(len(pids), xs.size, q2s.size)
result[0][0][0]
12.628715124538546
NOTE: To load all the PDF members at once, as in LHAPDF, one can simply do:
all_pdfs = NEOPDF.mkPDFs("NNPDF40_nnlo_as_01180")
2. Convert LHAPDF into NeoPDF Format ¶
One can convert PDF sets in the LHAPDF format into the new format used by NeoPDF.
from neopdf.converter import convert_lhapdf
# NeoPDF files should be appended with the suffix `neopdf.lz4`
convert_lhapdf(
pdf_name="NNPDF40_nnlo_as_01180", output_path="NNPDF40_nnlo_as_01180.neopdf.lz4"
)
The above code will generate a file called NNPDF40_nnlo_as_01180.neopdf.lz4
in the current directory. In order to load it, put it in the (default) NEOPDF_DATA_PATH
path. See the installation
section for more details. The set can be loaded in the same way as for standard LHAPDF sets.
neopdf = NEOPDF.mkPDF("NNPDF40_nnlo_as_01180.neopdf.lz4")
pdf.xfxQ2(id=21, x=1e-5, q2=1e3)
111.20984759980468
3. Create a NeoPDF Grid ¶
The following section shows how one can create a NeoPDF PDF grid from a given predictions/fits. As mentioned in the Design
section, a NeoPDF subgrid is a 5-dimensional array in $(A, \alpha_s, x, Q^2)$. In turn, each PDF member can contain multiple subgrids.
# Some imports
from neopdf.writer import compress
from neopdf.gridpdf import GridArray, SubGrid
from neopdf.metadata import InterpolatorType, SetType, MetaData, PhysicsParameters
# Construct the shared Metadata
num_members = 10
x_min = 1e-5
x_max = 1.0
q_min = 4.0
q_max = 1e5
flavors = [nf for nf in range(-4, 5)]
num_alphas = 6
alphas_vals = np.random.uniform(0.1, 0.2, num_alphas)
alphas_qvalues = np.geomspace(q_min, q_max, num_alphas)
# Define the Physical Parameters in the Metadata
physparams_kwargs = {
"flavor_scheme": "fixed",
"order_qcd": 2,
"alphas_order_qcd": 2,
"m_w": 80.3520,
"m_z": 91.1876,
"m_up": 0.0,
"m_down": 0.0,
"m_strange": 0.0,
"m_charm": 1.51,
"m_bottom": 4.92,
"m_top": 172.5,
}
phys_params = PhysicsParameters(**physparams_kwargs)
# Construct the PDF Metadata
metadata_kwargs = {
"set_desc": "Some Toy NeoPDF set",
"set_index": 123456,
"num_members": num_members,
"x_min": x_min,
"x_max": x_max,
"q_min": q_min,
"q_max": q_max,
"flavors": flavors,
"format": "neopdf",
"alphas_q_values": alphas_qvalues,
"alphas_vals": alphas_vals,
"polarised": False,
"set_type": SetType.SpaceLike,
"interpolator_type": InterpolatorType.LogBicubic,
"phys_params": phys_params,
}
metadata = MetaData(**metadata_kwargs)
# Construct the PDF Grid:
# - A subrid is represented by the `SubGrid` object
# - A member is represented by the `GridArray` object
nucleons = [1] # Proton Only
alphas_mZ = [0.118] # Only one value
kts = [0.0] # No dependence on the Transverse Momentum
grid_members = []
x_values = np.geomspace(x_min, x_max, 50)
# Subdivide the Q2 range into subgrids
q_mid = int(q_max / 3)
q2_sub1 = np.geomspace(q_min * q_min, q_mid * q_mid, 25)
q2_sub2 = np.geomspace(q_mid * q_mid, q_max * q_max, 25)
q2_values = [q2_sub1, q2_sub2]
for _ in range(num_members):
sub_grids = []
for q2_vals in q2_values:
grid_shape = (
len(nucleons),
len(alphas_mZ),
len(flavors),
len(kts),
x_values.size,
q2_vals.size,
)
grid = np.random.uniform(0, 1, grid_shape)
sub_grid = SubGrid(
xs=x_values,
q2s=q2_vals,
kts=kts,
nucleons=nucleons,
alphas=alphas_mZ,
grid=grid,
)
sub_grids.append(sub_grid)
grid_member = GridArray(pids=flavors, subgrids=sub_grids)
grid_members.append(grid_member)
# Write the compressed PDF set into disk
compress(grids=grid_members, metadata=metadata, path="TOY_NEOPDF.neopdf.lz4")
This will generate a file called TOY_NEOPDF.neopdf.lz4
in the current directory that one can then move into the (default) NEOPDF_DATA_PATH
. Alternatively, one can pass directly the path to where the the grid will be written in compress
via the path
argument.
4. Combining LHAPDF Nuclear Sets into a Single NeoPDF Grid ¶
NeoPDF
can combine multiple LHAPDF sets along a given parameter into a single NeoPDF
grid such that the parameter-dependence is intrinsically included. In the following example, we are going to combine the nuclear $A$ dependence, however combining multiple $\alpha_s$ variations also works the same.
NOTE: As a default, the bi-cubic log-interpolation is used, and because of how such an interpolation works, the cardinal of the $\left\{A\right\}$ set must at least be four (04).
from neopdf.converter import combine_lhapdf_npdfs
# List of nuclear PDFs to be combined - at least 4
npdfs_list = [
"nNNPDF30_nlo_as_0118_p",
"nNNPDF30_nlo_as_0118_A2_Z1",
"nNNPDF30_nlo_as_0118_A4_Z2",
"nNNPDF30_nlo_as_0118_A6_Z3",
"nNNPDF30_nlo_as_0118_A9_Z4",
]
combine_lhapdf_npdfs(
pdf_names=npdfs_list, output_path="nNNPDF30_nlo_as_0118_A1_A9.neopdf.lz4"
)
We can now load the combined nuclear PDF sets and check that we can now interpolate the $A$-dependence using the xfxQ2_ND
function.
nuclear_neopdfs = NEOPDF.mkPDF("nNNPDF30_nlo_as_0118_A1_A9.neopdf.lz4")
A_value = 5 # No nPDF available
nuclear_neopdfs.xfxQ2_ND(21, [A_value, 1e-5, 1e3])
134.7438112942475
If instead of the $A$ dependence the grid containts the $\alpha_s$ dependence, then one can just replace A_value
with alphas_value
. If the grid has a dependence on both, then syntax becomes:
A_value = 5
alphas_value = 0.117
nuclear_neopdfs.xfxQ2_ND(21, [A_value, alphas_value, 1e-5, 1e3]) # Note the Order