33 std::map<int, double>
Beta;
37 std::map<int, std::vector<double>>
KCS;
59 std::vector<double>
const& Thresholds,
60 double const& IntEps = 1e-5);
75 std::vector<double>
const& Thresholds,
76 double const& IntEps = 1e-5);
90 std::vector<double>
const& Thresholds,
91 double const& IntEps = 1e-5);
106 std::vector<double>
const& Thresholds,
107 double const& IntEps = 1e-5);
119 std::vector<double>
const& Thresholds,
120 double const& IntEps = 1e-5);
145 std::function<Set<Distribution>(
double const&,
double const&,
double const&)>
BuildTmdPDFs(std::map<int, TmdObjects>
const& TmdObj,
147 std::function<
double(
double const&)>
const& Alphas,
148 int const& PerturbativeOrder,
149 double const& Ci = 1,
150 double const& IntEps = 1e-7);
165 std::function<Set<Distribution>(
double const&,
double const&,
double const&)>
BuildTmdFFs(std::map<int, TmdObjects>
const& TmdObj,
167 std::function<
double(
double const&)>
const& Alphas,
168 int const& PerturbativeOrder,
169 double const& Ci = 1,
170 double const& IntEps = 1e-7);
187 std::function<double(
double const&,
double const&,
double const&)>
BuildTmdJet(std::map<int, TmdObjects>
const& TmdObj,
190 std::function<
double(
double const&)>
const& Alphas,
191 int const& PerturbativeOrder,
192 double const& CJ = 1,
193 double const& Ci = 1,
194 double const& IntEps = 1e-7);
206 std::function<Set<Distribution>(
double const&)>
MatchTmdPDFs(std::map<int, TmdObjects>
const& TmdObj,
208 std::function<
double(
double const&)>
const& Alphas,
209 int const& PerturbativeOrder,
210 double const& Ci = 1);
222 std::function<Set<Distribution>(
double const&)>
MatchTmdFFs(std::map<int, TmdObjects>
const& TmdObj,
224 std::function<
double(
double const&)>
const& Alphas,
225 int const& PerturbativeOrder,
226 double const& Ci = 1);
242 std::function<double(
double const&,
double const&)>
MatchTmdJet(std::map<int, TmdObjects>
const& TmdObj,
245 std::function<
double(
double const&)>
const& Alphas,
246 int const& PerturbativeOrder,
247 double const& CJ = 1,
248 double const& Ci = 1,
249 double const& IntEps = 1e-7);
262 std::function<
double(
double const&)>
const& Alphas,
263 int const& PerturbativeOrder,
264 double const& Ci = 1);
277 std::function<
double(
double const&)>
const& Alphas,
278 int const& PerturbativeOrder,
279 double const& Ci = 1);
294 std::function<std::vector<double>(
double const&,
double const&,
double const&)>
EvolutionFactors(std::map<int, TmdObjects>
const& TmdObj,
295 std::function<
double(
double const&)>
const& Alphas,
296 int const& PerturbativeOrder,
297 double const& Ci = 1,
298 double const& IntEps = 1e-7);
316 std::function<std::vector<double>(
double const&,
double const&,
double const&)>
EvolutionFactorsK(std::map<int, TmdObjects>
const& TmdObj,
317 std::function<
double(
double const&)>
const& Alphas,
318 int const& PerturbativeOrder,
319 double const& Ci = 1,
320 double const& IntEps = 1e-7);
334 std::function<double(
double const&,
double const&,
double const&)>
QuarkEvolutionFactor(std::map<int, TmdObjects>
const& TmdObj,
335 std::function<
double(
double const&)>
const& Alphas,
336 int const& PerturbativeOrder,
337 double const& Ci = 1,
338 double const& IntEps = 1e-7);
354 std::function<double(
double const&,
double const&,
double const&)>
QuarkEvolutionFactorxi(std::map<int, TmdObjects>
const& TmdObj,
355 std::function<
double(
double const&)>
const& Alphas,
356 int const& PerturbativeOrder,
357 double const& xi = 1,
358 double const& Ci = 1,
359 double const& IntEps = 1e-7);
373 std::function<double(
double const&,
double const&,
double const&)>
GluonEvolutionFactor(std::map<int, TmdObjects>
const& TmdObj,
374 std::function<
double(
double const&)>
const& Alphas,
375 int const& PerturbativeOrder,
376 double const& Ci = 1,
377 double const& IntEps = 1e-7);
391 double const& Alphas0,
393 double const& kappa0,
394 int const& PerturbativeOrder);
408 double const& Alphas0,
410 double const& kappa0,
411 int const& PerturbativeOrder);
425 std::function<double(
double const&,
double const&)>
CollinsSoperKernel(std::map<int, TmdObjects>
const& TmdObj,
426 std::function<
double(
double const&)>
const& Alphas,
427 int const& PerturbativeOrder,
428 double const& Ci = 1,
429 double const& IntEps = 1e-7);
440 std::function<double(
double const&)>
HardFactor(std::string
const& Process,
441 std::map<int, TmdObjects>
const& TmdObj,
442 std::function<
double(
double const&)>
const& Alphas,
443 int const& PerturbativeOrder,
444 double const& Cf = 1);
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 &)> 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< Set< Operator >(double const &)> 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::map< int, TmdObjects > InitializeTmdObjectsDYResScheme(Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
The InitializeTmdObjectsDYResScheme function precomputes the perturbative coefficients required for t...
std::map< int, TmdObjects > InitializeTmdObjectsSivers(Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
The InitializeTmdObjectsSivers function precomputes the perturbative coefficients required for the ev...
std::function< Set< Distribution >(double const &, double const &, double const &)> 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...
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< 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::function< Set< Distribution >(double const &)> 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< Set< Distribution >(double const &)> 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< std::vector< double >(double const &, double const &, double const &)> 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",...
std::function< Set< Distribution >(double const &, double const &, double const &)> 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 ...
std::function< double(double const &, double const &, double const &)> 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< double(double const &)> 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.
std::function< double(double const &, double const &)> 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< double(double const &, double const &)> 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 &, double const &, double const &)> 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-sca...
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::map< int, TmdObjects > InitializeTmdObjects(Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
The InitializeTmdObjects function precomputes the perturbative coefficients required for the evolutio...
std::map< int, TmdObjects > InitializeTmdObjectsg1(Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
The InitializeTmdObjects function precomputes the perturbative coefficients required for the evolutio...
std::function< Set< Operator >(double const &)> 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.
JetAlgorithm
Enumerator for the jet algoritms for the jet TMDs.
Definition tmdbuilder.h:21
@ KT
Definition tmdbuilder.h:21
@ CONE
Definition tmdbuilder.h:21
std::map< int, TmdObjects > InitializeTmdObjectsBM(Grid const &g, std::vector< double > const &Thresholds, double const &IntEps=1e-5)
The InitializeTmdObjectsBM function precomputes the perturbative coefficients required for the evolut...
std::function< double(double const &)> 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.
Structure that contains all precomputed quantities needed to perform the TMD evolution,...
Definition tmdbuilder.h:31
std::map< int, double > Beta
Definition tmdbuilder.h:33
std::map< int, double > GammaFg
Definition tmdbuilder.h:35
std::map< int, std::vector< Set< Operator > > > MatchingFunctionsFFs
Definition tmdbuilder.h:39
std::map< int, double > GammaFq
Definition tmdbuilder.h:34
std::map< int, double > GammaK
Definition tmdbuilder.h:36
std::map< int, std::vector< double > > KCS
Definition tmdbuilder.h:37
std::map< int, std::vector< Set< Operator > > > MatchingFunctionsPDFs
Definition tmdbuilder.h:38
std::map< std::string, std::map< int, double > > HardFactors
Definition tmdbuilder.h:40
double Threshold
Definition tmdbuilder.h:32