Nanga Parbat 1.1.0
A TMD fitting framework
createtmdgrid.h
Go to the documentation of this file.
1//
2// Author: Valerio Bertone: valerio.bertone@cern.ch
3//
4
5#pragma once
6
8
9#include <vector>
10#include <memory>
11#include <yaml-cpp/yaml.h>
12#include <apfel/apfelxx.h>
13
14namespace NangaParbat
15{
20 {
21 std::vector<double> Qg;
22 std::vector<double> xg;
23 std::vector<double> qToQg;
24 };
25
31 ThreeDGrid Inter3DGrid(std::string const& pf)
32 {
33 ThreeDGrid grid{};
34 if (pf == "pdf")
35 {
36 grid.Qg = std::vector<double>
37 {
38 1.000000e+00, 1.118034e+00, 1.224745e+00, 1.400000e+00, 1.581139e+00, 1.788854e+00,
39 2.000000e+00, 2.236068e+00, 2.529822e+00, 2.828427e+00, 3.162278e+00, 3.464102e+00,
40 4.750000e+00, 5.099020e+00, 6.324555e+00, 7.100000e+00, 8.000000e+00, 1.000000e+01,
41 1.118034e+01, 1.224745e+01, 1.400000e+01, 1.581139e+01, 1.788854e+01, 2.000000e+01,
42 2.236068e+01, 2.529822e+01, 2.828427e+01, 3.162278e+01, 3.464102e+01, 4.750000e+01,
43 5.099020e+01, 6.324555e+01, 7.100000e+01, 8.000000e+01, 1.000000e+02, 1.118034e+02,
44 1.224745e+02, 1.400000e+02, 1.581139e+02, 1.788854e+02, 2.000000e+02
45 };
46 grid.xg = std::vector<double>
47 {
48 1.000000e-05, 2.000000e-05, 4.000000e-05, 6.000000e-05, 8.000000e-05,
49 1.000000e-04, 2.000000e-04, 4.000000e-04, 6.000000e-04, 8.000000e-04,
50 1.000000e-03, 1.500000e-03, 2.000000e-03, 2.500000e-03, 3.000000e-03,
51 3.500000e-03, 4.000000e-03, 4.500000e-03, 5.000000e-03, 5.500000e-03,
52 6.000000e-03, 6.500000e-03, 7.000000e-03, 7.500000e-03, 8.000000e-03,
53 8.500000e-03, 9.000000e-03, 9.250000e-03, 9.500000e-03, 9.750000e-03,
54 1.000000e-02, 1.500000e-02, 2.000000e-02, 2.500000e-02, 3.000000e-02,
55 3.500000e-02, 4.000000e-02, 4.500000e-02, 5.000000e-02, 5.500000e-02,
56 6.000000e-02, 6.500000e-02, 7.000000e-02, 7.500000e-02, 8.000000e-02,
57 8.500000e-02, 9.000000e-02, 9.250000e-02, 9.500000e-02, 9.750000e-02,
58 1.000000e-01, 1.500000e-01, 2.000000e-01, 2.500000e-01, 3.000000e-01,
59 3.500000e-01, 4.000000e-01, 4.500000e-01, 5.000000e-01, 5.500000e-01,
60 6.000000e-01, 6.500000e-01, 7.000000e-01, 7.500000e-01, 8.000000e-01,
61 8.500000e-01, 9.000000e-01, 9.250000e-01, 9.500000e-01, 9.750000e-01,
62 1.000000e+00
63 };
64 grid.qToQg = std::vector<double>
65 {
66 0.0001, 0.0010, 0.0025, 0.0050, 0.0075, 0.0100, 0.0200, 0.0300, 0.0400,
67 0.0500, 0.0600, 0.0700, 0.0800, 0.0900, 0.1000, 0.1250, 0.1500, 0.1750,
68 0.2000, 0.2250, 0.2500, 0.2750, 0.3000, 0.3500, 0.4000, 0.4500, 0.5000,
69 0.5500, 0.6000, 0.6500, 0.7000, 0.8000, 0.9000, 1,
70 1.1000, 1.2000, 1.3000, 1.4000, 1.5000, 1.6000, 1.7000, 1.8000, 1.9000,
71 2.0010
72 };
73 }
74 else
75 {
76 grid.Qg = std::vector<double>
77 {
78 1.000000e+00, 1.080000e+00, 1.118034e+00, 1.170000e+00, 1.224745e+00, 1.300000e+00,
79 1.400000e+00, 1.581139e+00, 1.788854e+00, 2.000000e+00, 2.236068e+00, 2.529822e+00,
80 2.828427e+00, 3.162278e+00, 3.464102e+00, 4.750000e+00, 5.099020e+00, 6.324555e+00,
81 7.100000e+00, 8.000000e+00, 1.000000e+01, 1.118034e+01, 1.224745e+01, 1.400000e+01,
82 1.581139e+01, 1.788854e+01, 2.000000e+01, 2.236068e+01, 2.529822e+01, 2.828427e+01,
83 3.162278e+01, 3.464102e+01, 4.750000e+01, 5.099020e+01, 6.324555e+01, 7.100000e+01,
84 8.000000e+01, 1.000000e+02
85 };
86 grid.xg = std::vector<double>
87 {
88 1.000000e-01, 1.250000e-01, 1.500000e-01, 1.750000e-01, 2.000000e-01, 2.250000e-01,
89 2.500000e-01, 2.750000e-01, 3.000000e-01, 3.250000e-01, 3.500000e-01, 3.750000e-01,
90 4.000000e-01, 4.250000e-01, 4.500000e-01, 4.750000e-01, 5.000000e-01, 5.250000e-01,
91 5.500000e-01, 5.750000e-01, 6.000000e-01, 6.250000e-01, 6.500000e-01, 7.000000e-01,
92 7.500000e-01, 8.000000e-01, 8.500000e-01, 9.000000e-01, 9.500000e-01, 1.000000e+00
93 };
94 grid.qToQg = std::vector<double>
95 {
96 0.0001,
97 0.0010, 0.0025, 0.0050, 0.0075, 0.0100, 0.0200, 0.0300, 0.0400,
98 0.0500, 0.0600, 0.0700, 0.0800, 0.0900, 0.1000, 0.1100, 0.1200,
99 0.1300, 0.1400, 0.1500, 0.1600, 0.1700, 0.1800, 0.1900, 0.2000,
100 0.2200, 0.2400, 0.2600, 0.2800, 0.3000, 0.3200, 0.3400, 0.3600,
101 0.3800, 0.4000, 0.4500, 0.5000, 0.5500, 0.6000, 0.6500, 0.7000,
102 0.8000, 0.9000, 1,
103 1.1000, 1.2000, 1.3000, 1.4000, 1.5000, 1.6000, 1.7000, 1.8000,
104 1.9000, 2,
105 2.1000, 2.2000, 2.3000, 2.4000, 2.5000, 2.6000, 2.7000, 2.8000,
106 2.9000, 3,
107 3.2000, 3.4000, 3.6000, 3.8000, 4,
108 4.2000, 4.4000, 4.6000, 4.8000, 5
109 };
110 }
111 return grid;
112 };
113
117 const std::map<int, std::string> PtOrderMap{{0, "LL"}, {1, "NLL"}, {-1, "NLL'"}, {2, "NNLL"}, {-2, "NNLL'"}, {3, "NNNLL"}};
118
127 void ProduceTMDGrid(std::string const& ReportFolder, std::string const& Output, std::string const& distype = "pdf");
128
140 std::unique_ptr<YAML::Emitter> EmitTMDGrid(YAML::Node const& config,
141 std::string const& parameterisation,
142 std::vector<double> const& params,
143 std::string const& pf,
144 ThreeDGrid const& tdg);
145
166 std::unique_ptr<YAML::Emitter> EmitTMDInfo(YAML::Node const& config,
167 int const& NumMembers,
168 std::string const& pf,
169 ThreeDGrid const& tdg,
170 std::vector<int> const& Flavors = {-5, -4, -3, -2, -1, 1, 2, 3, 4, 5},
171 std::string const& SetDesc = "Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)",
172 std::string const& Authors = "M.A.P. Collaboration",
173 std::string const& Reference = "arXiv:xxxx.xxxx",
174 std::string const& SetIndex = "000000",
175 std::string const& SetName = "PV17nll",
176 std::string const& Format = "TMDlib2",
177 std::string const& DataVersion = "1",
178 std::string const& ErrorType = "Monte Carlo",
179 std::string const& FlavorScheme = "LHAPDF style");
180}
Definition: bstar.h:12
const std::map< int, std::string > PtOrderMap
Map of perturbative orders.
Definition: createtmdgrid.h:117
ThreeDGrid Inter3DGrid(std::string const &pf)
Function that returns ThreeDGrid object better tuned for PDFs or FFs, according to the input string.
Definition: createtmdgrid.h:31
void ProduceTMDGrid(std::string const &ReportFolder, std::string const &Output, std::string const &distype="pdf")
This function encapsulates and streamlines the production of an interpolation grid starting from the ...
std::unique_ptr< YAML::Emitter > EmitTMDGrid(YAML::Node const &config, std::string const &parameterisation, std::vector< double > const &params, std::string const &pf, ThreeDGrid const &tdg)
Function that produces the TMD interpolation grid in momentum space. This is supposed to resamble an ...
std::unique_ptr< YAML::Emitter > EmitTMDInfo(YAML::Node const &config, int const &NumMembers, std::string const &pf, ThreeDGrid const &tdg, std::vector< int > const &Flavors={-5, -4, -3, -2, -1, 1, 2, 3, 4, 5}, std::string const &SetDesc="Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)", std::string const &Authors="M.A.P. Collaboration", std::string const &Reference="arXiv:xxxx.xxxx", std::string const &SetIndex="000000", std::string const &SetName="PV17nll", std::string const &Format="TMDlib2", std::string const &DataVersion="1", std::string const &ErrorType="Monte Carlo", std::string const &FlavorScheme="LHAPDF style")
Function that produces the info file of the TMD set. This is suppose to resamble an LHAPDF info file ...
Structure that contains the 3D grid in Q, x, and qT.
Definition: createtmdgrid.h:20
std::vector< double > Qg
Grid in Q.
Definition: createtmdgrid.h:21
std::vector< double > qToQg
Grid in qT/Q.
Definition: createtmdgrid.h:23
std::vector< double > xg
Grid in x.
Definition: createtmdgrid.h:22