Back to home page

EIC code displayed by LXR

 
 

    


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

0001 //
0002 // APFEL++ 2017
0003 //
0004 // Author: Valerio Bertone: valerio.bertone@cern.ch
0005 //
0006 
0007 #pragma once
0008 
0009 #include "apfel/grid.h"
0010 #include "apfel/operator.h"
0011 #include "apfel/set.h"
0012 #include "apfel/dglapbuilder.h"
0013 #include "apfel/tabulateobject.h"
0014 #include "apfel/constants.h"
0015 
0016 namespace apfel
0017 {
0018   /**
0019    * @brief Structure that contains all precomputed quantities needed
0020    * to perform the TMD evolution, matching to the collinear PDFs, and
0021    * computation of cross sections, i.e. the perturbative coefficients
0022    * of matching functions, all anomalous dimensions, and hard
0023    * functions.
0024    */
0025   struct TmdObjects
0026   {
0027     double                                       Threshold;
0028     std::map<int, double>                        Beta;
0029     std::map<int, double>                        GammaFq;
0030     std::map<int, double>                        GammaFg;
0031     std::map<int, double>                        GammaK;
0032     std::map<int, std::vector<double>>           KCS;
0033     std::map<int, std::vector<Set<Operator>>>    MatchingFunctionsPDFs;
0034     std::map<int, std::vector<Set<Operator>>>    MatchingFunctionsFFs;
0035     std::map<std::string, std::map<int, double>> HardFactors;
0036   };
0037 
0038   /**
0039    * @name TMD object initializers
0040    * Collection of functions that initialise a TmdObjects structure
0041    * for the perturbartive evolution and the matching.
0042    */
0043   ///@{
0044   /**
0045    * @brief The InitializeTmdObjects function precomputes the
0046    * perturbative coefficients required for the evolution and matching
0047    * of TMD PDFs and FFs and store them into a 'TmdObjects' structure.
0048    * @param g: the x-space grid
0049    * @param Thresholds: the heavy quark thresholds
0050    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0051    * @return A map of TmdObject objects, one for each possible nf
0052    */
0053   std::map<int, TmdObjects> InitializeTmdObjects(Grid                const& g,
0054                                                  std::vector<double> const& Thresholds,
0055                                                  double              const& IntEps = 1e-5);
0056 
0057   /**
0058    * @brief The InitializeTmdObjectsDYResScheme function precomputes
0059    * the perturbative coefficients required for the evolution and
0060    * matching of TMD PDFs and FFs and store them into a 'TmdObjects'
0061    * structure. This function applies a resummation-scheme
0062    * transformation to produce the scheme often used in qT resummation
0063    * that has H = 1.
0064    * @param g: the x-space grid
0065    * @param Thresholds: the heavy quark thresholds
0066    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0067    * @return A map of TmdObject objects, one for each possible nf
0068    */
0069   std::map<int, TmdObjects> InitializeTmdObjectsDYResScheme(Grid                const& g,
0070                                                             std::vector<double> const& Thresholds,
0071                                                             double              const& IntEps = 1e-5);
0072 
0073   /**
0074    * @brief The InitializeTmdObjectsBM function precomputes the
0075    * perturbative coefficients required for the evolution and matching
0076    * of the (gluon) Boer-Mulders TMD PDF and store them into a
0077    * 'TmdObjects' structure. For now, quark and FF TMDs are not filled
0078    * in.
0079    * @param g: the x-space grid
0080    * @param Thresholds: the heavy quark thresholds
0081    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0082    * @return A map of TmdObject objects, one for each possible nf
0083    */
0084   std::map<int, TmdObjects> InitializeTmdObjectsBM(Grid                const& g,
0085                                                    std::vector<double> const& Thresholds,
0086                                                    double              const& IntEps = 1e-5);
0087 
0088   /**
0089    * @brief The InitializeTmdObjectsSivers function precomputes the
0090    * perturbative coefficients required for the evolution and matching
0091    * of the quark Sivers TMD PDF and store them into a 'TmdObjects'
0092    * structure. For now, gluon and FF TMDs (i.e. the Collins TMDs) are
0093    * not filled in. In addition, the matching is only present up to
0094    * one loop.
0095    * @param g: the x-space grid
0096    * @param Thresholds: the heavy quark thresholds
0097    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0098    * @return A map of TmdObject objects, one for each possible nf
0099    */
0100   std::map<int, TmdObjects> InitializeTmdObjectsSivers(Grid                const& g,
0101                                                        std::vector<double> const& Thresholds,
0102                                                        double              const& IntEps = 1e-5);
0103 
0104   /**
0105    * @brief The InitializeTmdObjects function precomputes the
0106    * perturbative coefficients required for the evolution and matching
0107    * of TMD g1 PDFs and FFs and store them into a 'TmdObjects'
0108    * structure. Matching functions are implemented only to one-loop
0109    * order.
0110    * @param g: the x-space grid
0111    * @param Thresholds: the heavy quark thresholds
0112    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0113    * @return A map of TmdObject objects, one for each possible nf
0114    */
0115   std::map<int, TmdObjects> InitializeTmdObjectsg1(Grid                const& g,
0116                                                    std::vector<double> const& Thresholds,
0117                                                    double              const& IntEps = 1e-5);
0118   ///@}
0119 
0120   /**
0121    * @name TMD builders
0122    * Collection of functions that build a TMD distributions (both PDFs
0123    * and FFs) as Set<Distribution>-valued functions. These functions
0124    * perform evolution and matching either separately or alltogether.
0125    * Also a function for the computation of the hard factors is
0126    * provided.
0127    */
0128   ///@{
0129   /**
0130    * @brief Function that returns the matched and evolved TMD PDFs in
0131    * b-space as functions of the final scale and rapidity.
0132    * @param TmdObj: the TMD objects
0133    * @param CollPDFs: the set of collinear PDFs to be matched
0134    * @param Alphas: the strong coupling function
0135    * @param PerturbativeOrder: the logarithmic perturbative order
0136    * @param Ci: the initial-scale variation factor (default: 1)
0137    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0138    * @return Set<Distribution>-valued function of the impact parameter
0139    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0140    * final rapidity scale &zeta; representing the evolved TMD PDFs
0141    */
0142   std::function<Set<Distribution>(double const&, double const&, double const&)> BuildTmdPDFs(std::map<int, TmdObjects>                       const& TmdObj,
0143                                                                                              std::function<Set<Distribution>(double const&)> const& CollPDFs,
0144                                                                                              std::function<double(double const&)>            const& Alphas,
0145                                                                                              int                                             const& PerturbativeOrder,
0146                                                                                              double                                          const& Ci = 1,
0147                                                                                              double                                          const& IntEps = 1e-7);
0148 
0149   /**
0150    * @brief Function that returns the matched and evolved TMD FFs in
0151    * b-space as functions of the final scale and rapidity.
0152    * @param TmdObj: the TMD objects
0153    * @param CollFFs: the set of collinear PDFs to be matched
0154    * @param Alphas: the strong coupling function
0155    * @param PerturbativeOrder: the logarithmic perturbative order
0156    * @param Ci: the initial-scale variation factor (default: 1)
0157    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0158    * @return Set<Distribution>-valued function of the impact parameter
0159    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0160    * final rapidity scale &zeta; representing the evolved TMD FFs
0161    */
0162   std::function<Set<Distribution>(double const&, double const&, double const&)> BuildTmdFFs(std::map<int, TmdObjects>                       const& TmdObj,
0163                                                                                             std::function<Set<Distribution>(double const&)> const& CollFFs,
0164                                                                                             std::function<double(double const&)>            const& Alphas,
0165                                                                                             int                                             const& PerturbativeOrder,
0166                                                                                             double                                          const& Ci = 1,
0167                                                                                             double                                          const& IntEps = 1e-7);
0168 
0169   /**
0170    * @brief Function that returns the matched TMD PDFs in b-space.
0171    * @param TmdObj: the TMD objects
0172    * @param CollPDFs: the set of collinear PDFs to be matched
0173    * @param Alphas: the strong coupling function
0174    * @param PerturbativeOrder: the logarithmic perturbative order
0175    * @param Ci: the initial-scale variation factor
0176    * @return Set<Distribution>-valued function of the impact parameter
0177    * b<SUB>T</SUB> representing the matched TMD PDFs
0178    */
0179   std::function<Set<Distribution>(double const&)> MatchTmdPDFs(std::map<int, TmdObjects>                       const& TmdObj,
0180                                                                std::function<Set<Distribution>(double const&)> const& CollPDFs,
0181                                                                std::function<double(double const&)>            const& Alphas,
0182                                                                int                                             const& PerturbativeOrder,
0183                                                                double                                          const& Ci = 1);
0184 
0185   /**
0186    * @brief Function that returns the matched TMD FFs in b-space.
0187    * @param TmdObj: the TMD objects
0188    * @param CollFFs: the set of collinear FFs to be matched
0189    * @param Alphas: the strong coupling function
0190    * @param PerturbativeOrder: the logarithmic perturbative order
0191    * @param Ci: the initial-scale variation factor
0192    * @return Set<Distribution>-valued function of the impact parameter
0193    * b<SUB>T</SUB> representing the matched TMD FFs
0194    */
0195   std::function<Set<Distribution>(double const&)> MatchTmdFFs(std::map<int, TmdObjects>                       const& TmdObj,
0196                                                               std::function<Set<Distribution>(double const&)> const& CollFFs,
0197                                                               std::function<double(double const&)>            const& Alphas,
0198                                                               int                                             const& PerturbativeOrder,
0199                                                               double                                          const& Ci = 1);
0200 
0201   /**
0202    * @brief Function that returns the mathing functions for the TMD PDFs.
0203    * @param TmdObj: the TMD objects
0204    * @param Alphas: the strong coupling function
0205    * @param PerturbativeOrder: the logarithmic perturbative order
0206    * @param Ci: the initial-scale variation factor
0207    * @return Set<Operator>-valued function of the scale mu
0208    * corresponding to the set of matching functions for PDFs in the
0209    * evolution basis.
0210    */
0211   std::function<Set<Operator>(double const&)> MatchingFunctionsPDFs(std::map<int, TmdObjects>            const& TmdObj,
0212                                                                     std::function<double(double const&)> const& Alphas,
0213                                                                     int                                  const& PerturbativeOrder,
0214                                                                     double                               const& Ci = 1);
0215 
0216   /**
0217    * @brief Function that returns the mathing functions for the TMD FFs.
0218    * @param TmdObj: the TMD objects
0219    * @param Alphas: the strong coupling function
0220    * @param PerturbativeOrder: the logarithmic perturbative order
0221    * @param Ci: the initial-scale variation factor
0222    * @return Set<Operator>-valued function of the scale mu
0223    * corresponding to the set of matching functions for FFs in the
0224    * evolution basis.
0225    */
0226   std::function<Set<Operator>(double const&)> MatchingFunctionsFFs(std::map<int, TmdObjects>            const& TmdObj,
0227                                                                    std::function<double(double const&)> const& Alphas,
0228                                                                    int                                  const& PerturbativeOrder,
0229                                                                    double                               const& Ci = 1);
0230 
0231   /**
0232    * @brief Function that returns the evolution factors for gluon and quarks.
0233    * @param TmdObj: the TMD objects
0234    * @param Alphas: the strong coupling function
0235    * @param PerturbativeOrder: the logarithmic perturbative order
0236    * @param Ci: the initial scale-variation factor (default: 1)
0237    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0238    * @return std::vector<double>-valued function of the impact parameter
0239    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0240    * final rapidity scale &zeta;. The 0-th component contains the
0241    * gluon evolution factor, the remaining 12, from 1 to 12, are all
0242    * equal and represent the quark evolution factors.
0243    */
0244   std::function<std::vector<double>(double const&, double const&, double const&)> EvolutionFactors(std::map<int, TmdObjects>            const& TmdObj,
0245                                                                                                    std::function<double(double const&)> const& Alphas,
0246                                                                                                    int                                  const& PerturbativeOrder,
0247                                                                                                    double                               const& Ci = 1,
0248                                                                                                    double                               const& IntEps = 1e-7);
0249 
0250   /**
0251    * @brief Function that returns the evolution factors for gluon and
0252    * quarks. As compared to "EvolutionFactors", this function isolates
0253    * the double logs into gammaK. This is reminiscent of the
0254    * qT-resummation typical way of computing the Sudakov form factor.
0255    * @param TmdObj: the TMD objects
0256    * @param Alphas: the strong coupling function
0257    * @param PerturbativeOrder: the logarithmic perturbative order
0258    * @param Ci: the initial scale-variation factor (default: 1)
0259    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0260    * @return std::vector<double>-valued function of the impact parameter
0261    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0262    * final rapidity scale &zeta;. The 0-th component contains the
0263    * gluon evolution factor, the remaining 12, from 1 to 12, are all
0264    * equal and represent the quark evolution factors.
0265    */
0266   std::function<std::vector<double>(double const&, double const&, double const&)> EvolutionFactorsK(std::map<int, TmdObjects>            const& TmdObj,
0267                                                                                                     std::function<double(double const&)> const& Alphas,
0268                                                                                                     int                                  const& PerturbativeOrder,
0269                                                                                                     double                               const& Ci = 1,
0270                                                                                                     double                               const& IntEps = 1e-7);
0271 
0272   /**
0273    * @brief Function that returns the evolution factor for quarks.
0274    * @param TmdObj: the TMD objects
0275    * @param Alphas: the strong coupling function
0276    * @param PerturbativeOrder: the logarithmic perturbative order
0277    * @param Ci: the initial scale-variation factor (default: 1)
0278    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0279    * @return double-valued function of the impact parameter
0280    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0281    * final rapidity scale &zeta;. It returns the quark evolution
0282    * factor.
0283    */
0284   std::function<double(double const&, double const&, double const&)> QuarkEvolutionFactor(std::map<int, TmdObjects>            const& TmdObj,
0285                                                                                           std::function<double(double const&)> const& Alphas,
0286                                                                                           int                                  const& PerturbativeOrder,
0287                                                                                           double                               const& Ci = 1,
0288                                                                                           double                               const& IntEps = 1e-7);
0289 
0290   /**
0291    * @brief Function that returns the evolution factor for quarks with
0292    * explicit dependence on the resummation-scale parameter.
0293    * @param TmdObj: the TMD objects
0294    * @param Alphas: the strong coupling function
0295    * @param PerturbativeOrder: the logarithmic perturbative order
0296    * @param xi: the resummation-scale parameter (default: 1)
0297    * @param Ci: the initial scale-variation factor (default: 1)
0298    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0299    * @return double-valued function of the impact parameter
0300    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0301    * final rapidity scale &zeta;. It returns the quark evolution
0302    * factor.
0303    */
0304   std::function<double(double const&, double const&, double const&)> QuarkEvolutionFactorxi(std::map<int, TmdObjects>            const& TmdObj,
0305                                                                                             std::function<double(double const&)> const& Alphas,
0306                                                                                             int                                  const& PerturbativeOrder,
0307                                                                                             double                               const& xi = 1,
0308                                                                                             double                               const& Ci = 1,
0309                                                                                             double                               const& IntEps = 1e-7);
0310 
0311   /**
0312    * @brief Function that returns the evolution factor for the gluon.
0313    * @param TmdObj: the TMD objects
0314    * @param Alphas: the strong coupling function
0315    * @param PerturbativeOrder: the logarithmic perturbative order
0316    * @param Ci: the initial scale-variation factor (default: 1)
0317    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0318    * @return double-valued function of the impact parameter
0319    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0320    * final rapidity scale &zeta;. It returns the gluon evolution
0321    * factor.
0322    */
0323   std::function<double(double const&, double const&, double const&)> GluonEvolutionFactor(std::map<int, TmdObjects>            const& TmdObj,
0324                                                                                           std::function<double(double const&)> const& Alphas,
0325                                                                                           int                                  const& PerturbativeOrder,
0326                                                                                           double                               const& Ci = 1,
0327                                                                                           double                               const& IntEps = 1e-7);
0328 
0329   /**
0330    * @brief Analytic evolution factor for the quark TMD
0331    * @param TmdObj: container of the anomalous dimensions with a fixed number of active flavours
0332    * @param mu0: the strong coupling reference scale
0333    * @param Alphas0: the value of the strong coupling and mu0
0334    * @param kappa: the resummation-scale parameter
0335    * @param kappa0: the ratio mu0 / M (i.e. the alphas reference scale over the hard scale)
0336    * @param PerturbativeOrder: the logarithmic perturbative accuracy
0337    * @return the quark Sudakov form factor as function of the impact parameter b
0338    */
0339   std::function<double(double const&)> QuarkAnalyticEvolutionFactor(TmdObjects const& TmdObj,
0340                                                                     double     const& mu0,
0341                                                                     double     const& Alphas0,
0342                                                                     double     const& kappa,
0343                                                                     double     const& kappa0,
0344                                                                     int        const& PerturbativeOrder);
0345 
0346   /**
0347    * @brief Analytic evolution factor for the gluon TMD
0348    * @param TmdObj: container of the anomalous dimensions with a fixed number of active flavours
0349    * @param mu0: the strong coupling reference scale
0350    * @param Alphas0: the value of the strong coupling and mu0
0351    * @param kappa: the resummation-scale parameter
0352    * @param kappa0: the ratio mu0 / M (i.e. the alphas reference scale over the hard scale)
0353    * @param PerturbativeOrder: the logarithmic perturbative accuracy
0354    * @return the gluon Sudakov form factor as function of the impact parameter b
0355    */
0356   std::function<double(double const&)> GluonAnalyticEvolutionFactor(TmdObjects const& TmdObj,
0357                                                                     double     const& mu0,
0358                                                                     double     const& Alphas0,
0359                                                                     double     const& kappa,
0360                                                                     double     const& kappa0,
0361                                                                     int        const& PerturbativeOrder);
0362 
0363   /**
0364    * @brief Function that returns the perturbative part of the
0365    * Collins-Soper kernel for quarks.
0366    * @param TmdObj: the TMD objects
0367    * @param Alphas: the strong coupling function
0368    * @param PerturbativeOrder: the logarithmic perturbative order
0369    * @param Ci: the initial scale-variation factor (default: 1)
0370    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0371    * @return double-valued function of the impact parameter
0372    * b<SUB>T</SUB> and of the the final renormalisation scale &mu;. It
0373    * returns perturbative part of the Collis-Soper kernel for quarks.
0374    */
0375   std::function<double(double const&, double const&)> CollinsSoperKernel(std::map<int, TmdObjects>            const& TmdObj,
0376                                                                          std::function<double(double const&)> const& Alphas,
0377                                                                          int                                  const& PerturbativeOrder,
0378                                                                          double                               const& Ci = 1,
0379                                                                          double                               const& IntEps = 1e-7);
0380 
0381   /**
0382    * @brief Function that returns the hard factor.
0383    * @param Process: the string corresponding to the process requested
0384    * @param TmdObj: the TMD objects
0385    * @param Alphas: the strong coupling function
0386    * @param PerturbativeOrder: the logarithmic perturbative order
0387    * @param Cf: the final scale-variation factor (default: 1)
0388    * @return double-valued function of the final renormalisation scale &mu;
0389    */
0390   std::function<double(double const&)> HardFactor(std::string                          const& Process,
0391                                                   std::map<int, TmdObjects>            const& TmdObj,
0392                                                   std::function<double(double const&)> const& Alphas,
0393                                                   int                                  const& PerturbativeOrder,
0394                                                   double                               const& Cf = 1);
0395   ///@}
0396 }