APFEL 4.8.0
A PDF evolution library in C++
Loading...
Searching...
No Matches
gtmdbuilder.h
Go to the documentation of this file.
1//
2// APFEL++ 2017
3//
4// Author: Valerio Bertone: valerio.bertone@cern.ch
5//
6
7#pragma once
8
9#include "apfel/grid.h"
10#include "apfel/operator.h"
11#include "apfel/set.h"
12#include "apfel/dglapbuilder.h"
14#include "apfel/constants.h"
15
16namespace apfel
17{
19 {
20 double Threshold;
21 double xi;
22 std::map<int, double> Beta;
23 std::map<int, double> GammaFq;
24 std::map<int, double> GammaFg;
25 std::map<int, double> GammaK;
26 std::map<int, std::vector<double>> KCS;
27 std::map<int, std::vector<Set<Operator>>> MatchingFunctions;
28 };
29
36
46 std::map<int, GtmdObjects> InitializeGtmdObjects(Grid const& g,
47 std::vector<double> const& Thresholds,
48 double const& xi,
49 double const& IntEps = 1e-5);
51
59
72 std::function<Set<Distribution>(double const&, double const&, double const&)> BuildGtmds(std::map<int, GtmdObjects> const& GtmdObj,
73 std::function<Set<Distribution>(double const&)> const& CollGPDs,
74 std::function<double(double const&)> const& Alphas,
75 int const& PerturbativeOrder,
76 double const& Ci = 1,
77 double const& IntEps = 1e-7);
78
89 std::function<Set<Distribution>(double const&)> MatchGtmds(std::map<int, GtmdObjects> const& GtmdObj,
90 std::function<Set<Distribution>(double const&)> const& CollGPDs,
91 std::function<double(double const&)> const& Alphas,
92 int const& PerturbativeOrder,
93 double const& Ci = 1);
94
105 std::function<Set<Operator>(double const&)> MatchingFunctions(std::map<int, GtmdObjects> const& GtmdObj,
106 std::function<double(double const&)> const& Alphas,
107 int const& PerturbativeOrder,
108 double const& Ci = 1);
110
123 std::function<std::vector<double>(double const&, double const&, double const&)> EvolutionFactors(std::map<int, GtmdObjects> const& GtmdObj,
124 std::function<double(double const&)> const& Alphas,
125 int const& PerturbativeOrder,
126 double const& Ci = 1,
127 double const& IntEps = 1e-7);
128
141 std::function<double(double const&, double const&, double const&)> QuarkEvolutionFactor(std::map<int, GtmdObjects> const& GtmdObj,
142 std::function<double(double const&)> const& Alphas,
143 int const& PerturbativeOrder,
144 double const& Ci = 1,
145 double const& IntEps = 1e-7);
146
159 std::function<double(double const&, double const&, double const&)> GluonEvolutionFactor(std::map<int, GtmdObjects> const& GtmdObj,
160 std::function<double(double const&)> const& Alphas,
161 int const& PerturbativeOrder,
162 double const& Ci = 1,
163 double const& IntEps = 1e-7);
166}
The Grid class defines ab object that is essentially a collection of "SubGrid" objects plus other glo...
Definition grid.h:22
The Set template class allocates a collection of objects of type T along the ConvolutionMap and provi...
Definition set.h:22
Namespace for all APFEL++ functions and classes.
Definition alphaqcd.h:14
std::function< double(double const &, double const &, double const &) 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::map< int, GtmdObjects > InitializeGtmdObjects(Grid const &g, std::vector< double > const &Thresholds, double const &xi, double const &IntEps=1e-5)
The InitializeGtmdObjects function precomputes the perturbative coefficients required for the evoluti...
std::function< Set< Distribution >(double const &) MatchGtmds)(std::map< int, GtmdObjects > const &GtmdObj, std::function< Set< Distribution >(double const &)> const &CollGPDs, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1)
Function that returns the matched TMD GPDs in b-space.
std::function< double(double const &, double const &, double const &) 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< Set< Operator >(double const &) MatchingFunctions)(std::map< int, GtmdObjects > const &GtmdObj, std::function< double(double const &)> const &Alphas, int const &PerturbativeOrder, double const &Ci=1)
Function that returns the mathing functions for the GTMDs.
std::function< std::vector< double >(double const &, double const &, double const &) 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< Set< Distribution >(double const &, double const &, double const &) BuildGtmds)(std::map< int, GtmdObjects > const &GtmdObj, std::function< Set< Distribution >(double const &)> const &CollGPDs, 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 GTMDs in b-space as functions of the final scale and ra...
Definition gtmdbuilder.h:19
std::map< int, double > GammaK
Definition gtmdbuilder.h:25
std::map< int, double > Beta
Definition gtmdbuilder.h:22
std::map< int, double > GammaFq
Definition gtmdbuilder.h:23
std::map< int, std::vector< double > > KCS
Definition gtmdbuilder.h:26
double Threshold
Definition gtmdbuilder.h:20
std::map< int, std::vector< Set< Operator > > > MatchingFunctions
Definition gtmdbuilder.h:27
std::map< int, double > GammaFg
Definition gtmdbuilder.h:24
double xi
Definition gtmdbuilder.h:21