PDFxTMDLib  1.0.0
Public Types | Public Member Functions | List of all members
PDFxTMD::PDFSet< Tag > Class Template Reference

Manages a set of Parton Distribution Functions (PDFs), providing tools for uncertainty and correlation analysis. More...

#include <PDFSet.h>

Public Types

using PDF_t = typename PDFType< Tag >::type
 The specific PDF interface type (ITMD or ICPDF) determined by the Tag. More...
 

Public Member Functions

 PDFSet (std::string pdfSetName, bool alternativeReplicaUncertainty=false)
 Constructs a PDFSet for a given PDF name. More...
 
 PDFSet ()
 Default constructor. More...
 
double alphasQ (double q) const
 Get the strong coupling constant alpha_s at a given scale Q. More...
 
double alphasQ2 (double q2) const
 Get the strong coupling constant alpha_s at a given squared scale Q^2. More...
 
PDF_toperator[] (int member)
 Access a specific PDF member from the set. More...
 
PDF_toperator[] (int member) const
 Access a specific PDF member from the set (const version). More...
 
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, TMDPDFTag>>>
void Uncertainty (PartonFlavor flavor, double x, double kt2, double mu2, double cl, PDFUncertainty &resUncertainty)
 Calculate the TMD uncertainty. More...
 
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, CollinearPDFTag>>>
void Uncertainty (PartonFlavor flavor, double x, double mu2, double cl, PDFUncertainty &resUncertainty)
 Calculate the collinear PDF uncertainty. More...
 
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, TMDPDFTag>>>
PDFUncertainty Uncertainty (PartonFlavor flavor, double x, double kt2, double mu2, double cl=NO_REQUESTED_CONFIDENCE_LEVEL)
 Calculate the TMD uncertainty and return the result. More...
 
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, CollinearPDFTag>>>
PDFUncertainty Uncertainty (PartonFlavor flavor, double x, double mu2, double cl=NO_REQUESTED_CONFIDENCE_LEVEL)
 Calculate the collinear PDF uncertainty and return the result. More...
 
void Uncertainty (const std::vector< double > &values, double cl, PDFUncertainty &resUncertainty)
 Calculate uncertainty from a pre-computed vector of PDF values. More...
 
PDFUncertainty Uncertainty (const std::vector< double > &values, double cl=NO_REQUESTED_CONFIDENCE_LEVEL)
 Calculate uncertainty from a pre-computed vector of PDF values and return the result. More...
 
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, CollinearPDFTag>>>
double Correlation (PartonFlavor flavorA, double xA, double mu2A, PartonFlavor flavorB, double xB, double mu2B)
 Calculate the correlation for collinear PDFs. More...
 
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, TMDPDFTag>>>
double Correlation (PartonFlavor flavorA, double xA, double kt2A, double mu2A, PartonFlavor flavorB, double xB, double kt2B, double mu2B)
 Calculate the correlation for TMDs. More...
 
double Correlation (const std::vector< double > &valuesA, const std::vector< double > &valuesB) const
 Calculate correlation from two pre-computed vectors of PDF values. More...
 
void CreatePDFSet (unsigned int setMember)
 Explicitly creates a single PDF member. More...
 
void CreateAllPDFSets ()
 Pre-loads all PDF members in the set. More...
 
void InitailizePDFSetName (std::string pdfSetName)
 Re-initializes the PDFSet with a new PDF set name. More...
 
size_t size () const
 Get the total number of members in this PDF set. More...
 
YamlStandardTMDInfo getStdPDFInfo () const
 Get the standard metadata for the PDF set. More...
 
YamlErrorInfo getPDFErrorInfo () const
 Get the error metadata for the PDF set. More...
 
ConfigWrapper info ()
 Get the full configuration metadata object for the set. More...
 
 PDFSet (const PDFSet &)=delete
 
PDFSetoperator= (const PDFSet &)=delete
 
 PDFSet (PDFSet &&)=default
 
PDFSetoperator= (PDFSet &&)=default
 

Detailed Description

template<typename Tag>
class PDFxTMD::PDFSet< Tag >

Manages a set of Parton Distribution Functions (PDFs), providing tools for uncertainty and correlation analysis.

This class acts as a container for all members of a specific PDF set (e.g., CT18, PB-NLO). It handles the loading of PDF data, calculation of uncertainties using various methods (replicas, Hessian), and evaluation of correlations between PDF values. It is a templated class that can be specialized for either TMD or Collinear PDFs.

Template Parameters
TagA tag to specify the PDF type (TMDPDFTag or CollinearPDFTag).

Member Typedef Documentation

◆ PDF_t

template<typename Tag >
using PDFxTMD::PDFSet< Tag >::PDF_t = typename PDFType<Tag>::type

The specific PDF interface type (ITMD or ICPDF) determined by the Tag.

Constructor & Destructor Documentation

◆ PDFSet() [1/4]

