Nanga Parbat 1.1.0
A TMD fitting framework
createstructgrid.h
Go to the documentation of this file.
1//
2// Author: Valerio Bertone: valerio.bertone@cern.ch
3// Chiara Bissolotti: chiara.bissolotti01@universitadipavia.it
4//
5
6#pragma once
7
9
10#include <vector>
11#include <memory>
12#include <yaml-cpp/yaml.h>
13#include <apfel/apfelxx.h>
14
15namespace NangaParbat
16{
20 struct FourDGrid
21 {
22 std::vector<double> Qg;
23 std::vector<double> xg;
24 std::vector<double> zg;
25 std::vector<double> qToQg;
26 };
27
33 FourDGrid Inter4DGrid(std::string const& pf)
34 {
35 FourDGrid grid{};
36 grid.Qg = std::vector<double>
37 {
38 1.000000e+00, 1.080000e+00, 1.118030e+00, 1.170000e+00, 1.224740e+00,
39 1.300000e+00, 1.400000e+00, 1.581140e+00, 1.788850e+00, 2.000000e+00,
40 2.236070e+00, 2.529820e+00, 2.828430e+00, 3.162280e+00, 3.464100e+00,
41 4.750000e+00, 5.099020e+00, 6.324560e+00, 7.100000e+00, 8.000000e+00,
42 1.000000e+01, 1.118030e+01, 1.224750e+01, 1.400000e+01, 1.581140e+01,
43 1.788850e+01, 2.000000e+01, 2.236070e+01, 2.529820e+01, 2.828430e+01,
44 3.162280e+01, 3.464100e+01, 4.750000e+01, 5.099020e+01, 6.324560e+01,
45 7.100000e+01, 8.000000e+01, 1.000000e+02
46 };
47 grid.xg = std::vector<double>
48 {
49 1.000000e-05, 2.000000e-05, 4.000000e-05, 6.000000e-05, 8.000000e-05,
50 1.000000e-04, 2.000000e-04, 4.000000e-04, 6.000000e-04, 8.000000e-04,
51 1.000000e-03, 1.500000e-03, 2.000000e-03, 2.500000e-03, 3.000000e-03,
52 3.500000e-03, 4.000000e-03, 4.500000e-03, 5.000000e-03, 5.500000e-03,
53 6.000000e-03, 6.500000e-03, 7.000000e-03, 7.500000e-03, 8.000000e-03,
54 8.500000e-03, 9.000000e-03, 9.250000e-03, 9.500000e-03, 9.750000e-03,
55 1.000000e-02, 1.500000e-02, 2.000000e-02, 2.500000e-02, 3.000000e-02,
56 3.500000e-02, 4.000000e-02, 4.500000e-02, 5.000000e-02, 5.500000e-02,
57 6.000000e-02, 6.500000e-02, 7.000000e-02, 7.500000e-02, 8.000000e-02,
58 8.500000e-02, 9.000000e-02, 9.250000e-02, 9.500000e-02, 9.750000e-02,
59 1.000000e-01, 1.500000e-01, 2.000000e-01, 2.500000e-01, 3.000000e-01,
60 3.500000e-01, 4.000000e-01, 4.500000e-01, 5.000000e-01, 5.500000e-01,
61 6.000000e-01, 6.500000e-01, 7.000000e-01, 7.500000e-01, 8.000000e-01
62 };
63 grid.zg = std::vector<double>
64 {
65 0.100, 0.125, 0.150, 0.175,
66 0.200, 0.225, 0.250, 0.275,
67 0.300, 0.325, 0.350, 0.375,
68 0.400, 0.425, 0.450, 0.475,
69 0.500, 0.525, 0.550, 0.575,
70 0.600, 0.625, 0.650, 0.700,
71 0.750, 0.80
72 };
73 grid.qToQg =
74 {
75 0.0001, 0.0010, 0.0025, 0.0050, 0.0075, 0.0100, 0.0200, 0.0300, 0.0400,
76 0.0500, 0.0600, 0.0700, 0.0800, 0.0900, 0.1000, 0.1100, 0.1200, 0.1300,
77 0.1400, 0.1500, 0.1600, 0.1700, 0.1800, 0.1900, 0.2000, 0.2200, 0.2400,
78 0.2600, 0.2800, 0.3000, 0.3200, 0.3400, 0.3600, 0.3800, 0.4000, 0.4500,
79 0.5000, 0.5500, 0.6000, 0.6500, 0.7000, 0.8000, 0.9000, 1
80 };
81
82 return grid;
83 };
84
85 // /**
86 // * @brief Map of perturbative orders
87 // */
88 // const std::map<int, std::string> PtOrderMapSF{{0, "LL"}, {1, "NLL"}, {-1, "NLL'"}, {2, "NNLL"}, {-2, "NNLL'"}, {3, "NNNLL"}};
89
101 void ProduceStructGrid(std::string const& GridsDirectory,
102 std::string const& GridTMDPDFfolder,
103 std::string const& GridTMDFFfolder,
104 std::string const& Output,
105 std::string const& repID = "none",
106 std::string const& structype = "FUUT");
107
121 std::unique_ptr<YAML::Emitter> EmitStructGrid(std::string const& GridsDirectory,
122 std::string const& GridTMDPDFfolder,
123 std::string const& GridTMDFFfolder,
124 int const& repnumber,
125 std::string const& pf,
126 FourDGrid const& fdg,
127 int const& qToQcut = 5);
128
140 std::unique_ptr<YAML::Emitter> EmitStructGridDirect(std::string const& FitDirectory,
141 int const& repnumber,
142 std::string const& pf,
143 FourDGrid const& fdg,
144 int const& qToQcut);
145
169 std::unique_ptr<YAML::Emitter> EmitStructInfo(std::string const& GridsDirectory,
170 std::string const& GridTMDPDFfolder,
171 std::string const& GridTMDFFfolder,
172 YAML::Node const& config,
173 int const& NumMembers,
174 std::string const& pf,
175 FourDGrid const& fdg,
176 std::vector<int> const& Flavors = {-5, -4, -3, -2, -1, 1, 2, 3, 4, 5},
177 std::string const& SetDesc = "Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)",
178 std::string const& Target = "proton",
179 std::string const& Hadron = "Pip",
180 std::string const& Authors = "A. Bacchetta, F. Delcarro, C. Pisano, M. Radici, A. Signori",
181 std::string const& Reference = "arXiv:1703.10157",
182 // std::string const& Authors = "A. Bacchetta, V. Bertone, C. Bissolotti, G. Bozzi, F. Delcarro, F. Piacenza, M. Radici",
183 // std::string const& Reference = "arXiv:1912.07550",
184 std::string const& SetIndex = "000000",
185 std::string const& Format = "TMDlib2",
186 std::string const& DataVersion = "1",
187 std::string const& ErrorType = "Monte Carlo");
188}
Definition: bstar.h:12
void ProduceStructGrid(std::string const &GridsDirectory, std::string const &GridTMDPDFfolder, std::string const &GridTMDFFfolder, std::string const &Output, std::string const &repID="none", std::string const &structype="FUUT")
This function encapsulates and streamlines the production of an interpolation grid for a structure fu...
std::unique_ptr< YAML::Emitter > EmitStructGrid(std::string const &GridsDirectory, std::string const &GridTMDPDFfolder, std::string const &GridTMDFFfolder, int const &repnumber, std::string const &pf, FourDGrid const &fdg, int const &qToQcut=5)
Function that produces the structure function interpolation grid in momentum space....
FourDGrid Inter4DGrid(std::string const &pf)
Function that returns FourDGrid object better tuned for F_UUT, according to the input string.
Definition: createstructgrid.h:33
std::unique_ptr< YAML::Emitter > EmitStructGridDirect(std::string const &FitDirectory, int const &repnumber, std::string const &pf, FourDGrid const &fdg, int const &qToQcut)
Function that produces the structure function interpolation grid in momentum space....
std::unique_ptr< YAML::Emitter > EmitStructInfo(std::string const &GridsDirectory, std::string const &GridTMDPDFfolder, std::string const &GridTMDFFfolder, YAML::Node const &config, int const &NumMembers, std::string const &pf, FourDGrid const &fdg, 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 &Target="proton", std::string const &Hadron="Pip", std::string const &Authors="A. Bacchetta, F. Delcarro, C. Pisano, M. Radici, A. Signori", std::string const &Reference="arXiv:1703.10157", std::string const &SetIndex="000000", std::string const &Format="TMDlib2", std::string const &DataVersion="1", std::string const &ErrorType="Monte Carlo")
Function that produces the info file of the TMD set. This is suppose to resamble an LHAPDF info file ...
Structure that contains the 4D grid in Q, x, z and qT.
Definition: createstructgrid.h:21
std::vector< double > zg
Grid in z.
Definition: createstructgrid.h:24
std::vector< double > xg
Grid in x.
Definition: createstructgrid.h:23
std::vector< double > qToQg
Grid in qT/Q.
Definition: createstructgrid.h:25
std::vector< double > Qg
Grid in Q.
Definition: createstructgrid.h:22