neopdf 0.2.0-alpha7

A modern, fast, and reliable PDF interpolation library
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
//! This module provides the high-level interface for working with PDF sets.
//!
//! It defines the [`PDF`] struct, which serves as the main entry point for accessing,
//! interpolating, and retrieving metadata from PDF sets. The module abstracts over different
//! PDF set formats (LHAPDF and NeoPDF) and provides convenient loader functions for both
//! single and multiple PDF members.
//!
//! # Main Features
//!
//! - Unified interface for loading and accessing PDF sets from different formats.
//! - Parallel loading of all PDF members for efficient batch operations.
//! - High-level interpolation methods for PDF values and strong coupling constant (`alpha_s`).
//! - Access to underlying grid data and metadata for advanced use cases.
//!
//! # Key Types
//!
//! - [`PDF`]: Represents a single PDF member, providing methods for interpolation and metadata access.
//! - [`PdfSet`]: Trait for abstracting over different PDF set backends.
//! - Loader functions: [`PDF::load`], [`PDF::load_pdfs`], and internal helpers for batch loading.
//!
//! See the documentation for [`PDF`] for more details on available methods and usage patterns.
use ndarray::{Array1, Array2};
use rayon::prelude::*;

use super::gridpdf::{ForcePositive, GridArray, GridPDF};
use super::metadata::MetaData;
use super::parser::{LhapdfSet, NeopdfSet};
use super::subgrid::{RangeParameters, SubGrid};

/// Trait for abstracting over different PDF set backends (e.g., LHAPDF, NeoPDF).
///
/// Provides a unified interface for accessing the number of members and retrieving individual
/// members as metadata and grid arrays.
trait PdfSet: Send + Sync {
    /// Returns the number of members in the PDF set.
    fn num_members(&self) -> usize;
    /// Retrieves the metadata and grid array for the specified member index.
    fn member(&self, idx: usize) -> (MetaData, GridArray);
}

impl PdfSet for LhapdfSet {
    fn num_members(&self) -> usize {
        self.info.num_members as usize
    }
    fn member(&self, idx: usize) -> (MetaData, GridArray) {
        self.member(idx)
    }
}

impl PdfSet for NeopdfSet {
    fn num_members(&self) -> usize {
        self.info.num_members as usize
    }
    fn member(&self, idx: usize) -> (MetaData, GridArray) {
        self.member(idx)
    }
}

/// Loads a single PDF member from a generic PDF set backend.
///
/// # Arguments
///
/// * `set` - The PDF set backend implementing [`PdfSet`].
/// * `member` - The index of the member to load.
///
/// # Returns
///
/// A [`PDF`] instance for the specified member.
fn pdfset_loader<T: PdfSet>(set: T, member: usize) -> PDF {
    let (info, knot_array) = set.member(member);
    PDF {
        grid_pdf: GridPDF::new(info, knot_array),
    }
}

/// Loads all PDF members from a generic PDF set backend in sequential.
///
/// # Arguments
///
/// * `set` - The PDF set backend implementing [`PdfSet`].
///
/// # Returns
///
/// A vector of [`PDF`] instances, one for each member in the set.
fn pdfsets_seq_loader<T: PdfSet + Send + Sync>(set: T) -> Vec<PDF> {
    (0..set.num_members())
        .map(|idx| {
            let (info, knot_array) = set.member(idx);
            PDF {
                grid_pdf: GridPDF::new(info, knot_array),
            }
        })
        .collect()
}

/// Loads all PDF members from a generic PDF set backend in parallel.
///
/// # Arguments
///
/// * `set` - The PDF set backend implementing [`PdfSet`].
///
/// # Returns
///
/// A vector of [`PDF`] instances, one for each member in the set.
fn pdfsets_par_loader<T: PdfSet + Send + Sync>(set: T) -> Vec<PDF> {
    (0..set.num_members())
        .into_par_iter()
        .map(|idx| {
            let (info, knot_array) = set.member(idx);
            PDF {
                grid_pdf: GridPDF::new(info, knot_array),
            }
        })
        .collect()
}

/// Represents a Parton Distribution Function (PDF) set.
///
/// This struct provides a high-level interface for accessing PDF data,
/// including interpolation and metadata retrieval. It encapsulates the
/// `GridPDF` struct, which handles the low-level grid operations.
pub struct PDF {
    grid_pdf: GridPDF,
}

impl PDF {
    /// Loads a given member of the PDF set.
    ///
    /// This function reads the `.info` file and the corresponding `.dat` member file
    /// to construct a `GridPDF` object, which is then wrapped in a `PDF` instance.
    ///
    /// # Arguments
    ///
    /// * `pdf_name` - The name of the PDF set (e.g., "NNPDF40_nnlo_as_01180").
    /// * `member` - The ID of the PDF member to load (0-indexed).
    ///
    /// # Returns
    ///
    /// A `PDF` instance representing the loaded PDF member.
    pub fn load(pdf_name: &str, member: usize) -> Self {
        if pdf_name.ends_with(".neopdf.lz4") {
            pdfset_loader(NeopdfSet::new(pdf_name), member)
        } else {
            pdfset_loader(LhapdfSet::new(pdf_name), member)
        }
    }