template<typename Tag >
PDFxTMD::PDFSet< Tag >::PDFSet ( std::string  pdfSetName,
bool  alternativeReplicaUncertainty = false 
)
inlineexplicit

Constructs a PDFSet for a given PDF name.

Parameters
pdfSetNameThe name of the PDF set to load.
alternativeReplicaUncertaintyIf true, use percentile strategy for replica uncertainties; otherwise, use standard deviation.

◆ PDFSet() [2/4]

template<typename Tag >
PDFxTMD::PDFSet< Tag >::PDFSet ( )
inline

Default constructor.

◆ PDFSet() [3/4]

template<typename Tag >
PDFxTMD::PDFSet< Tag >::PDFSet ( const PDFSet< Tag > &  )
delete

◆ PDFSet() [4/4]

template<typename Tag >
PDFxTMD::PDFSet< Tag >::PDFSet ( PDFSet< Tag > &&  )
default

Member Function Documentation

◆ alphasQ()

template<typename Tag >
double PDFxTMD::PDFSet< Tag >::alphasQ ( double  q) const
inline

Get the strong coupling constant alpha_s at a given scale Q.

Parameters
qThe momentum transfer scale Q in GeV.
Returns
The value of alpha_s(Q).

◆ alphasQ2()

template<typename Tag >
double PDFxTMD::PDFSet< Tag >::alphasQ2 ( double  q2) const
inline

Get the strong coupling constant alpha_s at a given squared scale Q^2.

Parameters
q2The squared momentum transfer scale Q^2 in GeV^2.
Returns
The value of alpha_s(Q^2).

◆ Correlation() [1/3]

template<typename Tag >
double PDFxTMD::PDFSet< Tag >::Correlation ( const std::vector< double > &  valuesA,
const std::vector< double > &  valuesB 
) const
inline

Calculate correlation from two pre-computed vectors of PDF values.

Parameters
valuesAA vector of PDF values for the first observable.
valuesBA vector of PDF values for the second observable.
Returns
The correlation coefficient.

◆ Correlation() [2/3]

template<typename Tag >
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, TMDPDFTag>>>
double PDFxTMD::PDFSet< Tag >::Correlation ( PartonFlavor  flavorA,
double  xA,
double  kt2A,
double  mu2A,
PartonFlavor  flavorB,
double  xB,
double  kt2B,
double  mu2B 
)
inline

Calculate the correlation for TMDs.

(Enabled only for TMDPDFTag)

Parameters
flavorAThe parton flavor for the first TMD.
xAThe momentum fraction 'x' for the first TMD.
kt2AThe transverse momentum 'kt2' for the first TMD.
mu2AThe scale 'mu2' for the first TMD.
flavorBThe parton flavor for the second TMD.
xBThe momentum fraction 'x' for the second TMD.
kt2BThe transverse momentum 'kt2' for the second TMD.
mu2BThe scale 'mu2' for the second TMD.
Returns
The correlation coefficient.

◆ Correlation() [3/3]

template<typename Tag >
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, CollinearPDFTag>>>
double PDFxTMD::PDFSet< Tag >::Correlation ( PartonFlavor  flavorA,
double  xA,
double  mu2A,
PartonFlavor  flavorB,
double  xB,
double  mu2B 
)
inline

Calculate the correlation for collinear PDFs.

(Enabled only for CollinearPDFTag)

Parameters
flavorAThe parton flavor for the first observable.
xAThe momentum fraction 'x' for the first observable.
mu2AThe scale 'mu2' for the first observable.
flavorBThe parton flavor for the second observable.
xBThe momentum fraction 'x' for the second observable.
mu2BThe scale 'mu2' for the second observable.
Returns
The correlation coefficient.

◆ CreateAllPDFSets()

template<typename Tag >
void PDFxTMD::PDFSet< Tag >::CreateAllPDFSets ( )
inline

Pre-loads all PDF members in the set.

◆ CreatePDFSet()

template<typename Tag >
void PDFxTMD::PDFSet< Tag >::CreatePDFSet ( unsigned int  setMember)
inline

Explicitly creates a single PDF member.

Parameters
setMemberThe index of the PDF member to create.

◆ getPDFErrorInfo()

template<typename Tag >
YamlErrorInfo PDFxTMD::PDFSet< Tag >::getPDFErrorInfo ( ) const
inline

Get the error metadata for the PDF set.

Returns
A YamlErrorInfo object.

◆ getStdPDFInfo()

template<typename Tag >
YamlStandardTMDInfo PDFxTMD::PDFSet< Tag >::getStdPDFInfo ( ) const
inline

Get the standard metadata for the PDF set.

Returns
A YamlStandardPDFInfo object.

◆ info()

template<typename Tag >
ConfigWrapper PDFxTMD::PDFSet< Tag >::info ( )
inline

Get the full configuration metadata object for the set.

Returns
A ConfigWrapper object.

◆ InitailizePDFSetName()

template<typename Tag >
void PDFxTMD::PDFSet< Tag >::InitailizePDFSetName ( std::string  pdfSetName)
inline

