API Documentation
-
namespace neopdf
Object Oriented interface to NeoPDF.
-
struct PhysicsParameters
- #include <NeoPDF.hpp>
C++ representation of PhysicsParameters.
Public Functions
-
inline NeoPDFPhysicsParameters to_c() const
-
inline NeoPDFPhysicsParameters to_c() const
-
struct MetaData
- #include <NeoPDF.hpp>
C++ representation of NeoPDFMetaData.
Subclassed by neopdf::MetaDataV2
Public Functions
-
inline NeoPDFMetaData to_c() const
Public Members
-
std::string set_desc
-
uint32_t set_index
-
uint32_t num_members
-
double x_min
-
double x_max
-
double q_min
-
double q_max
-
std::vector<int32_t> flavors
-
std::string format
-
std::vector<double> alphas_q_values
-
std::vector<double> alphas_vals
-
bool polarised
-
neopdf_set_type set_type
-
neopdf_interpolator_type interpolator_type
-
std::string error_type
-
int32_t hadron_pid
-
PhysicsParameters phys_params
-
inline NeoPDFMetaData to_c() const
-
struct MetaDataV2 : public neopdf::MetaData
- #include <NeoPDF.hpp>
C++ representation of NeoPDFMetaDataV2.
Public Functions
-
inline NeoPDFMetaDataV2 to_c_v2() const
-
inline NeoPDFMetaDataV2 to_c_v2() const
-
class NeoPDF
- #include <NeoPDF.hpp>
Base PDF class that instantiates the PDF object.
Public Functions
-
inline virtual ~NeoPDF()
Destructor.
-
inline NeoPDF(const std::string &pdf_name, size_t member = 0)
Constructor of the PDF object.
pdf_nameName of the PDF set.memberID number of the PDF member.
-
inline double x_min() const
Get the minimum value of the x-grid for the PDF.
-
inline double x_max() const
Get the maximum value of the x-grid for the PDF.
-
inline double q2_min() const
Get the minimum value of the Q2-grid for the PDF.
-
inline double q2_max() const
Get the maximum value of the Q2-grid for the PDF.
-
inline double xfxQ2(int pid, double x, double q2) const
Compute the
xfvalue for a given PID, x, and Q2.
-
inline double xfxQ2_ND(int pid, std::vector<double> params) const
Compute the
xfvalue for a generic set of parameters.
-
inline std::vector<double> xfxQ2_cheby_batch(int pid, const std::vector<std::vector<double>> &points) const
Compute the
xfvalue for a generic set of parameters using batch Chebyshev interpolation.
-
inline std::vector<double> xfxQ2_pids(const std::vector<int32_t> &pids, const std::vector<double> &points) const
Evaluate all requested flavors at a single kinematic point.
-
inline std::vector<double> xfxQ2s(const std::vector<int32_t> &pids, const std::vector<std::vector<double>> &points) const
Interpolate PDF values for multiple PIDs at multiple kinematic points.
- Returns:
Flat row-major vector of shape [pids.size(), points.size()].
-
inline double alphasQ2(double q2) const
Compute the value of
alphasat the Q2 value.
-
inline size_t num_pids() const
Get the number of PIDs.
-
inline std::vector<int32_t> pids() const
Get the PID representation of the PDF Grid.
-
inline size_t num_subgrids() const
Get the number of subgrids in the PDF Grid.
-
inline std::vector<double> param_range(NeopdfSubgridParams param) const
Get the minimum and maximum value for a given parameter.
-
inline std::vector<size_t> subgrids_shape_for_param(NeopdfSubgridParams param) const
Get the shape of the subgrids in the order of their index for a given parameter.
-
inline std::vector<double> subgrid_for_param(NeopdfSubgridParams param, size_t subgrid_index) const
Get the grid values of a parameter for a given subgrid.
-
inline void set_force_positive(neopdf_force_positive option)
Clip the interpolated values if they turned out negatives.
-
inline neopdf_force_positive is_force_positive() const
Returns the value of
ForcePositivedefining the PDF grid.
Public Static Functions
Protected Functions
-
inline NeoPDF(NeoPDFWrapper *pdf)
Constructor (protected to avoid direct instantiation).
-
NeoPDF() = delete
Deleted copy/move semantics.
Private Members
-
NeoPDFWrapper *raw
Underlying raw object.
Friends
- friend class NeoPDFs
-
inline virtual ~NeoPDF()
-
class NeoPDFs
- #include <NeoPDF.hpp>
Class to load and manage multiple PDF members.
Public Functions
-
inline NeoPDFs(const std::string &pdf_name)
Constructor that loads all PDF members for a given PDF set.
- Parameters:
pdf_name – Name of the PDF set.
-
inline size_t size() const
Get the number of loaded PDF members.
-
inline const NeoPDF &operator[](size_t index) const
Access a specific PDF member by index (const version).
-
inline const NeoPDF &at(size_t index) const
Access a specific PDF member by index with bounds checking (const version).
-
inline void set_force_positive_members(neopdf_force_positive option)
Clip the interpolated values if they turned out negatives for all members.
-
inline NeoPDFs(const std::string &pdf_name)
-
class NeoPDFLazy
- #include <NeoPDF.hpp>
Class for lazily loading PDF members from a .neopdf.lz4 file.
Public Functions
-
inline explicit NeoPDFLazy(const std::string &pdf_name)
Constructor that initializes the lazy iterator for a given PDF set.
- Parameters:
pdf_name – Name of the PDF set (must be a .neopdf.lz4 file).
- Throws:
std::runtime_error – if the iterator cannot be created.
-
inline ~NeoPDFLazy()
Destructor.
-
inline NeoPDFLazy(NeoPDFLazy &&other) noexcept
Move constructor.
-
inline NeoPDFLazy &operator=(NeoPDFLazy &&other) noexcept
Move assignment operator.
-
NeoPDFLazy(const NeoPDFLazy&) = delete
Deleted copy semantics.
-
NeoPDFLazy &operator=(const NeoPDFLazy&) = delete
Private Members
-
::NeoPDFLazyIterator *raw_iter
-
inline explicit NeoPDFLazy(const std::string &pdf_name)
-
class GridWriter
- #include <NeoPDF.hpp>
Class for writing NeoPDF grid data to a file.
Public Functions
-
inline GridWriter()
Constructor.
-
inline ~GridWriter()
Destructor.
-
inline void new_grid()
Starts a new grid for a new member.
-
inline void add_subgrid(const std::vector<double> &nucleons, const std::vector<double> &alphas, const std::vector<double> &kts, const std::vector<double> &xs, const std::vector<double> &q2s, const std::vector<double> &grid_data)
Adds a subgrid to the current grid.
- Parameters:
nucleons – Vector of nucleon numbers.
alphas – Vector of alpha_s values.
kts – Vector of kt values.
xs – Vector of x values.
q2s – Vector of Q2 values.
grid_data – Vector of grid data.
-
inline void add_subgrid_v2(const std::vector<double> &nucleons, const std::vector<double> &alphas, const std::vector<double> &xis, const std::vector<double> &deltas, const std::vector<double> &kts, const std::vector<double> &xs, const std::vector<double> &q2s, const std::vector<double> &grid_data)
Adds a subgrid to the current grid (v2 for 8D).
- Parameters:
nucleons – Vector of nucleon numbers.
alphas – Vector of alpha_s values.
xis – Vector of xi values.
deltas – Vector of delta values.
kts – Vector of kt values.
xs – Vector of x values.
q2s – Vector of Q2 values.
grid_data – Vector of grid data.
-
inline void push_grid(const std::vector<int32_t> &flavors)
Finalizes the current grid, sets its flavors, and adds it to the collection.
- Parameters:
flavors – Vector of flavor IDs.
-
inline void compress(const MetaData &metadata, const std::string &output_path)
Compresses the added grids and writes them to a file.
- Parameters:
metadata – The metadata for the PDF set.
output_path – The path to the output file.
-
inline void compress_v2(const MetaDataV2 &metadata, const std::string &output_path)
Compresses the added grids and writes them to a file using V2 metadata.
- Parameters:
metadata – The V2 metadata for the PDF set.
output_path – The path to the output file.
-
inline GridWriter()
-
struct PhysicsParameters
-
namespace NEOLHAPDF
LHAPDF compatibility for no-code migration.
Functions
-
inline void setVerbosity(int)
-
class PDF
- #include <NeoPDF.hpp>
Subclassed by NEOLHAPDF::GridPDF
-
inline void setVerbosity(int)