LTF
LTF_ROOTTools.h
1 //#include "include/LTF/LTF.h"
2 #include "LTF/LTF.h"
3 
4 #include <TROOT.h>
5 #include <TStyle.h>
6 #include <TLegend.h>
7 #include <TMatrixDfwd.h>
8 #include <TMatrixD.h>
9 #include <TMatrixTSym.h>
10 #include <TCanvas.h>
11 #include <TH1D.h>
12 #include <TSystem.h>
13 #include <TFile.h>
14 #include <string>
15 #include <fstream>
16 #include <TRandom3.h>
17 #include <TGraphErrors.h>
18 #include <TLatex.h>
19 #include <TF1.h>
20 #include <TF2.h>
21 #include <TLine.h>
22 #include <TRandom3.h>
23 #include <TGraph2DErrors.h>
24 #include <TGraph2D.h>
25 
26 using namespace std;
27 
29 public:
30 
31 // __________________________________________________________________________________ //
37 static
38 std::map < std::string, std::vector<double> > read_input_table2(std::string filename, int ncol ) {
39  std::map < std::string, std::vector<double> > ret;
40  std::vector<std::string> cols;
41  // open file for reading
42  std::ifstream istrm(filename.c_str(), std::ios::binary);
43  if (!istrm.is_open()) {
44  std::cout << "failed to open " << filename << std::endl;
45  exit(1);
46  }
47  for ( int c=0 ; c<ncol ; c++ ) {
48  std::string colname;
49  istrm >> colname;
50  ret[colname] = vector<double>();
51  cols.push_back(colname);
52  }
53  while ( istrm.good()) {
54  double value;
55  for ( int c=0 ; c<ncol ; c++ ) {
56  istrm >> value;
57  if ( !istrm.good() ) break;
58  std::string colname = cols[c];
59  ret[colname].push_back(value);
60  }
61  }
62  cout<<"Info. [read_input_table] Read "<<ret.size()<<" rows."<<endl;
63  return ret;
64 }
65 
66 
67 
68 // __________________________________________________________________________________ //
74 static
75 TGraphErrors* MakeTGraph(const Eigen::VectorXd& xvalues, int ibin, const Eigen::MatrixXd& Y,
76  const std::map<std::string,Eigen::MatrixXd >& VSysY = {}
77  ) {
78 
79  TGraphErrors* graph = new TGraphErrors();
80  for ( int i = 0 ; i<xvalues.size() ; i++ ) {
81  double yvalue = Y(ibin,i);
82  graph->SetPoint(i,xvalues(i),yvalue);
83  double ey2 = 0;
84  for ( auto [name,Vy] : VSysY ) {
85  ey2 += pow(Vy(ibin,i),2);
86  }
87  graph->SetPointError(i,0,sqrt(ey2));
88  }
89  return graph;
90 }
91 
92 
93 // __________________________________________________________________________________ //
99 static
100 TGraph2DErrors* MakeTGraph2D(const Eigen::VectorXd& xvalues, const Eigen::VectorXd& yvalues,
101  int ibin, const Eigen::MatrixXd& Y,
102  const std::map<std::string,Eigen::MatrixXd >& VSysY = {}
103  ) {
104 
105  TGraph2DErrors* graph = new TGraph2DErrors();
106  for ( int i = 0 ; i<xvalues.size() ; i++ ) {
107  double zvalue = Y(ibin,i);
108  graph->SetPoint(i,xvalues(i),yvalues(i),zvalue);
109  double ey2 = 0;
110  for ( auto [name,Vy] : VSysY ) {
111  ey2 += pow(Vy(ibin,i),2);
112  }
113  graph->SetPointError(i,0,0,sqrt(ey2));
114  }
115  return graph;
116 }
117 
118 
119 // __________________________________________________________________________________ //
126 static
127 TH1D* MakeHistogram(int nEvents, int seed, double mean, double sigma, vector<double> bins ) {
128  TRandom3 rn(seed);
129 
130  TH1D* hist = new TH1D("hist","hist",bins.size()-1, &bins[0]);
131  hist->Sumw2();
132  for ( int i=0 ; i<nEvents; i++ ) {
133  hist->Fill( rn.Gaus(mean,sigma), 1 );
134  }
135  return hist;
136 }
137 
138 
139 
140 
141 // __________________________________________________________________________________ //
147 static
148 TH1D* MakeHistogram(const Eigen::VectorXd& values, vector<double> bins ={}, const std::vector<std::pair<std::string,Eigen::MatrixXd > >& V = {} ) {
149  TH1D* hist = bins.empty() ?
150  new TH1D("hist","hist",values.size(),0,values.size() ) :
151  new TH1D("hist","hist",bins.size()-1, &bins[0]);
152  if ( bins.size() && values.size() != bins.size()-1 ) {cout<<"ERROR! binning and number of entries does not fit!"<<endl;exit(1);}
153  for ( size_t i = 0 ;i<bins.size()-1; i++ ) {
154  hist->SetBinContent(i+1, values(i));
155  if ( V.size() ) {
156  double e2 = 0;
157  for ( auto apair : V ) {
158  e2 += apair.second(i,i);
159  }
160  hist->SetBinError(i+1, sqrt(e2) );
161  }
162  }
163  return hist;
164 }
165 
166 
167 
168 // __________________________________________________________________________________ //
176 static
177 void plotLiTeFit(const LTF::LiTeFit& fit, const vector<double>& bins,
178  const string& yaxistitle = "value [unit]",
179  const string& referencename = "Reference value (#alpha) [unit]",
180  const string& observablename = "Observable [unit]"
181  ) {
182  gStyle->SetOptStat(0);
183  gSystem->mkdir("plots");
184  auto& M = fit.M;
185 
186  // sanity check
187  if ( M.cols() != 2 ) {cout<<"Error! only 1-dim plotting implemented."<<endl;exit(1);}
188  Eigen::VectorXd reference_values = M.col(1);
189  //TH1D* hist = new TH1D("hist","hist",bins.size()-1, &bins[0]);
190 
191  map<double,TH1D*> templates;
192  for ( int iref = 0 ; iref<reference_values.size() ; iref++ ) {
193  //double ref = reference_values(iref);
194  templates[iref] = MakeHistogram(fit.Y.col(iref),bins);
195  }
196 
197  TH1D* data = MakeHistogram(fit.Dt,bins,fit.Vs);
198  TH1D* TheoFit = MakeHistogram(fit.TheoFit,bins);
199 
200  TCanvas c1("c1","LTF plots",800,800);
201  c1.SetRightMargin(0.05);
202  c1.SetLeftMargin(0.15);
203  c1.SetTopMargin(0.08);
204  // c1.SetRightMargin(0.02);
205 
206  const char* ps_name = fit.GetLogNormal() ?
207  "LTFlog_plots.ps" :
208  "LTF_plots.ps";
209  c1.Print( (string(ps_name)+"[").c_str() );
210 
211  // ---------------------------------------------- //
212  // main plot
213  // ---------------------------------------------- //
214  for ( int iref = 0 ; iref<reference_values.size() ; iref++ ) {
215  templates[iref]->SetLineWidth(2);
216  if ( iref == 0 ) {
217  templates[0]->SetTitle(("Linear Template Fit;"+observablename+";"+yaxistitle).c_str());
218  templates[0]->SetLineColor(kRed+1);
219  if ( templates[0]->GetMaximum()>0 )templates[0]->SetMinimum(0);
220  if ( !fit.GetLogNormal() )
221  templates[0]->SetMaximum(templates[0]->GetMaximum()*1.7);
222  // else
223  // templates[0]->SetMaximum(templates[0]->GetMaximum()*3.);
224  templates[0]->SetLineWidth(3);
225  templates[0]->DrawClone("hist");
226  }
227  else if ( iref==reference_values.size()-1) {
228  if ( templates[0]->GetMaximum()>0 ) templates[iref]->SetFillColorAlpha(kBlue,0.15);
229  templates[iref]->SetLineColor(kBlue+2);
230  templates[iref]->SetLineWidth(3);
231  templates[iref]->Draw("histsame");
232  }
233  else {
234  templates[iref]->SetLineColor(iref+2);
235  templates[iref]->SetLineWidth(2);
236  templates[iref]->Draw("histsame");
237  }
238  }
239  if ( templates[0]->GetMaximum()>0 )templates[0]->SetFillColorAlpha(kRed,0.15);
240  templates[0]->Draw("histsame");
241 
242  data->SetMarkerStyle(20);
243  data->SetMarkerSize(1.4);
244  data->SetLineColor(kBlack);
245  data->Draw("e0same");
246 
247  TheoFit->SetLineColor(923);
248  TheoFit->SetLineWidth(4);
249  TheoFit->SetLineStyle(2);
250  TheoFit->Draw("histsame");
251 
252  TLegend legend(0.18,0.70,0.94,0.92,"","NDC");
253  legend.SetNColumns(3);
254  legend.SetFillStyle(0);
255  legend.SetBorderSize(0);
256  legend.AddEntry(data,"Data","E0P");
257  for ( int iref = 0 ; iref<reference_values.size() ; iref++ ) {
258  legend.AddEntry(templates[iref],Form("Template #alpha=%6.2f",reference_values[iref]),"FL");
259  }
260  legend.AddEntry(TheoFit,"Estimated best model","L");
261  legend.Draw();
262 
263  c1.Print(ps_name);
264  c1.Print("plots/LTF_plot.pdf");
265 
266 
267  // ---------------------------------------------- //
268  // print linear-functions in every bin
269  // ---------------------------------------------- //
270  c1.Clear();
271  c1.SetRightMargin(0.02);
272  c1.SetTopMargin(0.02);
273  c1.SetLeftMargin(0.16);
274  c1.SetBottomMargin(0.16);
275 
276  gPad->SetTicky(1);
277 
278  gStyle->SetLabelSize(0.05,"XYZ");
279  gStyle->SetTitleSize(0.05,"XYZ");
280  gStyle->SetTitleOffset(1.1,"X");
281  gStyle->SetTitleOffset(1.6,"Y");
282  gStyle->SetMarkerSize(2);
283 
284 
285  map < string, vector<double> > input_table = read_input_table2("data/CMS_data.txt",32);
286 
287  for ( int ibin = 0 ; ibin<fit.Dt.size() ; ibin++ ) {
288  //for ( int ibin = 0 ; ibin<0 ; ibin++ ) {
289  //for ( int ibin = 0 ; ibin<6 ; ibin++ ) {
290  //TGraphErrors* data = MakeTGraph(fit.ahat.row(0),fit.Dt.row(ibin));
291  TGraphErrors* gdata = new TGraphErrors();
292  gdata->SetPoint(0,fit.ahat(0),fit.Dt(ibin));
293  gdata->SetPointError(0,0,data->GetBinError(ibin+1));
294  gdata->SetMarkerStyle(20);
295 
296  TGraphErrors* graph = MakeTGraph(fit.M.col(1),ibin,fit.Y,fit.SysY);
297  graph->Fit("pol1","Q");
298  TF1* pol1w = (TF1*)graph->GetFunction("pol1")->Clone("pol1w");
299  pol1w->SetLineColor(922);
300  pol1w->SetLineStyle(3);
301  pol1w->SetLineWidth(2);
302 
303  TF1* f2= new TF1("f2","[0]+[1]*pow(x,2)", -FLT_MIN,FLT_MAX );
304  graph->Fit(f2,"QW");
305  TF1* f1= new TF1("f1","[0]+[1]*pow(x,[2])", -FLT_MIN,FLT_MAX );
306 
307 
308  graph->Fit("pol1","QW"); // "W": Ignore all point errors when fitting a TGraphErrors
309  bool UseLTFOutput = true;
310  if ( UseLTFOutput ) {
311  Eigen::VectorXd ltfpol1param =(fit.Y*fit.Mc().transpose()).row(ibin);
312  //cout<<"M+*Y:"<<endl<< ltfpol1param <<endl;
313  graph->GetFunction("pol1")->SetParameter(0,ltfpol1param(0));
314  graph->GetFunction("pol1")->SetParameter(1,ltfpol1param(1));
315  f1->SetParameter(0,ltfpol1param(0));
316  f1->SetParameter(1,ltfpol1param(1));
317  f1->SetParameter(2,fit.Gamma[0]);
318  }
319  if ( fit.Gamma[0]!=1 ) cout<<"Warning! Plotting with gamma factor !=1 not correctly implmeneted!"<<endl;
320 
321  graph->GetFunction("pol1")->SetLineWidth(3);
322 
323  graph->Fit("pol2","QW");
324  graph->GetFunction("pol2")->SetLineWidth(83);
325 
326  graph->SetMarkerStyle(47);
327  graph->SetMarkerColor(kRed+2);
328  graph->SetLineColor(kRed+2);
329  if ( fit.GetLogNormal() )
330  graph->SetTitle((";"+referencename+";log("+yaxistitle+")").c_str()); // log(value/unit)
331  else
332  graph->SetTitle((";"+referencename+";"+yaxistitle).c_str());
333 
334 
335  TF1* f1log = NULL;
336  if ( fit.GetLogNormal() ) {
337  //f1log = new TF1("pol1log","[0]+[1]*exp(x)", -FLT_MIN,FLT_MAX );
338  f1log = new TF1("pol1log","log([0]+[1]*pow(x,1))", -FLT_MIN,FLT_MAX );
339 
340  TGraph* gexp = new TGraph();
341  for ( int i = 0 ; i<graph->GetN() ; i++ )
342  gexp->SetPoint(i,graph->GetX()[i], exp(graph->GetY()[i]));
343  gexp->Fit("pol1","QW");
344 
345  Eigen::VectorXd ltfpol1param =(fit.Y*fit.Mc().transpose()).row(ibin);
346  graph->GetFunction("pol1")->SetParameter(0,ltfpol1param(0));
347 
348  f1log->SetParameter(0, gexp->GetFunction("pol1")->GetParameter(0));
349  f1log->SetParameter(1, gexp->GetFunction("pol1")->GetParameter(1));
350  if ( fit.Gamma[0]!=1 ) cout<<"Warning! Plotting with gamma factor !=1 not correctly implmeneted!"<<endl;
351 
352  // //f1 = new TF1("pol1","[0]+[1]*x", -FLT_MIN,FLT_MAX );
353  // graph->Fit(f1log,"QW"); // "W": Ignore all point errors when fitting a TGraphErrors
354  f1log->SetLineColor(kBlue+1);
355  f1log->SetLineStyle(7);
356  f1log->SetLineWidth(2);
357  }
358 
359 
360  if ( gdata->GetY()[0] > 0 )
361  graph->SetMinimum(0);
362 
363  graph->SetMaximum( max(gdata->GetY()[0],max(graph->GetY()[0],graph->GetY()[graph->GetN()-1]))*1.2);
364 
365  graph->Draw("ap");
366  if ( fit.GetLogNormal() )
367  f1log->Draw("Lsame");
368  else
369  pol1w->Draw("Lsame");
370  //graph->GetFunction("pol1")->Draw("Lsame");
371  if ( reference_values.size()+1<= 6 )
372  graph->GetHistogram()->GetXaxis()->SetNdivisions(graph->GetN()+1);
373  else
374  graph->GetHistogram()->GetXaxis()->SetNdivisions(int(graph->GetN()/2)+1+200);
375  //f2->Draw("same");
376 
377  gdata->Draw("P");
378 
379  if ( ibin==0 ) {
380  double xmin = fit.GetLogNormal() ? 0.36 : 0.45;
381  TLegend legend(xmin,0.19,0.96,0.47,"","NDC");
382  //legend.SetNColumns(3);
383  legend.SetFillStyle(0);
384  legend.SetBorderSize(0);
385  legend.SetTextSize(0.05);
386  legend.AddEntry(data,"Data","E0P");
387  legend.AddEntry(graph,"Templates","PE0");
388  if ( fit.GetLogNormal() ) {
389  legend.AddEntry(graph->GetFunction("pol1"),"Linear log(model)","L");
390  legend.AddEntry(f1log,"#scale[0.9]{Linearized model #scale[0.7]{(unused)}}","L");
391  }
392  else {
393  legend.AddEntry(graph->GetFunction("pol1"),"Linearized model","L");
394  legend.AddEntry(pol1w,"Weighted fit #scale[0.7]{(unused)}","L");
395  }
396  legend.DrawClone();
397  }
398 
399  TLatex text;
400  text.SetNDC();
401  text.SetTextFont(42);
402  text.SetTextAlign(11);
403  text.SetTextSize(0.05);
404  //text.DrawLatex(0.20,0.20,Form("Bin %d",ibin));
405  // text.DrawLatex(0.20,0.30,"CMS inclusive jets");
406  // text.SetTextSize(0.04);
407  // text.DrawLatex(0.20,0.25,Form("%3.1f_{ }<_{ }|y|_{ }<_{ }%3.1f",input_table["ylow"][ibin],input_table["yhigh"][ibin]));
408  // text.DrawLatex(0.20,0.20,Form("%3.0f_{ }<_{ }p_{T}_{ }<_{ }%3.0f_{ }GeV",input_table["ptlow"][ibin],input_table["pthigh"][ibin]));
409 
410  text.SetTextSize(0.04);
411  text.DrawLatex(0.20,0.93,Form("%3.1f_{ }<_{ }|y|_{ }<_{ }%3.1f",input_table["ylow"][ibin],input_table["yhigh"][ibin]));
412  text.DrawLatex(0.20,0.88,Form("%3.0f_{ }<_{ }p_{T}_{ }<_{ }%3.0f_{ }GeV",input_table["ptlow"][ibin],input_table["pthigh"][ibin]));
413 
414  cout<<input_table["ylow"][ibin]<<"\t"
415  <<input_table["yhigh"][ibin]<<"\t"
416  <<input_table["ptlow"][ibin]<<"\t"
417  <<input_table["pthigh"][ibin]<<endl;
418 
419  c1.Print(ps_name);
420  if ( fit.GetLogNormal() )
421  c1.Print( Form("plots/LTFlog_bin_%02d.pdf",ibin));
422  else
423  c1.Print( Form("plots/LTF_bin_%02d.pdf",ibin));
424 
425  }
426 
427 
428  // ---------------------------------------------- //
429  // chisq plot
430  // ---------------------------------------------- //
431  TGraph* gChi2 = new TGraph();
432  int ndf = (fit.Dt.rows()-(fit.M.cols()-1));
433  for ( int itmpl = 0 ; itmpl<fit.chisq_y.size() ; itmpl++ ) {
434  gChi2->SetPoint(itmpl, reference_values[itmpl], fit.chisq_y(itmpl)/ndf);
435  }
436  gChi2->Fit("pol2","QW");
437  gChi2->SetMarkerStyle(20);
438 
439  TGraph* gChi2LTF = new TGraph();
440  gChi2LTF->SetPoint(0, fit.ahat(0), fit.chisq/ndf);
441  gChi2LTF->SetMarkerStyle(29);
442  gChi2LTF->SetMarkerSize(3.1);
443  gChi2LTF->SetMarkerColor(kViolet+2);
444 
445  TGraph* gChi2chk = new TGraph();
446  gChi2chk->SetPoint(0, fit.achk(0), fit.achk_chisq/ndf);
447  gChi2chk->SetMarkerSize(2.2);
448  gChi2chk->SetMarkerStyle(24);
449  gChi2chk->SetMarkerColor(kRed);
450 
451  //gChi2LTF->Print("all");
452 
453  //gChi2->SetTitle(";#alpha_{0} [unit];#chi^{2}/ndf");
454  gChi2->SetTitle((";"+referencename+";#chi^{2}/ndf").c_str());
455  gChi2->SetMinimum(0.75);
456  gChi2->SetMaximum(1.2);
457  gChi2->Draw("ap");
458  if ( reference_values.size()+1<= 6 )
459  gChi2->GetHistogram()->SetNdivisions(reference_values.size()+1+500,"X");
460  else
461  gChi2->GetHistogram()->SetNdivisions(int(reference_values.size()/2)+1+500,"X");
462 
463  TLine line;
464  line.SetLineColor(920);
465  line.SetLineStyle(3);
466  line.DrawLine(
467  gChi2->GetHistogram()->GetXaxis()->GetXmin(),1.,
468  gChi2->GetHistogram()->GetXaxis()->GetXmax(), 1);
469 
470  line.DrawLine(
471  fit.achk_chisq/ndf+1,1.,
472  fit.achk_chisq/ndf+1,1.);
473 
474  gChi2chk->Draw("Psame");
475  gChi2LTF->Draw("Psame");
476 
477  {
478  TLegend legend(0.18,0.70,0.66,0.97,"","NDC");
479  //legend.SetNColumns(3);
480  legend.SetFillStyle(0);
481  legend.SetBorderSize(0);
482  legend.SetTextSize(0.045);
483  legend.AddEntry(gChi2LTF,"#hat#chi^{2} of the linear template fit","P");
484  legend.AddEntry(gChi2, "#chi^{2}_{#font[12]{j}} of the individual templates","P");
485  legend.AddEntry(gChi2->GetFunction("pol2"),"Parabola","L");
486  legend.AddEntry(gChi2chk,"Minimum of #chi^{2} parabola #scale[0.8]{(#check#chi^{2})}","P"); // (#check#chi^{2})
487  legend.DrawClone();
488  }
489 
490  c1.Print(ps_name);
491  c1.Print( "plots/LTF_chi2.pdf");
492 
493  c1.Print( (string(ps_name)+"]").c_str() );
494 
495 }
496 
497 
498 // __________________________________________________________________________________ //
506 static
507 void plotLiTeFit_2D(const LTF::LiTeFit& fit, const vector<double> bins ){
508 
509  auto& M = fit.M;
510  if ( M.cols() != 3 ) {cout<<"Error! only 2-dim plotting implemented."<<endl;exit(1);}
511  Eigen::VectorXd reference_values1 = M.col(1);
512  Eigen::VectorXd reference_values2 = M.col(2);
513 
514  map<double,TH1D*> templates;
515  set<double> r1,r2;
516  for ( int iref = 0 ; iref<reference_values1.size() ; iref++ ) {
517  //double ref = reference_values1(iref);
518  templates[iref] = MakeHistogram(fit.Y.col(iref),bins);
519  r1.insert(reference_values1(iref));
520  r2.insert(reference_values2(iref));
521  }
522 
523  if ( r1.size()<=1 ) {
524  cout<<"ERROR! at least two distinct reference points for dimension 1 must be given"<<endl;
525  exit(1);
526  }
527  if ( r2.size()<=1 ) {
528  cout<<"ERROR! at least two distinct reference points for dimension 2 must be given"<<endl;
529  exit(1);
530  }
531 
532  double xmin = 2.* (*r1.begin()) - (*(++r1.begin()))*0.9999;
533  double xmax = 2.* (*r1.rbegin()) - (*(++r1.rbegin()))*1.00001;
534  double ymin = 2.* (*r2.begin()) - (*(++r2.begin()))*0.9999;
535  double ymax = 2.* (*r2.rbegin()) - (*(++r2.rbegin()))*1.00001;
536 
537  TH1D* data = MakeHistogram(fit.Dt,bins,fit.Vs);
538  TH1D* TheoFit = MakeHistogram(fit.TheoFit,bins);
539 
540  TCanvas c1("c1","LTF plots",800,800);
541  c1.SetRightMargin(0.05);
542  c1.SetLeftMargin(0.15);
543  c1.SetTopMargin(0.08);
544  // c1.SetRightMargin(0.02);
545 
546  const char* ps_name = "LTF2D_plots.ps";
547  c1.Print( (string(ps_name)+"[").c_str() );
548 
549  // ---------------------------------------------- //
550  // main plot
551  // ---------------------------------------------- //
552  for ( int iref = 0 ; iref<reference_values1.size() ; iref++ ) {
553  templates[iref]->SetLineWidth(2);
554  if ( iref == 0 ) {
555  templates[0]->SetTitle("Linear Template Fit;Observable [unit]; Value [unit]");
556  templates[0]->SetLineColor(kRed+2);
557  templates[0]->SetMinimum(0);
558  templates[0]->SetMaximum(templates[0]->GetMaximum()*1.7);
559  templates[0]->SetLineWidth(3);
560  templates[0]->DrawClone("hist");
561  }
562  else if ( iref==reference_values1.size()-1) {
563  templates[iref]->SetFillColorAlpha(kBlue,0.15);
564  templates[iref]->SetLineColor(kBlue+2);
565  templates[iref]->SetLineWidth(3);
566  templates[iref]->Draw("histsame");
567  }
568  else {
569  int color = iref+1;
570  if (color >= 10 ) color=(color-10)*2+28;
571  templates[iref]->SetLineColor(color);
572  //templates[iref]->SetFillColorAlpha(color,0.08);
573  templates[iref]->Draw("histsame");
574  }
575  }
576  templates[0]->SetFillColorAlpha(kRed,0.15);
577  templates[0]->Draw("histsame");
578 
579  data->SetMarkerStyle(20);
580  data->SetMarkerSize(1.4);
581  data->SetLineColor(kBlack);
582  data->Draw("e0same");
583 
584  TheoFit->SetLineColor(923);
585  TheoFit->SetLineWidth(4);
586  TheoFit->SetLineStyle(2);
587  TheoFit->Draw("histsame");
588 
589  TLegend legend(0.18,0.70,0.94,0.92,"","NDC");
590  legend.SetNColumns(3);
591  legend.SetFillStyle(0);
592  legend.SetBorderSize(0);
593  legend.AddEntry(data,"Data","E0P");
594  for ( int iref = 0 ; iref<reference_values1.size() ; iref++ ) {
595  legend.AddEntry(templates[iref],Form("Tpl. #alpha_{0}=%5.1f, #alpha_{1}=%3.1f",reference_values1[iref],reference_values2[iref]),"FL");
596  }
597  legend.AddEntry(TheoFit,"Estimated best model","L");
598  legend.Draw();
599 
600  c1.Print(ps_name);
601  c1.Print("plots/LTF2D_plot.pdf");
602 
603 
604  // ---------------------------------------------- //
605  // print linear-functions in every bin
606  // ---------------------------------------------- //
607  c1.Clear();
608  c1.SetRightMargin(0.02);
609  c1.SetTopMargin(0.02);
610 
611  c1.SetLeftMargin(0.12);
612  c1.SetBottomMargin(0.12);
613 
614  for ( int ibin = 0 ; ibin<fit.Dt.size() ; ibin++ ) {
615  //TGraphErrors* data = MakeTGraph(fit.ahat.row(0),fit.Dt.row(ibin));
616  TGraph2DErrors* gdata = new TGraph2DErrors();
617  gdata->SetPoint(0, fit.ahat(0), fit.ahat(1), fit.Dt(ibin));
618  gdata->SetPointError(0,0,0,data->GetBinError(ibin+1));
619  gdata->SetMarkerStyle(20);
620  gdata->SetMarkerSize(2);
621 
622 
623  TGraph2DErrors* graph = MakeTGraph2D(fit.M.col(1),fit.M.col(2),ibin,fit.Y,fit.SysY);
624  TF2* f2 = new TF2(Form("Linear Template Fit (2-dim., bin %d)",ibin),"[0]+[1]*x+[2]*y",
625  xmin,xmax,ymin,ymax
626  );
627  graph->Fit(f2,"QW");
628 
629 
630 
631  TGraph2D* gtheo = new TGraph2D();
632  gtheo->SetPoint(0, fit.ahat(0), fit.ahat(1), f2->Eval(fit.ahat(0), fit.ahat(1)));
633  gtheo->SetMarkerStyle(20);
634  gtheo->SetMarkerSize(0.4);
635  gtheo->SetMarkerColor(kBlack);
636  for ( int it = 0 ; it<fit.M.col(1).size(); it++ ) {
637  gtheo->SetPoint(it+1, reference_values1(it), reference_values2(it),
638  f2->Eval(reference_values1(it), reference_values2(it)));
639  }
640 
641  f2->SetLineColor(922);
642  f2->SetLineStyle(3);
643  f2->SetLineWidth(1);
644  f2->SetMinimum(0);
645  f2->Draw("");
646  f2->GetHistogram()->GetXaxis()->SetNdivisions(r1.size()+1);
647  f2->GetHistogram()->GetYaxis()->SetNdivisions(r2.size()+1);
648  f2->GetHistogram()->GetXaxis()->SetTitleOffset(2.2);
649  f2->GetHistogram()->GetYaxis()->SetTitleOffset(2.2);
650  f2->GetHistogram()->GetZaxis()->SetTitleOffset(1.8);
651  f2->GetHistogram()->SetTitle(";Reference value 0 (#alpha_{0}) [unit] ;Reference value 1 (#alpha_{1}) [unit];Value [unit]");
652 
653  graph->SetMarkerStyle(25);
654  //graph->SetMarkerSize(2);
655  graph->SetMarkerColor(kRed+2);
656  graph->SetLineColor(kRed+2);
657  graph->SetTitle(";Reference value (#alpha) [unit];Value [unit]");
658 
659  //graph->SetMinimum(0);
660  //graph->SetMaximum( max(gdata->GetY()[0],max(graph->GetY()[0],graph->GetY()[graph->GetN()-1]))*1.2);
661  //graph->Draw("ap");
662  //pol1w->Draw("Lsame");
663  //graph->GetFunction("pol1")->Draw("Lsame");
664  //graph->GetHistogram()->GetXaxis()->SetNdivisions(graph->GetN()+1);
665 
666  f2->Draw("surf1");
667  gtheo->Draw("P0 same");
668  graph->Draw("err p0 same");
669  gdata->Draw("err P0 same");
670 
671 
672  TF2* f2leg = (TF2*)f2->Clone("f2leg");
673  f2leg->SetLineWidth(3);
674  gtheo->SetMarkerStyle(24);
675  gtheo->SetMarkerSize(0.4);
676  gdata->SetMarkerStyle(24);
677  graph->SetMarkerStyle(24);
678  if ( ibin==1 ) {
679  TLegend legend(0.42,0.18,0.80,0.43,"","NDC");
680  //legend.SetNColumns(3);
681  legend.SetFillStyle(0);
682  legend.SetBorderSize(0);
683  legend.SetTextSize(0.03);
684  legend.AddEntry(graph,"Templates","PE0");
685  legend.AddEntry(f2leg,"Linearized model","L");
686  legend.AddEntry(gtheo,"Projections onto model","P");
687  legend.AddEntry(gdata,"Data","E0P");
688  legend.DrawClone();
689  }
690 
691  TLatex text;
692  text.SetNDC();
693  text.SetTextFont(42);
694  text.SetTextAlign(13);
695  text.SetTextSize(0.05);
696  //text.DrawLatex(0.75,0.25,Form("Bin %d",ibin));
697  text.DrawLatex(0.02,0.97,Form("Bin %d",ibin));
698 
699  c1.Print(ps_name);
700  c1.Print( Form("plots/LTF2D_bin_%02d.pdf",ibin));
701 
702  }
703 
704  c1.Print( (string(ps_name)+"]").c_str() );
705 
706 }
707 
708 
709 
710 
711 
712 
713 // __________________________________________________________________________________ //
721 static
722 void plotLiTeFitPol2Test(const LTF::LiTeFit& fit, const vector<double>& bins,
723  const string& yaxistitle = "value [unit]",
724  const string& referencename = "Reference value (#alpha) [unit]",
725  const string& observablename = "Observable [unit]"
726  ) {
727  gStyle->SetOptStat(0);
728  gSystem->mkdir("plots");
729  auto& M = fit.M;
730 
731  // sanity check
732  if ( M.cols() != 2 ) {cout<<"Error! only 1-dim plotting implemented."<<endl;exit(1);}
733  Eigen::VectorXd reference_values = M.col(1);
734  //TH1D* hist = new TH1D("hist","hist",bins.size()-1, &bins[0]);
735 
736  map<double,TH1D*> templates;
737  for ( int iref = 0 ; iref<reference_values.size() ; iref++ ) {
738  //double ref = reference_values(iref);
739  templates[iref] = MakeHistogram(fit.Y.col(iref),bins);
740  }
741 
742  TH1D* data = MakeHistogram(fit.Dt,bins,fit.Vs);
743  TH1D* TheoFit = MakeHistogram(fit.TheoFit,bins);
744 
745  TCanvas c1("c1","LTF plots",800,800);
746  c1.SetRightMargin(0.05);
747  c1.SetLeftMargin(0.15);
748  c1.SetTopMargin(0.08);
749  // c1.SetRightMargin(0.02);
750 
751  const char* ps_name = fit.GetLogNormal() ?
752  "LTFlog_plots.ps" :
753  "LTF_plots.ps";
754  c1.Print( (string(ps_name)+"[").c_str() );
755 
756  // ---------------------------------------------- //
757  // main plot
758  // ---------------------------------------------- //
759  for ( int iref = 0 ; iref<reference_values.size() ; iref++ ) {
760  templates[iref]->SetLineWidth(2);
761  if ( iref == 0 ) {
762  templates[0]->SetTitle(("Linear Template Fit;"+observablename+";"+yaxistitle).c_str());
763  templates[0]->SetLineColor(kRed+1);
764  if ( templates[0]->GetMaximum()>0 )templates[0]->SetMinimum(0);
765  if ( !fit.GetLogNormal() )
766  templates[0]->SetMaximum(templates[0]->GetMaximum()*1.7);
767  // else
768  // templates[0]->SetMaximum(templates[0]->GetMaximum()*3.);
769  templates[0]->SetLineWidth(3);
770  templates[0]->DrawClone("hist");
771  }
772  else if ( iref==reference_values.size()-1) {
773  if ( templates[0]->GetMaximum()>0 ) templates[iref]->SetFillColorAlpha(kBlue,0.15);
774  templates[iref]->SetLineColor(kBlue+2);
775  templates[iref]->SetLineWidth(3);
776  templates[iref]->Draw("histsame");
777  }
778  else {
779  templates[iref]->SetLineColor(iref+2);
780  templates[iref]->SetLineWidth(2);
781  templates[iref]->Draw("histsame");
782  }
783  }
784  if ( templates[0]->GetMaximum()>0 )templates[0]->SetFillColorAlpha(kRed,0.15);
785  templates[0]->Draw("histsame");
786 
787  data->SetMarkerStyle(20);
788  data->SetMarkerSize(1.4);
789  data->SetLineColor(kBlack);
790  data->Draw("e0same");
791 
792  TheoFit->SetLineColor(923);
793  TheoFit->SetLineWidth(4);
794  TheoFit->SetLineStyle(2);
795  TheoFit->Draw("histsame");
796 
797  TLegend legend(0.18,0.70,0.94,0.92,"","NDC");
798  legend.SetNColumns(3);
799  legend.SetFillStyle(0);
800  legend.SetBorderSize(0);
801  legend.AddEntry(data,"Data","E0P");
802  for ( int iref = 0 ; iref<reference_values.size() ; iref++ ) {
803  legend.AddEntry(templates[iref],Form("Template #alpha=%6.2f",reference_values[iref]),"FL");
804  }
805  legend.AddEntry(TheoFit,"Estimated best model","L");
806  legend.Draw();
807 
808  c1.Print(ps_name);
809  c1.Print("plots/LTF_plot.pdf");
810 
811 
812  // ---------------------------------------------- //
813  // print linear-functions in every bin
814  // ---------------------------------------------- //
815  c1.Clear();
816  c1.SetRightMargin(0.02);
817  c1.SetTopMargin(0.02);
818  c1.SetLeftMargin(0.16);
819  c1.SetBottomMargin(0.16);
820 
821  gPad->SetTicky(1);
822 
823  gStyle->SetLabelSize(0.05,"XYZ");
824  gStyle->SetTitleSize(0.05,"XYZ");
825  gStyle->SetTitleOffset(1.1,"X");
826  gStyle->SetTitleOffset(1.6,"Y");
827  gStyle->SetMarkerSize(2);
828 
829 
830  map < string, vector<double> > input_table = read_input_table2("data/CMS_data.txt",32);
831 
832  for ( int ibin = 0 ; ibin<fit.Dt.size() ; ibin++ ) {
833  //for ( int ibin = 0 ; ibin<0 ; ibin++ ) {
834  //for ( int ibin = 0 ; ibin<6 ; ibin++ ) {
835  //TGraphErrors* data = MakeTGraph(fit.ahat.row(0),fit.Dt.row(ibin));
836  TGraphErrors* gdata = new TGraphErrors();
837  gdata->SetPoint(0,fit.ahat(0),fit.Dt(ibin));
838  gdata->SetPointError(0,0,data->GetBinError(ibin+1));
839  gdata->SetMarkerStyle(20);
840 
841  TGraphErrors* graph = MakeTGraph(fit.M.col(1),ibin,fit.Y,fit.SysY);
842  graph->Fit("pol1","QW"); // in this function it is "QW" (unweighted
843  TF1* pol1 = (TF1*)graph->GetFunction("pol1")->Clone("pol1");
844  pol1->SetLineColor(kRed);
845  //pol1w->SetLineStyle(3);
846  pol1->SetLineWidth(2);
847 
848  TF1* f2= new TF1("f2","[0]+[1]*pow(x,2)", -FLT_MIN,FLT_MAX );
849  graph->Fit(f2,"QW");
850  TF1* f1= new TF1("f1","[0]+[1]*pow(x,[2])", -FLT_MIN,FLT_MAX );
851 
852 
853  graph->Fit("pol1","QW"); // "W": Ignore all point errors when fitting a TGraphErrors
854  bool UseLTFOutput = true;
855  if ( UseLTFOutput ) {
856  Eigen::VectorXd ltfpol1param =(fit.Y*fit.Mc().transpose()).row(ibin);
857  //cout<<"M+*Y:"<<endl<< ltfpol1param <<endl;
858  graph->GetFunction("pol1")->SetParameter(0,ltfpol1param(0));
859  graph->GetFunction("pol1")->SetParameter(1,ltfpol1param(1));
860  f1->SetParameter(0,ltfpol1param(0));
861  f1->SetParameter(1,ltfpol1param(1));
862  f1->SetParameter(2,fit.Gamma[0]);
863  }
864  if ( fit.Gamma[0]!=1 ) cout<<"Warning! Plotting with gamma factor !=1 not correctly implmeneted!"<<endl;
865 
866  graph->GetFunction("pol1")->SetLineWidth(3);
867 
868  graph->Fit("pol2","QW");
869  graph->GetFunction("pol2")->SetLineWidth(83);
870 
871  graph->SetMarkerStyle(47);
872  graph->SetMarkerColor(kRed+2);
873  graph->SetLineColor(kRed+2);
874  if ( fit.GetLogNormal() )
875  graph->SetTitle((";"+referencename+";log("+yaxistitle+")").c_str()); // log(value/unit)
876  else
877  graph->SetTitle((";"+referencename+";"+yaxistitle).c_str());
878 
879 
880  TF1* f1log = NULL;
881  if ( fit.GetLogNormal() ) {
882  //f1log = new TF1("pol1log","[0]+[1]*exp(x)", -FLT_MIN,FLT_MAX );
883  f1log = new TF1("pol1log","log([0]+[1]*pow(x,1))", -FLT_MIN,FLT_MAX );
884 
885  TGraph* gexp = new TGraph();
886  for ( int i = 0 ; i<graph->GetN() ; i++ )
887  gexp->SetPoint(i,graph->GetX()[i], exp(graph->GetY()[i]));
888  gexp->Fit("pol1","QW");
889 
890  Eigen::VectorXd ltfpol1param =(fit.Y*fit.Mc().transpose()).row(ibin);
891  graph->GetFunction("pol1")->SetParameter(0,ltfpol1param(0));
892 
893  f1log->SetParameter(0, gexp->GetFunction("pol1")->GetParameter(0));
894  f1log->SetParameter(1, gexp->GetFunction("pol1")->GetParameter(1));
895  if ( fit.Gamma[0]!=1 ) cout<<"Warning! Plotting with gamma factor !=1 not correctly implmeneted!"<<endl;
896 
897  // //f1 = new TF1("pol1","[0]+[1]*x", -FLT_MIN,FLT_MAX );
898  // graph->Fit(f1log,"QW"); // "W": Ignore all point errors when fitting a TGraphErrors
899  f1log->SetLineColor(kBlue+1);
900  f1log->SetLineStyle(7);
901  f1log->SetLineWidth(2);
902  }
903 
904 
905  if ( gdata->GetY()[0] > 0 )
906  graph->SetMinimum(0);
907 
908  graph->SetMaximum( max(gdata->GetY()[0],max(graph->GetY()[0],graph->GetY()[graph->GetN()-1]))*1.2);
909  graph->Fit("pol2","QW");
910  graph->GetFunction("pol2")->SetLineColor(kBlue);;
911 
912  TF1* tang = new TF1("tang","[0]+[1]*x", -FLT_MIN,FLT_MAX );
913  double f0 = graph->GetFunction("pol2")->Eval(gdata->GetX()[0]);
914  double fp0 = graph->GetFunction("pol2")->GetParameter(1) + 2.*graph->GetFunction("pol2")->GetParameter(2)*gdata->GetX()[0];
915  tang->SetParameter(0, f0 - fp0*gdata->GetX()[0] );
916  tang->SetParameter(1, fp0);
917  tang->SetLineWidth(2);
918  tang->SetLineColor(kTeal+4);
919  tang->SetLineStyle(7);
920 
921  pol1->SetLineStyle(3);
922  graph->Draw("ap");
923  if ( fit.GetLogNormal() )
924  f1log->Draw("Lsame");
925  else
926  pol1->Draw("Lsame");
927  //graph->GetFunction("pol1")->Draw("Lsame");
928  if ( reference_values.size()+1<= 8 )
929  graph->GetHistogram()->GetXaxis()->SetNdivisions(graph->GetN()+1);
930  else
931  graph->GetHistogram()->GetXaxis()->SetNdivisions(int(graph->GetN()/2)+1+200);
932  //f2->Draw("same");
933 
934  tang->Draw("Lsame");
935  gdata->Draw("P");
936 
937  if ( ibin==0 ) {
938  double xmin = fit.GetLogNormal() ? 0.36 : 0.45;
939  //TLegend legend(xmin,0.19,0.96,0.47,"","NDC");
940  TLegend legend(0.19,0.19,0.56,0.47,"","NDC");
941  //legend.SetNColumns(3);
942  legend.SetFillStyle(0);
943  legend.SetBorderSize(0);
944  legend.SetTextSize(0.045);
945  legend.AddEntry(data,"Data","E0P");
946  legend.AddEntry(graph,"Templates","PE0");
947  if ( fit.GetLogNormal() ) {
948  legend.AddEntry(graph->GetFunction("pol1"),"Linear log(model)","L");
949  legend.AddEntry(f1log,"#scale[0.9]{Linearized model #scale[0.7]{(unused)}}","L");
950  }
951  else {
952  legend.AddEntry(graph->GetFunction("pol2"),"Second-order model","L");
953  legend.AddEntry(tang,"Linear model","L");
954  legend.AddEntry(pol1,"Linear Template Fit","L");
955  }
956  legend.DrawClone();
957  }
958 
959  TLatex text;
960  text.SetNDC();
961  text.SetTextFont(42);
962  text.SetTextAlign(11);
963  text.SetTextSize(0.045);
964  //text.DrawLatex(0.20,0.20,Form("Bin %d",ibin));
965 
966  text.SetTextAlign(31);
967  text.DrawLatex(0.955,0.20,Form("Bin %d",ibin));
968 
969  // text.DrawLatex(0.20,0.30,"CMS inclusive jets");
970  // text.SetTextSize(0.04);
971  // text.DrawLatex(0.20,0.25,Form("%3.1f_{ }<_{ }|y|_{ }<_{ }%3.1f",input_table["ylow"][ibin],input_table["yhigh"][ibin]));
972  // text.DrawLatex(0.20,0.20,Form("%3.0f_{ }<_{ }p_{T}_{ }<_{ }%3.0f_{ }GeV",input_table["ptlow"][ibin],input_table["pthigh"][ibin]));
973 
974  // text.SetTextSize(0.04);
975  // text.DrawLatex(0.20,0.93,Form("%3.1f_{ }<_{ }|y|_{ }<_{ }%3.1f",input_table["ylow"][ibin],input_table["yhigh"][ibin]));
976  // text.DrawLatex(0.20,0.88,Form("%3.0f_{ }<_{ }p_{T}_{ }<_{ }%3.0f_{ }GeV",input_table["ptlow"][ibin],input_table["pthigh"][ibin]));
977 
978  cout<<input_table["ylow"][ibin]<<"\t"
979  <<input_table["yhigh"][ibin]<<"\t"
980  <<input_table["ptlow"][ibin]<<"\t"
981  <<input_table["pthigh"][ibin]<<endl;
982 
983  c1.Print(ps_name);
984  if ( fit.GetLogNormal() )
985  c1.Print( Form("plots/LTFlog_bin_%02d.pdf",ibin));
986  else
987  c1.Print( Form("plots/LTF_bin_%02d.pdf",ibin));
988 
989  }
990 
991 
992  // ---------------------------------------------- //
993  // chisq plot
994  // ---------------------------------------------- //
995  TGraph* gChi2 = new TGraph();
996  int ndf = (fit.Dt.rows()-(fit.M.cols()-1));
997  for ( int itmpl = 0 ; itmpl<fit.chisq_y.size() ; itmpl++ ) {
998  gChi2->SetPoint(itmpl, reference_values[itmpl], fit.chisq_y(itmpl)/ndf);
999  }
1000  gChi2->Fit("pol2","QW");
1001  gChi2->SetMarkerStyle(20);
1002 
1003  TGraph* gChi2LTF = new TGraph();
1004  gChi2LTF->SetPoint(0, fit.ahat(0), fit.chisq/ndf);
1005  gChi2LTF->SetMarkerStyle(29);
1006  gChi2LTF->SetMarkerSize(3.1);
1007  gChi2LTF->SetMarkerColor(kViolet+2);
1008 
1009  TGraph* gChi2chk = new TGraph();
1010  gChi2chk->SetPoint(0, fit.achk(0), fit.achk_chisq/ndf);
1011  gChi2chk->SetMarkerSize(2.2);
1012  gChi2chk->SetMarkerStyle(24);
1013  gChi2chk->SetMarkerColor(kRed);
1014 
1015  //gChi2LTF->Print("all");
1016 
1017  //gChi2->SetTitle(";#alpha_{0} [unit];#chi^{2}/ndf");
1018  gChi2->SetTitle((";"+referencename+";#chi^{2}/ndf").c_str());
1019  gChi2->SetMinimum(0.);
1020  gChi2->SetMaximum(20.2);
1021  gChi2->Draw("apc");
1022  if ( reference_values.size()+1<= 8 )
1023  gChi2->GetHistogram()->SetNdivisions(reference_values.size()+1+500,"X");
1024  else
1025  gChi2->GetHistogram()->SetNdivisions(int(reference_values.size()/2)+1+500,"X");
1026 
1027  TLine line;
1028  line.SetLineColor(920);
1029  line.SetLineStyle(3);
1030  line.DrawLine(
1031  gChi2->GetHistogram()->GetXaxis()->GetXmin(),1.,
1032  gChi2->GetHistogram()->GetXaxis()->GetXmax(), 1);
1033 
1034  line.DrawLine(
1035  fit.achk_chisq/ndf+1,1.,
1036  fit.achk_chisq/ndf+1,1.);
1037 
1038  gChi2chk->Draw("Psame");
1039  gChi2LTF->Draw("Psame");
1040 
1041  {
1042  TLegend legend(0.18,0.75,0.56,0.97,"","NDC");
1043  //legend.SetNColumns(3);
1044  legend.SetFillStyle(0);
1045  legend.SetBorderSize(0);
1046  legend.SetTextSize(0.045);
1047  legend.AddEntry(gChi2LTF,"#hat#chi^{2} of the iterative template fit","P");
1048  legend.AddEntry(gChi2, "#chi^{2}_{#font[12]{j}} of the individual templates","PL");
1049  legend.AddEntry(gChi2->GetFunction("pol2"),"Parabola","L");
1050  //legend.AddEntry(gChi2chk,"Minimum of #chi^{2} parabola #scale[0.8]{(#check#chi^{2})}","P"); // (#check#chi^{2})
1051  legend.DrawClone();
1052  }
1053 
1054  c1.Print(ps_name);
1055  c1.Print( "plots/LTF_chi2.pdf");
1056 
1057  c1.Print( (string(ps_name)+"]").c_str() );
1058 
1059 }
1060 
1061 
1062 
1063 };
static void plotLiTeFitPol2Test(const LTF::LiTeFit &fit, const vector< double > &bins, const string &yaxistitle="value [unit]", const string &referencename="Reference value (#alpha) [unit]", const string &observablename="Observable [unit]")
Definition: LTF_ROOTTools.h:722
Eigen::MatrixXd Mc(int nPow=1, int nInfrc=0) const
calculated matrix M^+ from matrix M
Definition: LTF.cxx:389
Eigen::MatrixXd Y
Matrix Y (templates)
Definition: LTF.h:83
Eigen::VectorXd Dt
data distribution
Definition: LTF.h:81
double achk_chisq
chisq_min from polynomial fit to chi2's of templates
Definition: LTF.h:123
Definition: LTF_ROOTTools.h:28
Eigen::MatrixXd M
Matrix M (reference values)
Definition: LTF.h:82
Eigen::VectorXd achk
result from polynomial fit to chi2's of templates
Definition: LTF.h:121
bool GetLogNormal() const
return boolean if fit was performed with relative uncertainties
Definition: LTF.h:140
Eigen::VectorXd ahat
results (a0, a1, ..., epsilon1, epsilon2,...)
Definition: LTF.h:107
static TGraph2DErrors * MakeTGraph2D(const Eigen::VectorXd &xvalues, const Eigen::VectorXd &yvalues, int ibin, const Eigen::MatrixXd &Y, const std::map< std::string, Eigen::MatrixXd > &VSysY={})
Definition: LTF_ROOTTools.h:100
static TH1D * MakeHistogram(const Eigen::VectorXd &values, vector< double > bins={}, const std::vector< std::pair< std::string, Eigen::MatrixXd > > &V={})
Definition: LTF_ROOTTools.h:148
Small helper class to collect all results, All input parameters are set by class LTF The user can acc...
Definition: LTF.h:65
static TH1D * MakeHistogram(int nEvents, int seed, double mean, double sigma, vector< double > bins)
Definition: LTF_ROOTTools.h:127
std::map< std::string, Eigen::MatrixXd > SysY
uncertainty of the templates
Definition: LTF.h:91
std::vector< std::pair< std::string, Eigen::MatrixXd > > Vs
(uncor./Cov) uncertainties
Definition: LTF.h:86
static TGraphErrors * MakeTGraph(const Eigen::VectorXd &xvalues, int ibin, const Eigen::MatrixXd &Y, const std::map< std::string, Eigen::MatrixXd > &VSysY={})
Definition: LTF_ROOTTools.h:75
static void plotLiTeFit(const LTF::LiTeFit &fit, const vector< double > &bins, const string &yaxistitle="value [unit]", const string &referencename="Reference value (#alpha) [unit]", const string &observablename="Observable [unit]")
Definition: LTF_ROOTTools.h:177
static std::map< std::string, std::vector< double > > read_input_table2(std::string filename, int ncol)
Definition: LTF_ROOTTools.h:38
Eigen::VectorXd chisq_y
chi^2 for each template
Definition: LTF.h:105
Eigen::VectorXd TheoFit
best-fit theory predictions
Definition: LTF.h:112
vector< double > Gamma
Use non-linear Gamma-factor.
Definition: LTF.h:94
double chisq
chi^2
Definition: LTF.h:102
static void plotLiTeFit_2D(const LTF::LiTeFit &fit, const vector< double > bins)
Definition: LTF_ROOTTools.h:507