Re-initializes the PDFSet with a new PDF set name.

Parameters
pdfSetNameThe name of the new PDF set.

◆ operator=() [1/2]

template<typename Tag >
PDFSet& PDFxTMD::PDFSet< Tag >::operator= ( const PDFSet< Tag > &  )
delete

◆ operator=() [2/2]

template<typename Tag >
PDFSet& PDFxTMD::PDFSet< Tag >::operator= ( PDFSet< Tag > &&  )
default

◆ operator[]() [1/2]

template<typename Tag >
PDF_t* PDFxTMD::PDFSet< Tag >::operator[] ( int  member)
inline

Access a specific PDF member from the set.

Parameters
memberThe member index (0 is the central value).
Returns
A pointer to the PDF object. Creates the object if it doesn't exist.

◆ operator[]() [2/2]

template<typename Tag >
PDF_t* PDFxTMD::PDFSet< Tag >::operator[] ( int  member) const
inline

Access a specific PDF member from the set (const version).

Parameters
memberThe member index.
Returns
A const pointer to the PDF object, or nullptr if it has not been created.

◆ size()

template<typename Tag >
size_t PDFxTMD::PDFSet< Tag >::size ( ) const
inline

Get the total number of members in this PDF set.

Returns
The number of members.

◆ Uncertainty() [1/6]

template<typename Tag >
void PDFxTMD::PDFSet< Tag >::Uncertainty ( const std::vector< double > &  values,
double  cl,
PDFUncertainty resUncertainty 
)
inline

Calculate uncertainty from a pre-computed vector of PDF values.

Parameters
valuesA vector of PDF values from all members of the set.
clThe desired confidence level in percent.
resUncertaintyThe output PDFUncertainty object.

◆ Uncertainty() [2/6]

template<typename Tag >
PDFUncertainty PDFxTMD::PDFSet< Tag >::Uncertainty ( const std::vector< double > &  values,
double  cl = NO_REQUESTED_CONFIDENCE_LEVEL 
)
inline

Calculate uncertainty from a pre-computed vector of PDF values and return the result.

Parameters
valuesA vector of PDF values from all members of the set.
clThe desired confidence level in percent (default is the set's native CL).
Returns
A PDFUncertainty object.

◆ Uncertainty() [3/6]

template<typename Tag >
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, TMDPDFTag>>>
void PDFxTMD::PDFSet< Tag >::Uncertainty ( PartonFlavor  flavor,
double  x,
double  kt2,
double  mu2,
double  cl,
PDFUncertainty resUncertainty 
)
inline

Calculate the TMD uncertainty.

(Enabled only for TMDPDFTag)

Parameters
flavorThe parton flavor.
xThe momentum fraction.
kt2The squared transverse momentum.
mu2The squared factorization scale.
clThe desired confidence level in percent.
resUncertaintyThe output PDFUncertainty object.

◆ Uncertainty() [4/6]

template<typename Tag >
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, TMDPDFTag>>>
PDFUncertainty PDFxTMD::PDFSet< Tag >::Uncertainty ( PartonFlavor  flavor,
double  x,
double  kt2,
double  mu2,
double  cl = NO_REQUESTED_CONFIDENCE_LEVEL 
)
inline

Calculate the TMD uncertainty and return the result.

(Enabled only for TMDPDFTag)

Parameters
flavorThe parton flavor.
xThe momentum fraction.
kt2The squared transverse momentum.
mu2The squared factorization scale.
clThe desired confidence level in percent (default is the set's native CL).
Returns
A PDFUncertainty object with the results.

◆ Uncertainty() [5/6]

template<typename Tag >
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, CollinearPDFTag>>>
void PDFxTMD::PDFSet< Tag >::Uncertainty ( PartonFlavor  flavor,
double  x,
double  mu2,
double  cl,
PDFUncertainty resUncertainty 
)
inline

Calculate the collinear PDF uncertainty.

(Enabled only for CollinearPDFTag)

Parameters
flavorThe parton flavor.
xThe momentum fraction.
mu2The squared factorization scale.
clThe desired confidence level in percent.
resUncertaintyThe output PDFUncertainty object.

◆ Uncertainty() [6/6]

template<typename Tag >
template<typename T = Tag, typename = std::enable_if_t<std::is_same_v<T, CollinearPDFTag>>>
PDFUncertainty PDFxTMD::PDFSet< Tag >::Uncertainty ( PartonFlavor  flavor,
double  x,
double  mu2,
double  cl = NO_REQUESTED_CONFIDENCE_LEVEL 
)
inline

Calculate the collinear PDF uncertainty and return the result.

(Enabled only for CollinearPDFTag)

Parameters
flavorThe parton flavor.
xThe momentum fraction.
mu2The squared factorization scale.
clThe desired confidence level in percent (default is the set's native CL).
Returns
A PDFUncertainty object with the results.

The documentation for this class was generated from the following file: