Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-06-02 08:17:12

0001 //
0002 // APFEL++ 2017
0003 //
0004 // Author: Valerio Bertone: valerio.bertone@cern.ch
0005 //
0006 
0007 #pragma once
0008 
0009 #include "apfel/dglapbuilder.h"
0010 
0011 namespace apfel
0012 {
0013   /**
0014    * @name GPD object initializers
0015    * Collection of functions that initialise DglapObjects structure
0016    * for the different kinds of GPD evolution currently available.
0017    */
0018   ///@{
0019   /**
0020    * @brief The InitializeGPDObjects function precomputes the
0021    * perturbative coefficients of unpolarised GPD evolution kernels
0022    * and store them into a 'DglapObjects' structure. GPDs are assumed
0023    * to be continuous over heavy-quark thresholds.
0024    * @param g: the x-space grid
0025    * @param Thresholds: the heavy quark thresholds
0026    * @param xi: value of the skewness
0027    * @param OpEvol: the switch for the computation of the evolution operator (default: false)
0028    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0029    * @return A map of DglapObject objects, one for each possible nf
0030    */
0031   std::map<int, DglapObjects> InitializeGpdObjects(Grid                const& g,
0032                                                    std::vector<double> const& Thresholds,
0033                                                    double              const& xi,
0034                                                    bool                const& OpEvol = false,
0035                                                    double              const& IntEps = 1e-5);
0036 
0037   /**
0038    * @brief The InitializeGPDObjectsPol function precomputes the
0039    * perturbative coefficients of polarised GPD evolution kernels
0040    * and store them into a 'DglapObjects' structure. GPDs are assumed
0041    * to be continuous over heavy-quark thresholds.
0042    * @param g: the x-space grid
0043    * @param Thresholds: the heavy quark thresholds
0044    * @param xi: value of the skewness
0045    * @param OpEvol: the switch for the computation of the evolution operator (default: false)
0046    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0047    * @return A map of DglapObject objects, one for each possible nf
0048    */
0049   std::map<int, DglapObjects> InitializeGpdObjectsPol(Grid                const& g,
0050                                                       std::vector<double> const& Thresholds,
0051                                                       double              const& xi,
0052                                                       bool                const& OpEvol = false,
0053                                                       double              const& IntEps = 1e-5);
0054 
0055   /**
0056    * @brief The InitializeGPDObjectsTrans function precomputes the
0057    * perturbative coefficients of transversely polarised GPD evolution
0058    * kernels and store them into a 'DglapObjects' structure. GPDs are
0059    * assumed to be continuous over heavy-quark thresholds.
0060    * @param g: the x-space grid
0061    * @param Thresholds: the heavy quark thresholds
0062    * @param xi: value of the skewness
0063    * @param OpEvol: the switch for the computation of the evolution operator (default: false)
0064    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0065    * @return A map of DglapObject objects, one for each possible nf
0066    */
0067   std::map<int, DglapObjects> InitializeGpdObjectsTrans(Grid                const& g,
0068                                                         std::vector<double> const& Thresholds,
0069                                                         double              const& xi,
0070                                                         bool                const& OpEvol = false,
0071                                                         double              const& IntEps = 1e-5);
0072   ///@}
0073 }