    /// Loads all members of a PDF set in parallel.
    ///
    /// This function reads the `.info` file and all `.dat` member files
    /// to construct a `Vec<PDF>`, with each `PDF` instance representing a member
    /// of the set. The loading is performed in parallel.
    ///
    /// # Arguments
    ///
    /// * `pdf_name` - The name of the PDF set.
    ///
    /// # Returns
    ///
    /// A `Vec<PDF>` where each element is a `PDF` instance for a member of the set.
    pub fn load_pdfs(pdf_name: &str) -> Vec<PDF> {
        if pdf_name.ends_with(".neopdf.lz4") {
            pdfsets_par_loader(NeopdfSet::new(pdf_name))
        } else {
            pdfsets_par_loader(LhapdfSet::new(pdf_name))
        }
    }

    /// Loads all members of a PDF set in sequential.
    ///
    /// This function reads the `.info` file and all `.dat` member files
    /// to construct a `Vec<PDF>`, with each `PDF` instance representing a member
    /// of the set. The loading is performed in parallel.
    ///
    /// # Arguments
    ///
    /// * `pdf_name` - The name of the PDF set.
    ///
    /// # Returns
    ///
    /// A `Vec<PDF>` where each element is a `PDF` instance for a member of the set.
    pub fn load_pdfs_seq(pdf_name: &str) -> Vec<PDF> {
        if pdf_name.ends_with(".neopdf.lz4") {
            pdfsets_seq_loader(NeopdfSet::new(pdf_name))
        } else {
            pdfsets_seq_loader(LhapdfSet::new(pdf_name))
        }
    }

    /// Creates an iterator that loads PDF members lazily.
    ///
    /// This function is suitable for `.neopdf.lz4` files, which support lazy loading.
    /// It returns an iterator that yields `PDF` instances on demand, which is useful
    /// for reducing memory consumption when working with large PDF sets.
    ///
    /// # Arguments
    ///
    /// * `pdf_name` - The name of the PDF set (must end with `.neopdf.lz4`).
    ///
    /// # Returns
    ///
    /// An iterator over `Result<PDF, Box<dyn std::error::Error>>`.
    pub fn load_pdfs_lazy(
        pdf_name: &str,
    ) -> impl Iterator<Item = Result<PDF, Box<dyn std::error::Error>>> {
        assert!(
            pdf_name.ends_with(".neopdf.lz4"),
            "Lazy loading is only supported for .neopdf.lz4 files"
        );

        let iter_lazy = NeopdfSet::new(pdf_name).into_lazy_iterators();

        iter_lazy.map(|grid_array_with_metadata_result| {
            grid_array_with_metadata_result.map(|grid_array_with_metadata| {
                let info = (*grid_array_with_metadata.metadata).clone();
                let knot_array = grid_array_with_metadata.grid;
                PDF {
                    grid_pdf: GridPDF::new(info, knot_array),
                }
            })
        })
    }

    /// Clip the negative values for the `PDF` object.
    ///
    /// # Arguments
    ///
    /// * `option` - The method used to clip negative values.
    pub fn set_force_positive(&mut self, option: ForcePositive) {
        self.grid_pdf.set_force_positive(option);
    }

    /// Clip the negative values for all the `PDF` objects.
    ///
    /// # Arguments
    ///
    /// * `pdfs` - A `Vec<PDF>` where each element is a `PDF` instance.
    /// * `option` - The method used to clip negative values.
    pub fn set_force_positive_members(pdfs: &mut [PDF], option: ForcePositive) {
        for pdf in pdfs {
            pdf.set_force_positive(option.clone());
        }
    }

    /// Returns the clipping method used for a single `PDF` object.
    ///
    /// # Returns
    ///
    /// The clipping method given as a `ForcePositive` object.
    pub fn is_force_positive(&self) -> &ForcePositive {
        self.grid_pdf
            .force_positive
            .as_ref()
            .unwrap_or(&ForcePositive::NoClipping)
    }

    /// Interpolates the PDF value (xf) for a given nucleon, alphas, flavor, x, and Q2.
    ///
    /// Abstraction to the `GridPDF::xfxq2` method.
    ///
    /// # Arguments
    ///
    /// * `id` - The flavor ID (PDG ID).
    /// * `points` - A slice containing the collection of points to interpolate on.
    ///
    /// # Returns
    ///
    /// The interpolated PDF value `xf(nuclone, alphas, flavor, x, Q^2)`.
    pub fn xfxq2(&self, pid: i32, points: &[f64]) -> f64 {
        self.grid_pdf.xfxq2(pid, points).unwrap()
    }

