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 //         Simone Rodini: rodini.simone.luigi@gmail.com
0006 //
0007 
0008 #pragma once
0009 
0010 #include "apfel/grid.h"
0011 #include "apfel/operator.h"
0012 #include "apfel/set.h"
0013 #include "apfel/dglapbuilder.h"
0014 #include "apfel/tabulateobject.h"
0015 #include "apfel/constants.h"
0016 
0017 namespace apfel
0018 {
0019   struct GtmdObjects
0020   {
0021     double                                    Threshold;
0022     double                                    xi;
0023     std::map<int, double>                     Beta;
0024     std::map<int, double>                     GammaFq;
0025     std::map<int, double>                     GammaFg;
0026     std::map<int, double>                     GammaK;
0027     std::map<int, std::vector<double>>        KCS;
0028     std::map<int, std::vector<Set<Operator>>> MatchingFunctions;
0029   };
0030 
0031   /**
0032    * @name GTMD object initializers
0033    * Collection of functions that initialise GtmdObjects structure for
0034    * the perturbartive evolution and matching currently available.
0035    */
0036   ///@{
0037   /**
0038    * @brief The InitializeGtmdObjectsEvenUU function precomputes the
0039    * perturbative coefficients required for the evolution and matching
0040    * of the parity even part of twist-2 GTMDs for the UU polarisation channel
0041    * and store them into a 'GtmdObjects' structure.
0042    * @param g: the x-space grid
0043    * @param Thresholds: the heavy quark thresholds
0044    * @param xi: value of the skewness
0045    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0046    * @return A map of GtmdObjects, one for each possible nf
0047    */
0048   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenUU(Grid                const& g,
0049                                                          std::vector<double> const& Thresholds,
0050                                                          double              const& xi,
0051                                                          double              const& IntEps = 1e-5);
0052 
0053   /**
0054    * @brief The InitializeGtmdObjectsOddUU function precomputes the
0055    * perturbative coefficients required for the evolution and matching
0056    * of the parity odd part of twist-2 GTMDs for the UU polarisation channel
0057    * and store them into a 'GtmdObjects' structure.
0058    * @param g: the x-space grid
0059    * @param Thresholds: the heavy quark thresholds
0060    * @param xi: value of the skewness
0061    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0062    * @return A map of GtmdObjects, one for each possible nf
0063    */
0064   std::map<int, GtmdObjects> InitializeGtmdObjectsOddUU(Grid                const& g,
0065                                                         std::vector<double> const& Thresholds,
0066                                                         double              const& xi,
0067                                                         double              const& IntEps = 1e-5);
0068 
0069   /**
0070    * @brief The InitializeGtmdObjectsEvenLL function precomputes the
0071    * perturbative coefficients required for the evolution and matching
0072    * of the parity even part of twist-2 GTMDs for the LL polarisation channel
0073    * and store them into a 'GtmdObjects' structure.
0074    * @param g: the x-space grid
0075    * @param Thresholds: the heavy quark thresholds
0076    * @param xi: value of the skewness
0077    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0078    * @return A map of GtmdObjects, one for each possible nf
0079    */
0080   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenLL(Grid                const& g,
0081                                                          std::vector<double> const& Thresholds,
0082                                                          double              const& xi,
0083                                                          double              const& IntEps = 1e-5);
0084 
0085   /**
0086    * @brief The InitializeGtmdObjectsOddLL function precomputes the
0087    * perturbative coefficients required for the evolution and matching
0088    * of the parity odd part of twist-2 GTMDs for the LL polarisation channel
0089    * and store them into a 'GtmdObjects' structure.
0090    * @param g: the x-space grid
0091    * @param Thresholds: the heavy quark thresholds
0092    * @param xi: value of the skewness
0093    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0094    * @return A map of GtmdObjects, one for each possible nf
0095    */
0096   std::map<int, GtmdObjects> InitializeGtmdObjectsOddLL(Grid                const& g,
0097                                                         std::vector<double> const& Thresholds,
0098                                                         double              const& xi,
0099                                                         double              const& IntEps = 1e-5);
0100 
0101   /**
0102    * @brief The InitializeGtmdObjectsEvenTT function precomputes the
0103    * perturbative coefficients required for the evolution and matching
0104    * of the parity even part of twist-2 GTMDs for the TT polarisation channel
0105    * and store them into a 'GtmdObjects' structure.
0106    * @param g: the x-space grid
0107    * @param Thresholds: the heavy quark thresholds
0108    * @param xi: value of the skewness
0109    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0110    * @return A map of GtmdObjects, one for each possible nf
0111    */
0112   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenTT(Grid                const& g,
0113                                                          std::vector<double> const& Thresholds,
0114                                                          double              const& xi,
0115                                                          double              const& IntEps = 1e-5);
0116 
0117   /**
0118    * @brief The InitializeGtmdObjectsOddTT function precomputes the
0119    * perturbative coefficients required for the evolution and matching
0120    * of the parity odd part of twist-2 GTMDs for the TT polarisation channel
0121    * and store them into a 'GtmdObjects' structure.
0122    * @param g: the x-space grid
0123    * @param Thresholds: the heavy quark thresholds
0124    * @param xi: value of the skewness
0125    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0126    * @return A map of GtmdObjects, one for each possible nf
0127    */
0128   std::map<int, GtmdObjects> InitializeGtmdObjectsOddTT(Grid                const& g,
0129                                                         std::vector<double> const& Thresholds,
0130                                                         double              const& xi,
0131                                                         double              const& IntEps = 1e-5);
0132 
0133   /**
0134    * @brief The InitializeGtmdObjectsEvenUT function precomputes the
0135    * perturbative coefficients required for the evolution and matching
0136    * of the parity even part of twist-2 GTMDs for the UT polarisation channel
0137    * and store them into a 'GtmdObjects' structure.
0138    * @param g: the x-space grid
0139    * @param Thresholds: the heavy quark thresholds
0140    * @param xi: value of the skewness
0141    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0142    * @return A map of GtmdObjects, one for each possible nf
0143    */
0144   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenUT(Grid                const& g,
0145                                                          std::vector<double> const& Thresholds,
0146                                                          double              const& xi,
0147                                                          double              const& IntEps = 1e-5);
0148 
0149   /**
0150    * @brief The InitializeGtmdObjectsOddUT function precomputes the
0151    * perturbative coefficients required for the evolution and matching
0152    * of the parity odd part of twist-2 GTMDs for the UT polarisation channel
0153    * and store them into a 'GtmdObjects' structure.
0154    * @param g: the x-space grid
0155    * @param Thresholds: the heavy quark thresholds
0156    * @param xi: value of the skewness
0157    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0158    * @return A map of GtmdObjects, one for each possible nf
0159    */
0160   std::map<int, GtmdObjects> InitializeGtmdObjectsOddUT(Grid                const& g,
0161                                                         std::vector<double> const& Thresholds,
0162                                                         double              const& xi,
0163                                                         double              const& IntEps = 1e-5);
0164 
0165   /**
0166    * @brief The InitializeGtmdObjectsEvenLT function precomputes the
0167    * perturbative coefficients required for the evolution and matching
0168    * of the parity even part of twist-2 GTMDs for the LT polarisation channel
0169    * and store them into a 'GtmdObjects' structure.
0170    * @param g: the x-space grid
0171    * @param Thresholds: the heavy quark thresholds
0172    * @param xi: value of the skewness
0173    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0174    * @return A map of GtmdObjects, one for each possible nf
0175    */
0176   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenLT(Grid                const& g,
0177                                                          std::vector<double> const& Thresholds,
0178                                                          double              const& xi,
0179                                                          double              const& IntEps = 1e-5);
0180 
0181   /**
0182    * @brief The InitializeGtmdObjectsOddLT function precomputes the
0183    * perturbative coefficients required for the evolution and matching
0184    * of the parity odd part of twist-2 GTMDs for the LT polarisation channel
0185    * and store them into a 'GtmdObjects' structure.
0186    * @param g: the x-space grid
0187    * @param Thresholds: the heavy quark thresholds
0188    * @param xi: value of the skewness
0189    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0190    * @return A map of GtmdObjects, one for each possible nf
0191    */
0192   std::map<int, GtmdObjects> InitializeGtmdObjectsOddLT(Grid                const& g,
0193                                                         std::vector<double> const& Thresholds,
0194                                                         double              const& xi,
0195                                                         double              const& IntEps = 1e-5);
0196 
0197   /**
0198    * @brief The InitializeGtmdObjectsEvenTU function precomputes the
0199    * perturbative coefficients required for the evolution and matching
0200    * of the parity even part of twist-2 GTMDs for the TU polarisation channel
0201    * and store them into a 'GtmdObjects' structure.
0202    * @param g: the x-space grid
0203    * @param Thresholds: the heavy quark thresholds
0204    * @param xi: value of the skewness
0205    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0206    * @return A map of GtmdObjects, one for each possible nf
0207    */
0208   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenTU(Grid                const& g,
0209                                                          std::vector<double> const& Thresholds,
0210                                                          double              const& xi,
0211                                                          double              const& IntEps = 1e-5);
0212 
0213   /**
0214    * @brief The InitializeGtmdObjectsOddTU function precomputes the
0215    * perturbative coefficients required for the evolution and matching
0216    * of the parity odd part of twist-2 GTMDs for the TU polarisation channel
0217    * and store them into a 'GtmdObjects' structure.
0218    * @param g: the x-space grid
0219    * @param Thresholds: the heavy quark thresholds
0220    * @param xi: value of the skewness
0221    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0222    * @return A map of GtmdObjects, one for each possible nf
0223    */
0224   std::map<int, GtmdObjects> InitializeGtmdObjectsOddTU(Grid                const& g,
0225                                                         std::vector<double> const& Thresholds,
0226                                                         double              const& xi,
0227                                                         double              const& IntEps = 1e-5);
0228 
0229   /**
0230    * @brief The InitializeGtmdObjectsEvenTL function precomputes the
0231    * perturbative coefficients required for the evolution and matching
0232    * of the parity even part of twist-2 GTMDs for the TL polarisation channel
0233    * and store them into a 'GtmdObjects' structure.
0234    * @param g: the x-space grid
0235    * @param Thresholds: the heavy quark thresholds
0236    * @param xi: value of the skewness
0237    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0238    * @return A map of GtmdObjects, one for each possible nf
0239    */
0240   std::map<int, GtmdObjects> InitializeGtmdObjectsEvenTL(Grid                const& g,
0241                                                          std::vector<double> const& Thresholds,
0242                                                          double              const& xi,
0243                                                          double              const& IntEps = 1e-5);
0244 
0245   /**
0246    * @brief The InitializeGtmdObjectsOddTL function precomputes the
0247    * perturbative coefficients required for the evolution and matching
0248    * of the parity odd part of twist-2 GTMDs for the TL polarisation channel
0249    * and store them into a 'GtmdObjects' structure.
0250    * @param g: the x-space grid
0251    * @param Thresholds: the heavy quark thresholds
0252    * @param xi: value of the skewness
0253    * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>)
0254    * @return A map of GtmdObjects, one for each possible nf
0255    */
0256   std::map<int, GtmdObjects> InitializeGtmdObjectsOddTL(Grid                const& g,
0257                                                         std::vector<double> const& Thresholds,
0258                                                         double              const& xi,
0259                                                         double              const& IntEps = 1e-5);
0260   ///@}
0261 
0262   /**
0263    * @name GTMD builders
0264    * Collection of functions that build a GTMD distributions as
0265    * Set<Distribution>-valued functions. These functions perform
0266    * evolution and matching either separately or alltogether.
0267    */
0268   ///@{
0269   /**
0270    * @brief Function that returns the matched and evolved GTMDs in
0271    * b-space as functions of the final scale and rapidity.
0272    * @param GtmdObj: the GTMD objects
0273    * @param CollGPDs: the set of collinear GPDs to be matched
0274    * @param Alphas: the strong coupling function
0275    * @param PerturbativeOrder: the perturbative order
0276    * @param Ci: the initial-scale variation factor (default: 1)
0277    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0278    * @return Set<Distribution>-valued function of the impact parameter
0279    * b<SUB>T</SUB>, the final renormalisation scale &mu;, and the
0280    * final rapidity scale &zeta; representing the evolved GTMDs.
0281    */
0282   std::function<Set<Distribution>(double const&, double const&, double const&)> BuildGtmds(std::map<int, GtmdObjects>                      const& GtmdObj,
0283                                                                                            std::function<Set<Distribution>(double const&)> const& CollGPDs,
0284                                                                                            std::function<double(double const&)>            const& Alphas,
0285                                                                                            int                                             const& PerturbativeOrder,
0286                                                                                            double                                          const& Ci = 1,
0287                                                                                            double                                          const& IntEps = 1e-7);
0288 
0289   /**
0290    * @brief Function that returns the matched TMD GPDs in b-space.
0291    * @param GtmdObj: the GTMD objects
0292    * @param CollGPDs: the set of collinear GPDs to be matched
0293    * @param Alphas: the strong coupling function
0294    * @param PerturbativeOrder: the perturbative order
0295    * @param Ci: the initial-scale variation factor
0296    * @return Set<Distribution>-valued function of the impact parameter
0297    * b<SUB>T</SUB> representing the matched GTMDs.
0298    */
0299   std::function<Set<Distribution>(double const&)> MatchGtmds(std::map<int, GtmdObjects>                      const& GtmdObj,
0300                                                              std::function<Set<Distribution>(double const&)> const& CollGPDs,
0301                                                              std::function<double(double const&)>            const& Alphas,
0302                                                              int                                             const& PerturbativeOrder,
0303                                                              double                                          const& Ci = 1);
0304 
0305   /**
0306    * @brief Function that returns the mathing functions for the GTMDs.
0307    * @param GtmdObj: the GTMD objects
0308    * @param Alphas: the strong coupling function
0309    * @param PerturbativeOrder: the perturbative order
0310    * @param Ci: the initial-scale variation factor
0311    * @return Set<Operator>-valued function of the scale mu
0312    * corresponding to the set of matching functions for GPDs in the
0313    * evolution basis.
0314    */
0315   std::function<Set<Operator>(double const&)> MatchingFunctions(std::map<int, GtmdObjects>           const& GtmdObj,
0316                                                                 std::function<double(double const&)> const& Alphas,
0317                                                                 int                                  const& PerturbativeOrder,
0318                                                                 double                               const& Ci = 1);
0319   /// @cond UNNECESSARY
0320   /**
0321    * @brief Function that returns the evolution factors for gluon and quarks.
0322    * @param GtmdObj: the GTMD objects
0323    * @param Alphas: the strong coupling function
0324    * @param PerturbativeOrder: the perturbative order
0325    * @param Ci: the initial scale-variation factor (default: 1)
0326    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0327    * @return std::vector<double>-valued function of the longitudinal
0328    * momentum fraction x, the impact parameter b<SUB>T</SUB>, the
0329    * final renormalisation scale &mu;, and the final rapidity scale
0330    * &zeta;. The 0-th component contains the gluon evolution factor,
0331    * the remaining 12, from 1 to 12, are all equal and represent the
0332    * quark evolution factors.
0333    */
0334   std::function<std::vector<double>(double const&, double const&, double const&, double const&)> EvolutionFactors(std::map<int, GtmdObjects>           const& GtmdObj,
0335                                                                                                                   std::function<double(double const&)> const& Alphas,
0336                                                                                                                   int                                  const& PerturbativeOrder,
0337                                                                                                                   double                               const& Ci = 1,
0338                                                                                                                   double                               const& IntEps = 1e-7);
0339 
0340   /**
0341    * @brief Function that returns the evolution factor for quarks.
0342    * @param GtmdObj: the GTMD objects
0343    * @param Alphas: the strong coupling function
0344    * @param PerturbativeOrder: the perturbative order
0345    * @param Ci: the initial scale-variation factor (default: 1)
0346    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0347    * @return double-valued function of the longitudinal momentum
0348    * fraction x, the impact parameter b<SUB>T</SUB>, the final
0349    * renormalisation scale &mu;, and the final rapidity scale
0350    * &zeta;. It returns the quark evolution factor.
0351    */
0352   std::function<double(double const&, double const&, double const&, double const&)> QuarkEvolutionFactor(std::map<int, GtmdObjects>           const& GtmdObj,
0353                                                                                                          std::function<double(double const&)> const& Alphas,
0354                                                                                                          int                                  const& PerturbativeOrder,
0355                                                                                                          double                               const& Ci = 1,
0356                                                                                                          double                               const& IntEps = 1e-7);
0357 
0358   /**
0359    * @brief Function that returns the evolution factor for the gluon.
0360    * @param GtmdObj: the GTMD objects
0361    * @param Alphas: the strong coupling function
0362    * @param PerturbativeOrder: the perturbative order
0363    * @param Ci: the initial scale-variation factor (default: 1)
0364    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0365    * @return double-valued function of the longitudinal momentum
0366    * fraction x, the impact parameter b<SUB>T</SUB>, the final
0367    * renormalisation scale &mu;, and the final rapidity scale
0368    * &zeta;. It returns the gluon evolution factor.
0369    */
0370   std::function<double(double const&, double const&, double const&, double const&)> GluonEvolutionFactor(std::map<int, GtmdObjects>           const& GtmdObj,
0371                                                                                                          std::function<double(double const&)> const& Alphas,
0372                                                                                                          int                                  const& PerturbativeOrder,
0373                                                                                                          double                               const& Ci = 1,
0374                                                                                                          double                               const& IntEps = 1e-7);
0375 
0376   /**
0377    * @brief Function that returns the perturbative part of the
0378    * Collins-Soper kernel for quarks.
0379    * @param GtmdObj: the GTMD objects
0380    * @param Alphas: the strong coupling function
0381    * @param PerturbativeOrder: the logarithmic perturbative order
0382    * @param Ci: the initial scale-variation factor (default: 1)
0383    * @param IntEps: the integration accuracy (default: 10<SUP>-7</SUP>)
0384    * @return double-valued function of the impact parameter
0385    * b<SUB>T</SUB> and of the the final renormalisation scale &mu;. It
0386    * returns perturbative part of the Collis-Soper kernel for quarks.
0387    */
0388   std::function<double(double const&, double const&)> CollinsSoperKernel(std::map<int, GtmdObjects>           const& GtmdObj,
0389                                                                          std::function<double(double const&)> const& Alphas,
0390                                                                          int                                  const& PerturbativeOrder,
0391                                                                          double                               const& Ci = 1,
0392                                                                          double                               const& IntEps = 1e-7);
0393   /// @endcond
0394   ///@}
0395 }