PDFxTMDLib  1.0.0
Public Member Functions | List of all members
PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator > Class Template Reference

#include <GenericPDF.h>

Public Member Functions

 GenericPDF (const std::string &pdfName, int setNumber)
 
 ~GenericPDF ()=default
 
double pdf (PartonFlavor flavor, double x, double mu2) const
 Retrieves the collinear PDF value for a specific parton flavor. More...
 
void pdf (double x, double mu2, std::array< double, DEFAULT_TOTAL_PDFS > &output) const
 Evaluate the array of Collinear PDF values for {tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t}. More...
 
 GenericPDF (GenericPDF &&other) noexcept
 
GenericPDFoperator= (GenericPDF &&other) noexcept
 
 GenericPDF (const GenericPDF &other)
 
GenericPDFoperator= (const GenericPDF &other)
 
double tmd (PartonFlavor flavor, double x, double kt2, double mu2)
 Evaluates the TMD PDF value for a specific parton flavor. More...
 
void tmd (double x, double kt2, double mu2, std::array< double, DEFAULT_TOTAL_PDFS > &output)
 Evaluates the vector of TMD PDF values for {tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t}. More...
 
YamlStandardTMDInfo getStdPDFInfo () const
 Retrieves the standard PDF info. More...
 

Constructor & Destructor Documentation

◆ GenericPDF() [1/3]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::GenericPDF ( const std::string &  pdfName,
int  setNumber 
)
inline

◆ ~GenericPDF()

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::~GenericPDF ( )
default

◆ GenericPDF() [2/3]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::GenericPDF ( GenericPDF< Tag, Reader, Interpolator, Extrapolator > &&  other)
inlinenoexcept

◆ GenericPDF() [3/3]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::GenericPDF ( const GenericPDF< Tag, Reader, Interpolator, Extrapolator > &  other)
inline

Member Function Documentation

◆ getStdPDFInfo()

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
YamlStandardTMDInfo PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::getStdPDFInfo ( ) const
inline

Retrieves the standard PDF info.

Returns
YamlStandardTMDInfo The standard PDF info

◆ operator=() [1/2]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
GenericPDF& PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::operator= ( const GenericPDF< Tag, Reader, Interpolator, Extrapolator > &  other)
inline

◆ operator=() [2/2]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
GenericPDF& PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::operator= ( GenericPDF< Tag, Reader, Interpolator, Extrapolator > &&  other)
inlinenoexcept

◆ pdf() [1/2]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
void PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::pdf ( double  x,
double  mu2,
std::array< double, DEFAULT_TOTAL_PDFS > &  output 
) const
inline

Evaluate the array of Collinear PDF values for {tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t}.

Parameters
xBjorken x variable (momentum fraction)
mu2Factorization scale squared
outputThe array to store the Collinear PDF values for all flavors.
Returns
std::array<double, 13> The vector of Collinear PDF values
Exceptions
std::logic_errorIf called on a PDF type that doesn't support Collinear PDF

◆ pdf() [2/2]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
double PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::pdf ( PartonFlavor  flavor,
double  x,
double  mu2 
) const
inline

Retrieves the collinear PDF value for a specific parton flavor.

Parameters
flavorThe parton flavor
xBjorken x variable (momentum fraction)
mu2Factorization scale squared

◆ tmd() [1/2]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
void PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::tmd ( double  x,
double  kt2,
double  mu2,
std::array< double, DEFAULT_TOTAL_PDFS > &  output 
)
inline

Evaluates the vector of TMD PDF values for {tbar, bbar, cbar, sbar, ubar, dbar, g, d, u, s, c, b, t}.

Parameters
xBjorken x variable (momentum fraction)
kt2Transverse momentum squared
mu2Factorization scale squared
Returns
std::vector<double> The vector of TMD PDF values
Exceptions
std::logic_errorIf called on a PDF type that doesn't support TMD

◆ tmd() [2/2]

template<typename Tag , typename Reader = typename DefaultPDFImplementations<Tag>::Reader, typename Interpolator = typename DefaultPDFImplementations<Tag>::Interpolator, typename Extrapolator = typename DefaultPDFImplementations<Tag>::Extrapolator>
double PDFxTMD::GenericPDF< Tag, Reader, Interpolator, Extrapolator >::tmd ( PartonFlavor  flavor,
double  x,
double  kt2,
double  mu2 
)
inline

Evaluates the TMD PDF value for a specific parton flavor.

Parameters
flavorThe parton flavor
xBjorken x variable (momentum fraction)
kt2Transverse momentum squared

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