Nanga Parbat 1.1.0
A TMD fitting framework
datahandler.h
Go to the documentation of this file.
1//
2// Author: Valerio Bertone: valerio.bertone@cern.ch
3//
4
5#pragma once
6
7#include <string>
8#include <vector>
9#include <utility>
10
11#include <apfel/apfelxx.h>
12#include <yaml-cpp/yaml.h>
13#include <gsl/gsl_rng.h>
14
15namespace NangaParbat
16{
23 {
24 public:
25
29 enum Process: int {UnknownProcess = -1, DY = 0, SIDIS = 1, SIA = 2, DIA = 3};
30
35
41 {
43 bool empty() const;
44 int ndata;
45 double Vs;
46 std::vector<double> qTv;
47 std::vector<std::pair<double, double>> qTmap;
48 std::vector<double> qTfact;
49 std::pair<double, double> var1b;
50 std::pair<double, double> var2b;
51 std::pair<double, double> var3b;
52 bool IntqT;
53 bool Intv1;
54 bool Intv2;
55 bool Intv3;
56 bool PSRed;
57 double pTMin;
58 std::pair<double, double> etaRange;
59 };
60
66 struct Binning
67 {
69 double zmin;
70 double zmax;
71 double zav;
72 bool Intz;
73 double xmin;
74 double xmax;
75 double xav;
76 bool Intx;
77 double Qmin;
78 double Qmax;
79 double Qav;
80 bool IntQ;
81 double ymin;
82 double ymax;
83 double yav;
84 bool Inty;
85 };
86
91
100 DataHandler(std::string const& name, YAML::Node const& datafile, gsl_rng* rng = nullptr, int const& fluctuation = 0, std::vector<double> const& t0 = {});
101
107 void FluctuateData(gsl_rng *rng, int const &fluctuation);
108
116 void SetMeans(std::vector<double> const& means, gsl_rng* rng = nullptr, int const& fluctuation = 0);
117
121 std::string GetName() const { return _name; };
122
126 YAML::Node GetDataFile() const { return _datafile; };
127
131 Process GetProcess() const { return _proc; };
132
136 Observable GetObservable() const { return _obs; };
137
144 double GetTargetIsoscalarity() const { return _targetiso; };
145
150 std::string GetHadron() const { return _hadron; };
151
156 int GetCharge() const { return _charge; };
157
162 std::vector<apfel::QuarkFlavour> GetTagging() const { return _tagging; };
163
168 double GetPrefactor() const { return _prefact; };
169
173 Kinematics GetKinematics() const { return _kin; };
174
178 std::vector<double> GetMeanValues() const { return _means; };
179
183 std::vector<double> GetFluctutatedData() const { return _fluctuations; };
184
189 std::vector<double> GetUncorrelatedUnc() const { return _uncor; };
190
195 std::vector<std::vector<double>> GetAddCorrelatedUnc() const { return _corra; };
196
201 std::vector<std::vector<double>> GetMultCorrelatedUnc() const { return _corrm; };
202
207 std::vector<std::vector<double>> GetCorrelatedUnc() const { return _corr; };
208
213 apfel::matrix<double> GetCovarianceMatrix() const { return _covmat; };
214
219 apfel::matrix<double> GetCholeskyDecomposition() const { return _CholL; };
220
224 std::vector<double> GetT0() const { return _t0; };
225
229 std::map<std::string, std::string> GetLabels() const { return _labels; };
230
237 std::vector<Binning> GetBinning() const { return _bins; }
238
239 protected:
240
241 std::string _name;
242 YAML::Node _datafile;
245 double _targetiso;
246 std::string _hadron;
247 double _charge;
248 std::vector<apfel::QuarkFlavour> _tagging;
249 double _prefact;
251 std::vector<double> _means;
252 std::vector<double> _uncor;
253 std::vector<std::vector<double>> _corra;
254 std::vector<std::vector<double>> _corrm;
255 std::vector<std::vector<double>> _corr;
256 apfel::matrix<double> _covmat;
257 apfel::matrix<double> _CholL;
258 std::map<std::string, std::string> _labels;
259 std::vector<double> _fluctuations;
260 std::vector<double> _t0;
261 std::vector<Binning> _bins;
262
263 friend std::ostream& operator << (std::ostream& os, DataHandler const& DH);
264 };
265
266 std::ostream& operator << (std::ostream &os, DataHandler const& DH);
267}
The "DataHandler" class provides a common interface to all datasets. It provides methods to get kinem...
Definition: datahandler.h:23
std::vector< Binning > _bins
Vector of bins (currently used only for the FF_SIDIS project)
Definition: datahandler.h:261
YAML::Node _datafile
Datafile in YAML.
Definition: datahandler.h:242
std::vector< std::vector< double > > GetAddCorrelatedUnc() const
Function that returns the additive correlated systematic uncertainties.
Definition: datahandler.h:195
Process _proc
The process.
Definition: datahandler.h:243
std::vector< Binning > GetBinning() const
Get vector of bins. This is currently used only for the FF_SIDIS project.
Definition: datahandler.h:237
double _charge
Charge of the identified final state.
Definition: datahandler.h:247
Observable _obs
The observable.
Definition: datahandler.h:244
DataHandler(std::string const &name, YAML::Node const &datafile, gsl_rng *rng=nullptr, int const &fluctuation=0, std::vector< double > const &t0={})
The "DataHandler" constructor.
std::vector< double > _t0
Vector of t0-predictions.
Definition: datahandler.h:260
void FluctuateData(gsl_rng *rng, int const &fluctuation)
Function that fluctuates data.
int GetCharge() const
Function that returns the charge of the identified final state.
Definition: datahandler.h:156
std::string GetName() const
Function that returns the name of the dataset.
Definition: datahandler.h:121
std::vector< double > _means
Vector of central values.
Definition: datahandler.h:251
std::string GetHadron() const
Function that returns the possible identified hadron species in the final state.
Definition: datahandler.h:150
Kinematics GetKinematics() const
Function that returns the kinematic object.
Definition: datahandler.h:173
std::map< std::string, std::string > GetLabels() const
Function that returns the plotting labels.
Definition: datahandler.h:229
std::string _hadron
Hadron species identified in the final state.
Definition: datahandler.h:246
std::vector< double > GetMeanValues() const
Function that returns the mean values.
Definition: datahandler.h:178
double GetTargetIsoscalarity() const
Function that returns the target isoscalarity.
Definition: datahandler.h:144
double _prefact
Possible overall prefactor to multiply the theoretical predictions.
Definition: datahandler.h:249
std::map< std::string, std::string > _labels
Labels used for plotting.
Definition: datahandler.h:258
Kinematics _kin
Kinematics block.
Definition: datahandler.h:250
YAML::Node GetDataFile() const
Function that returns the datafile in YAML format.
Definition: datahandler.h:126
Observable GetObservable() const
Function that returns the observable code.
Definition: datahandler.h:136
apfel::matrix< double > _CholL
Cholesky decomposition of the covariance matrix.
Definition: datahandler.h:257
Process GetProcess() const
Function that returns the process code.
Definition: datahandler.h:131
std::vector< std::vector< double > > _corra
Additive correlated uncertainties.
Definition: datahandler.h:253
friend std::ostream & operator<<(std::ostream &os, DataHandler const &DH)
std::vector< std::vector< double > > _corrm
Multiplicative correlated uncertainties.
Definition: datahandler.h:254
std::vector< std::vector< double > > GetCorrelatedUnc() const
Function that returns the all the correlated systematic uncertainties (additive first and multiplicat...
Definition: datahandler.h:207
double GetPrefactor() const
Function that returns any possible constant prefactor to be used to multiply the theoretical predicti...
Definition: datahandler.h:168
apfel::matrix< double > _covmat
Covariance matrix.
Definition: datahandler.h:256
std::string _name
Name of the dataset.
Definition: datahandler.h:241
std::vector< apfel::QuarkFlavour > GetTagging() const
Function that returns the quark-tagged compoments. Zero corresponds to total.
Definition: datahandler.h:162
std::vector< double > GetFluctutatedData() const
Function that returns the fluctuated data.
Definition: datahandler.h:183
std::vector< std::vector< double > > GetMultCorrelatedUnc() const
Function that returns the multiplicative correlated systematic uncertainties.
Definition: datahandler.h:201
apfel::matrix< double > GetCholeskyDecomposition() const
Function that returns the Cholesky decomposition of the covariance matrix.
Definition: datahandler.h:219
std::vector< double > GetT0() const
Function that returns the set of t0 predictions.
Definition: datahandler.h:224
std::vector< apfel::QuarkFlavour > _tagging
Possible quark-tagged components.
Definition: datahandler.h:248
double _targetiso
Isoscalarity of the target.
Definition: datahandler.h:245
Process
The process enumerator.
Definition: datahandler.h:29
@ SIDIS
Definition: datahandler.h:29
@ SIA
Definition: datahandler.h:29
@ DY
Definition: datahandler.h:29
@ DIA
Definition: datahandler.h:29
@ UnknownProcess
Definition: datahandler.h:29
std::vector< std::vector< double > > _corr
All correlated uncertainties.
Definition: datahandler.h:255
std::vector< double > GetUncorrelatedUnc() const
Function that returns the sum in quadrature of the uncorrelated uncertainties.
Definition: datahandler.h:189
std::vector< double > _fluctuations
Vector of fluctuated data.
Definition: datahandler.h:259
apfel::matrix< double > GetCovarianceMatrix() const
Function that returns the covariance matrix of the correlated uncertainties.
Definition: datahandler.h:213
DataHandler(DataHandler const &DH)
The "DataHandler" copy constructor.
void SetMeans(std::vector< double > const &means, gsl_rng *rng=nullptr, int const &fluctuation=0)
Function that sets the data central values replacing that introduced in the constructor.
std::vector< double > _uncor
Vector of uncorrelated uncertainties.
Definition: datahandler.h:252
Observable
The observable enumerator.
Definition: datahandler.h:34
@ opposite_sign_ratio
Definition: datahandler.h:34
@ dsigma_dxdQdz
Definition: datahandler.h:34
@ UnknownObservable
Definition: datahandler.h:34
@ dsigma_dxdydz
Definition: datahandler.h:34
@ multiplicity
Definition: datahandler.h:34
@ F_uut
Definition: datahandler.h:34
Definition: bstar.h:12
YAML::Emitter & operator<<(YAML::Emitter &os, ChiSquare const &chi2)
Method which prints ChiSquare feautures with cout <<.
Structure containing the single bin information. This is currently used only for the FF_SIDIS project...
Definition: datahandler.h:67
double Qav
Definition: datahandler.h:79
double xmax
Definition: datahandler.h:74
bool Intz
Definition: datahandler.h:72
double xav
Definition: datahandler.h:75
double zmin
Definition: datahandler.h:69
bool IntQ
Definition: datahandler.h:80
bool Inty
Definition: datahandler.h:84
double zav
Definition: datahandler.h:71
double Qmin
Definition: datahandler.h:77
double zmax
Definition: datahandler.h:70
double ymax
Definition: datahandler.h:82
double xmin
Definition: datahandler.h:73
double Qmax
Definition: datahandler.h:78
double ymin
Definition: datahandler.h:81
bool Intx
Definition: datahandler.h:76
double yav
Definition: datahandler.h:83
Structure containing the kinematic information of one single data set.
Definition: datahandler.h:41
double Vs
Center of mass energy.
Definition: datahandler.h:45
std::pair< double, double > etaRange
Allowed range in eta of the final-state leptons.
Definition: datahandler.h:58
std::vector< std::pair< double, double > > qTmap
Map of qT bounds to associate to the single bins.
Definition: datahandler.h:47
double pTMin
Minimum pT of the final-state leptons.
Definition: datahandler.h:57
std::pair< double, double > var2b
Variable 2 integration bounds.
Definition: datahandler.h:50
bool Intv1
Whether the bins variable 1 are integrated over.
Definition: datahandler.h:53
std::vector< double > qTv
Vector of qT values.
Definition: datahandler.h:46
bool PSRed
Whether there is a final-state PS reduction.
Definition: datahandler.h:56
std::pair< double, double > var1b
Variable 1 integration bounds.
Definition: datahandler.h:49
bool Intv2
Whether the bins variable 2 are integrated over.
Definition: datahandler.h:54
int ndata
Number of data points.
Definition: datahandler.h:44
bool IntqT
Whether the bins in qTv are integrated over.
Definition: datahandler.h:52
std::vector< double > qTfact
Possible bin-by-bin prefactors to multiply the theoretical predictions.
Definition: datahandler.h:48
bool Intv3
Whether the bins variable 3 are integrated over.
Definition: datahandler.h:55
std::pair< double, double > var3b
Variable 3 integration bounds.
Definition: datahandler.h:51