    /// Interpolates the PDF value (xf) for multiple nucleons, alphas, flavors, xs, and Q2s.
    ///
    /// Abstraction to the `GridPDF::xfxq2s` method.
    ///
    /// # Arguments
    ///
    /// * `ids` - A vector of flavor IDs.
    /// * `slice_points` - A slice containing the collection of knots to interpolate on.
    ///   A knot is a collection of points containing `(nucleon, alphas, x, Q2)`.
    ///
    /// # Returns
    ///
    /// A 2D array of interpolated PDF values with shape `[flavors, N_knots]`.
    pub fn xfxq2s(&self, pids: Vec<i32>, slice_points: &[&[f64]]) -> Array2<f64> {
        self.grid_pdf.xfxq2s(pids, slice_points)
    }

    /// Interpolates the PDF value (xf) for multiple points using Chebyshev batch interpolation.
    ///
    /// Abstraction to the `GridPDF::xfxq2_cheby_batch` method.
    ///
    /// # Arguments
    ///
    /// * `pid` - The flavor ID.
    /// * `points` - A slice containing the collection of knots to interpolate on.
    ///   A knot is a collection of points containing `(nucleon, alphas, x, Q2)`.
    ///
    /// # Returns
    ///
    /// A `Vec<f64>` of interpolated PDF values.
    pub fn xfxq2_cheby_batch(&self, pid: i32, points: &[&[f64]]) -> Vec<f64> {
        self.grid_pdf.xfxq2_cheby_batch(pid, points).unwrap()
    }

    /// Interpolates the strong coupling constant `alpha_s` for a given Q2.
    ///
    /// Abstraction to the `GridPDF::alphas_q2` method.
    ///
    /// # Arguments
    ///
    /// * `q2` - The squared energy scale.
    ///
    /// # Returns
    ///
    /// The interpolated `alpha_s` value.
    pub fn alphas_q2(&self, q2: f64) -> f64 {
        self.grid_pdf.alphas_q2(q2)
    }

    /// Returns a reference to the PDF metadata.
    ///
    /// Abstraction to the `GridPDF::info` method.
    ///
    /// # Returns
    ///
    /// A `MetaData` struct containing information about the PDF set.
    pub fn metadata(&self) -> &MetaData {
        self.grid_pdf.metadata()
    }

    /// Returns the number of subgrids in the PDF set.
    ///
    /// # Returns
    ///
    /// The number of subgrids.
    pub fn num_subgrids(&self) -> usize {
        self.grid_pdf.knot_array.subgrids.len()
    }

    /// Returns a reference to the subgrid at the given index.
    ///
    /// # Arguments
    ///
    /// * `index` - The index of the subgrid.
    ///
    /// # Returns
    ///
    /// A reference to the `SubGrid`.
    pub fn subgrid(&self, index: usize) -> &SubGrid {
        &self.grid_pdf.knot_array.subgrids[index]
    }

    /// Returns references to all the subgrid at the given index.
    ///
    /// # Returns
    ///
    /// A reference to all the `SubGrid`.
    pub fn subgrids(&self) -> &Vec<SubGrid> {
        &self.grid_pdf.knot_array.subgrids
    }

    /// Returns the flavor PIDS of the PDG Grid.
    ///
    /// # Returns
    ///
    /// PID representation of the PDF.
    pub fn pids(&self) -> &Array1<i32> {
        &self.grid_pdf.knot_array.pids
    }

    /// Retrieves the ranges for the parameters.
    ///
    /// Abstraction to the `GridPDF::param_ranges` method.
    ///
    /// # Returns
    ///
    /// The minimum and maximum values for the parameters (x, q2, ...).
    pub fn param_ranges(&self) -> RangeParameters {
        self.grid_pdf.param_ranges()
    }

    /// Retrieves the PDF value (xf) at a specific knot point in the grid.
    ///
    /// Abstraction to the `GridArray::xf_from_index` method. This method does not
    /// perform any interpolation.
    ///
    /// # Arguments
    ///
    /// * `i_nucleons` - The index of the nucleon.
    /// * `i_alphas` - The index of the alpha_s value.
    /// * `i_kt` - The index of the `kT` value.
    /// * `ix` - The index of the x-value.
    /// * `iq2` - The index of the Q2-value.
    /// * `id` - The flavor ID.
    /// * `subgrid_id` - The ID of the subgrid.
    ///
    /// # Returns
    ///
    /// The PDF value at the specified knot.
    #[allow(clippy::too_many_arguments)]
    pub fn xf_from_index(
        &self,
        i_nucleons: usize,
        i_alphas: usize,
        i_kt: usize,
        ix: usize,
        iq2: usize,
        id: i32,
        subgrid_id: usize,
    ) -> f64 {
        self.grid_pdf
            .knot_array
            .xf_from_index(i_nucleons, i_alphas, i_kt, ix, iq2, id, subgrid_id)
    }
}