APFEL 4.8.0
A PDF evolution library in C++
Loading...
Searching...
No Matches
Classes | Namespaces | Enumerations
tmdbuilder.h File Reference
#include "apfel/grid.h"
#include "apfel/operator.h"
#include "apfel/set.h"
#include "apfel/dglapbuilder.h"
#include "apfel/tabulateobject.h"
#include "apfel/constants.h"

Go to the source code of this file.

Classes

struct  apfel::TmdObjects
 Structure that contains all precomputed quantities needed to perform the TMD evolution, matching to the collinear PDFs, and computation of cross sections, i.e. the perturbative coefficients of matching functions, all anomalous dimensions, and hard functions. More...
 

Namespaces

namespace  apfel
 Namespace for all APFEL++ functions and classes.
 

Enumerations

enum  apfel::JetAlgorithm : int { apfel::CONE = 0 , apfel::KT = 1 }
 Enumerator for the jet algoritms for the jet TMDs. More...
 

Functions

TMD object initializers

Collection of functions that initialise a TmdObjects structure for the perturbartive evolution and the matching.

std::map< int, TmdObjectsapfel::InitializeTmdObjects (Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
 The InitializeTmdObjects function precomputes the perturbative coefficients required for the evolution and matching of TMD PDFs and FFs and store them into a 'TmdObjects' structure.
 
std::map< int, TmdObjectsapfel::InitializeTmdObjectsDYResScheme (Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
 The InitializeTmdObjectsDYResScheme function precomputes the perturbative coefficients required for the evolution and matching of TMD PDFs and FFs and store them into a 'TmdObjects' structure. This function applies a resummation-scheme transformation to produce the scheme often used in qT resummation that has H = 1.
 
std::map< int, TmdObjectsapfel::InitializeTmdObjectsBM (Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
 The InitializeTmdObjectsBM function precomputes the perturbative coefficients required for the evolution and matching of the (gluon) Boer-Mulders TMD PDF and store them into a 'TmdObjects' structure. For now, quark and FF TMDs are not filled in.
 
std::map< int, TmdObjectsapfel::InitializeTmdObjectsSivers (Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
 The InitializeTmdObjectsSivers function precomputes the perturbative coefficients required for the evolution and matching of the quark Sivers TMD PDF and store them into a 'TmdObjects' structure. For now, gluon and FF TMDs (i.e. the Collins TMDs) are not filled in. In addition, the matching is only present up to one loop.
 
std::map< int, TmdObjectsapfel::InitializeTmdObjectsg1 (Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
 The InitializeTmdObjects function precomputes the perturbative coefficients required for the evolution and matching of TMD g1 PDFs and store them into a 'TmdObjects' structure.
 

Variables

TMD builders

Collection of functions that build a TMD distributions (both PDFs and FFs) as Set<Distribution>-valued functions. These functions perform evolution and matching either separately or alltogether. Also a function for the computation of the hard factors is provided.

std::function< Set< Distribution >(double const &, double const &, double const &) apfel::BuildTmdPDFs )(std::map< int, TmdObjects > const &TmdObj, std::function< Set< Distribution >(double const &)> const &CollPDFs, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the matched and evolved TMD PDFs in b-space as functions of the final scale and rapidity.
 
std::function< Set< Distribution >(double const &, double const &, double const &) apfel::BuildTmdFFs )(std::map< int, TmdObjects > const &TmdObj, std::function< Set< Distribution >(double const &)> const &CollFFs, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the matched and evolved TMD FFs in b-space as functions of the final scale and rapidity.
 
std::function< double(double const &, double const &, double const &) apfel::BuildTmdJet )(std::map< int, TmdObjects > const &TmdObj, JetAlgorithm const &JetAlgo, double const &JetR, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &CJ=1, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the TMD of a jet in b-space as functions of the final scale and rapidity.
 
std::function< Set< Distribution >(double const &) apfel::MatchTmdPDFs )(std::map< int, TmdObjects > const &TmdObj, std::function< Set< Distribution >(double const &)> const &CollPDFs, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1)
 Function that returns the matched TMD PDFs in b-space.
 
std::function< Set< Distribution >(double const &) apfel::MatchTmdFFs )(std::map< int, TmdObjects > const &TmdObj, std::function< Set< Distribution >(double const &)> const &CollFFs, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1)
 Function that returns the matched TMD FFs in b-space.
 
std::function< double(double const &, double const &) apfel::MatchTmdJet )(std::map< int, TmdObjects > const &TmdObj, JetAlgorithm const &JetAlgo, double const &tR, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &CJ=1, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the jet TMD in b-space at the initial scale.
 
std::function< Set< Operator >(double const &) apfel::MatchingFunctionsPDFs )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1)
 Function that returns the mathing functions for the TMD PDFs.
 
std::function< Set< Operator >(double const &) apfel::MatchingFunctionsFFs )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1)
 Function that returns the mathing functions for the TMD FFs.
 
std::function< std::vector< double >(double const &, double const &, double const &) apfel::EvolutionFactors )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the evolution factors for gluon and quarks.
 
std::function< std::vector< double >(double const &, double const &, double const &) apfel::EvolutionFactorsK )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the evolution factors for gluon and quarks. As compared to "EvolutionFactors", this function isolates the double logs into gammaK. This is reminiscent of the qT-resummation typical way of computing the Sudakov form factor.
 
std::function< double(double const &, double const &, double const &) apfel::QuarkEvolutionFactor )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the evolution factor for quarks.
 
std::function< double(double const &, double const &, double const &) apfel::QuarkEvolutionFactorxi )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &xi=1, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the evolution factor for quarks with explicit dependence on the resummation-scale parameter.
 
std::function< double(double const &, double const &, double const &) apfel::GluonEvolutionFactor )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the evolution factor for the gluon.
 
std::function< double(double const &) apfel::QuarkAnalyticEvolutionFactor )(TmdObjects const &TmdObj, double const &mu0, double const &Alphas0, double const &kappa, double const &kappa0, int const &PerturbativeOrder)
 Analytic evolution factor for the quark TMD.
 
std::function< double(double const &) apfel::GluonAnalyticEvolutionFactor )(TmdObjects const &TmdObj, double const &mu0, double const &Alphas0, double const &kappa, double const &kappa0, int const &PerturbativeOrder)
 Analytic evolution factor for the gluon TMD.
 
std::function< double(double const &, double const &) apfel::CollinsSoperKernel )(std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1, double const &IntEps=1e-7)
 Function that returns the perturbative part of the Collins-Soper kernel for quarks.
 
std::function< double(double const &) apfel::HardFactor )(std::string const &Process, std::map< int, TmdObjects > const &TmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Cf=1)
 Function that returns the hard factor.