Logo ROOT   6.08/07
Reference Guide
testUnfold3.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_unfold
3 /// \notebook
4 /// Simple Test program for the class TUnfoldDensity
5 ///
6 /// 1-dimensional unfolding with background subtraction
7 ///
8 /// the goal is to unfold the underlying "true" distribution of a variable Pt
9 ///
10 /// the reconstructed Pt is measured in 24 bins from 4 to 28
11 /// the generator-level Pt is unfolded into 10 bins from 6 to 26
12 /// plus underflow bin from 0 to 6
13 /// plus overflow bin above 26
14 /// there are two background sources
15 /// bgr1 and bgr2
16 /// the signal has a finite trigger efficiency at a threshold of 8 GeV
17 ///
18 /// one type of systematic error is studied, where the signal parameters are
19 /// changed
20 ///
21 /// Finally, the unfolding is compared to a "bin-by-bin" correction method
22 ///
23 /// History:
24 /// - Version 17.0, changeto use the class TUnfoldDensity
25 /// - Version 16.1, parallel to changes in TUnfold
26 /// - Version 16.0, parallel to changes in TUnfold
27 /// - Version 15, simple example including background subtraction
28 ///
29 /// \macro_image
30 /// \macro_output
31 /// \macro_code
32 ///
33 /// \author Stefan Schmitt, DESY
34 
35 #include <TMath.h>
36 #include <TCanvas.h>
37 #include <TRandom3.h>
38 #include <TFitter.h>
39 #include <TF1.h>
40 #include <TStyle.h>
41 #include <TVector.h>
42 #include <TGraph.h>
43 
44 #include "TUnfoldDensity.h"
45 
46 using namespace std;
47 
48 
49 TRandom *rnd=0;
50 
51 Double_t GenerateEvent(const Double_t *parm,
52  const Double_t *triggerParm,
53  Double_t *intLumi,
54  Bool_t *triggerFlag,
55  Double_t *ptGen,Int_t *iType)
56 {
57  // generate an event
58  // input:
59  // parameters for the event generator
60  // return value:
61  // reconstructed Pt
62  // output to pointers:
63  // integrated luminosity
64  // several variables only accessible on generator level
65  //
66  // the parm array defines the physical parameters
67  // parm[0]: background source 1 fraction
68  // parm[1]: background source 2 fraction
69  // parm[2]: lower limit of generated Pt distribution
70  // parm[3]: upper limit of generated Pt distribution
71  // parm[4]: exponent for Pt distribution signal
72  // parm[5]: exponent for Pt distribution background 1
73  // parm[6]: exponent for Pt distribution background 2
74  // parm[7]: resolution parameter a goes with sqrt(Pt)
75  // parm[8]: resolution parameter b goes with Pt
76  // triggerParm[0]: trigger threshold turn-on position
77  // triggerParm[1]: trigger threshold turn-on width
78  // triggerParm[2]: trigger efficiency for high pt
79  //
80  // intLumi is advanced bu 1 for each *generated* event
81  // for data, several events may be generated, until one passes the trigger
82  //
83  // some generator-level quantities are also returned:
84  // triggerFlag: whether the event passed the trigger threshold
85  // ptGen: the generated pt
86  // iType: which type of process was simulated
87  //
88  // the "triggerFlag" also has another meaning:
89  // if(triggerFlag==0) only triggered events are returned
90  //
91  // Usage to generate data events
92  // ptObs=GenerateEvent(parm,triggerParm,0,0,0)
93  //
94  // Usage to generate MC events
95  // ptGen=GenerateEvent(parm,triggerParm,&triggerFlag,&ptGen,&iType);
96  //
97  Double_t ptObs;
98  Bool_t isTriggered=kFALSE;
99  do {
100  Int_t itype;
101  Double_t ptgen;
102  Double_t f=rnd->Rndm();
103  // decide whether this is background or signal
104  itype=0;
105  if(f<parm[0]) itype=1;
106  else if(f<parm[0]+parm[1]) itype=2;
107  // generate Pt according to distribution pow(ptgen,a)
108  // get exponent
109  Double_t a=parm[4+itype];
110  Double_t a1=a+1.0;
111  Double_t t=rnd->Rndm();
112  if(a1 == 0.0) {
113  Double_t x0=TMath::Log(parm[2]);
114  ptgen=TMath::Exp(t*(TMath::Log(parm[3])-x0)+x0);
115  } else {
116  Double_t x0=pow(parm[2],a1);
117  ptgen=pow(t*(pow(parm[3],a1)-x0)+x0,1./a1);
118  }
119  if(iType) *iType=itype;
120  if(ptGen) *ptGen=ptgen;
121 
122  // smearing in Pt with large asymmetric tail
123  Double_t sigma=
124  TMath::Sqrt(parm[7]*parm[7]*ptgen+parm[8]*parm[8]*ptgen*ptgen);
125  ptObs=rnd->BreitWigner(ptgen,sigma);
126 
127  // decide whether event was triggered
128  Double_t triggerProb =
129  triggerParm[2]/(1.+TMath::Exp((triggerParm[0]-ptObs)/triggerParm[1]));
130  isTriggered= rnd->Rndm()<triggerProb;
131  (*intLumi) ++;
132  } while((!triggerFlag) && (!isTriggered));
133  // return trigger decision
134  if(triggerFlag) *triggerFlag=isTriggered;
135  return ptObs;
136 }
137 
138 //int main(int argc, char *argv[])
139 void testUnfold3()
140 {
141  // switch on histogram errors
143 
144  // show fit result
145  gStyle->SetOptFit(1111);
146 
147  // random generator
148  rnd=new TRandom3();
149 
150  // data and MC luminosities
151  Double_t const lumiData= 30000;
152  Double_t const lumiMC =1000000;
153 
154  //========================
155  // Step 1: define binning, book histograms
156 
157  // reconstructed pt (fine binning)
158  Int_t const nDet=24;
159  Double_t const xminDet=4.0;
160  Double_t const xmaxDet=28.0;
161 
162  // generated pt (coarse binning)
163  Int_t const nGen=10;
164  Double_t const xminGen= 6.0;
165  Double_t const xmaxGen=26.0;
166 
167  //==================================================================
168  // book histograms
169  // (1) unfolding input: binning scheme is fine for detector,coarse for gen
170  // histUnfoldInput : reconstructed data, binning for unfolding
171  // histUnfoldMatrix : generated vs reconstructed distribution
172  // histUnfoldBgr1 : background source1, as predicted by MC
173  // histUnfoldBgr2 : background source2, as predicted by MC
174  // for systematic studies
175  // histUnfoldMatrixSys : generated vs reconstructed with different signal
176  //
177  // (2) histograms required for bin-by-bin method
178  // histDetDATAbbb : reconstructed data for bin-by-bin
179  // histDetMCbbb : reconstructed MC for bin-by-bin
180  // histDetMCBgr1bbb : reconstructed MC bgr1 for bin-by-bin
181  // histDetMCBgr2bbb : reconstructed MC bgr2 for bin-by-bin
182  // histDetMCBgrPtbbb : reconstructed MC bgr from low/high pt for bin-by-bin
183  // histGenMCbbb : generated MC truth
184  // for systematic studies
185  // histDetMCSysbbb : reconstructed with changed trigger
186  // histGenMCSysbbb : generated MC truth
187  //
188  // (3) monitoring and control
189  // histGenData : data truth for bias tests
190  // histDetMC : MC prediction
191 
192  // (1) create histograms required for unfolding
193  TH1D *histUnfoldInput=
194  new TH1D("unfolding input rec",";ptrec",nDet,xminDet,xmaxDet);
195  TH2D *histUnfoldMatrix=
196  new TH2D("unfolding matrix",";ptgen;ptrec",
197  nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
198  TH1D *histUnfoldBgr1=
199  new TH1D("unfolding bgr1 rec",";ptrec",nDet,xminDet,xmaxDet);
200  TH1D *histUnfoldBgr2=
201  new TH1D("unfolding bgr2 rec",";ptrec",nDet,xminDet,xmaxDet);
202  TH2D *histUnfoldMatrixSys=
203  new TH2D("unfolding matrix sys",";ptgen;ptrec",
204  nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
205 
206  // (2) histograms required for bin-by-bin
207  TH1D *histBbbInput=
208  new TH1D("bbb input rec",";ptrec",nGen,xminGen,xmaxGen);
209  TH1D *histBbbSignalRec=
210  new TH1D("bbb signal rec",";ptrec",nGen,xminGen,xmaxGen);
211  TH1D *histBbbBgr1=
212  new TH1D("bbb bgr1 rec",";ptrec",nGen,xminGen,xmaxGen);
213  TH1D *histBbbBgr2=
214  new TH1D("bbb bgr2 rec",";ptrec",nGen,xminGen,xmaxGen);
215  TH1D *histBbbBgrPt=
216  new TH1D("bbb bgrPt rec",";ptrec",nGen,xminGen,xmaxGen);
217  TH1D *histBbbSignalGen=
218  new TH1D("bbb signal gen",";ptgen",nGen,xminGen,xmaxGen);
219  TH1D *histBbbSignalRecSys=
220  new TH1D("bbb signal sys rec",";ptrec",nGen,xminGen,xmaxGen);
221  TH1D *histBbbBgrPtSys=
222  new TH1D("bbb bgrPt sys rec",";ptrec",nGen,xminGen,xmaxGen);
223  TH1D *histBbbSignalGenSys=
224  new TH1D("bbb signal sys gen",";ptgen",nGen,xminGen,xmaxGen);
225 
226  // (3) control histograms
227  TH1D *histDataTruth=
228  new TH1D("DATA truth gen",";ptgen",nGen,xminGen,xmaxGen);
229  TH1D *histDetMC=
230  new TH1D("MC prediction rec",";ptrec",nDet,xminDet,xmaxDet);
231  // ==============================================================
232  // Step 2: generate data distributions
233  //
234  // data parameters: in real life these are unknown
235  static Double_t parm_DATA[]={
236  0.05, // fraction of background 1 (on generator level)
237  0.05, // fraction of background 2 (on generator level)
238  3.5, // lower Pt cut (generator level)
239  100.,// upper Pt cut (generator level)
240  -3.0,// signal exponent
241  0.1, // background 1 exponent
242  -0.5, // background 2 exponent
243  0.2, // energy resolution a term
244  0.01, // energy resolution b term
245  };
246  static Double_t triggerParm_DATA[]={8.0, // threshold is 8 GeV
247  4.0, // width is 4 GeV
248  0.95 // high Pt efficiency os 95%
249  };
250 
251  Double_t intLumi=0.0;
252  while(intLumi<lumiData) {
253  Int_t iTypeGen;
254  Bool_t isTriggered;
255  Double_t ptGen;
256  Double_t ptObs=GenerateEvent(parm_DATA,triggerParm_DATA,&intLumi,
257  &isTriggered,&ptGen,&iTypeGen);
258  if(isTriggered) {
259  // (1) histogram for unfolding
260  histUnfoldInput->Fill(ptObs);
261 
262  // (2) histogram for bin-by-bin
263  histBbbInput->Fill(ptObs);
264  }
265  // (3) monitoring
266  if(iTypeGen==0) histDataTruth->Fill(ptGen);
267  }
268 
269  // ==============================================================
270  // Step 3: generate default MC distributions
271  //
272  // MC parameters
273  // default settings
274  static Double_t parm_MC[]={
275  0.05, // fraction of background 1 (on generator level)
276  0.05, // fraction of background 2 (on generator level)
277  3.5, // lower Pt cut (generator level)
278  100.,// upper Pt cut (generator level)
279  -4.0,// signal exponent !!! steeper than in data
280  // to illustrate bin-by-bin bias
281  0.1, // background 1 exponent
282  -0.5, // background 2 exponent
283  0.2, // energy resolution a term
284  0.01, // energy resolution b term
285  };
286  static Double_t triggerParm_MC[]={8.0, // threshold is 8 GeV
287  4.0, // width is 4 GeV
288  0.95 // high Pt efficiency is 95%
289  };
290 
291  // weighting factor to accomodate for the different luminosity in data and MC
292  Double_t lumiWeight = lumiData/lumiMC;
293  intLumi=0.0;
294  while(intLumi<lumiMC) {
295  Int_t iTypeGen;
296  Bool_t isTriggered;
297  Double_t ptGen;
298  Double_t ptObs=GenerateEvent(parm_MC,triggerParm_MC,&intLumi,&isTriggered,
299  &ptGen,&iTypeGen);
300  if(!isTriggered) ptObs=0.0;
301 
302  // (1) distribution required for unfolding
303 
304  if(iTypeGen==0) {
305  histUnfoldMatrix->Fill(ptGen,ptObs,lumiWeight);
306  } else if(iTypeGen==1) {
307  histUnfoldBgr1->Fill(ptObs,lumiWeight);
308  } else if(iTypeGen==2) {
309  histUnfoldBgr2->Fill(ptObs,lumiWeight);
310  }
311 
312  // (2) distributions for bbb
313  if(iTypeGen==0) {
314  if((ptGen>=xminGen)&&(ptGen<xmaxGen)) {
315  histBbbSignalRec->Fill(ptObs,lumiWeight);
316  } else {
317  histBbbBgrPt->Fill(ptObs,lumiWeight);
318  }
319  histBbbSignalGen->Fill(ptGen,lumiWeight);
320  } else if(iTypeGen==1) {
321  histBbbBgr1->Fill(ptObs,lumiWeight);
322  } else if(iTypeGen==2) {
323  histBbbBgr2->Fill(ptObs,lumiWeight);
324  }
325 
326  // (3) control distribution
327  histDetMC ->Fill(ptObs,lumiWeight);
328  }
329 
330  // ==============================================================
331  // Step 4: generate MC distributions for systematic study
332  // example: study dependence on initial signal shape
333  // -> BGR distributions do not change
334  static Double_t parm_MC_SYS[]={
335  0.05, // fraction of background: unchanged
336  0.05, // fraction of background: unchanged
337  3.5, // lower Pt cut (generator level)
338  100.,// upper Pt cut (generator level)
339  -2.0, // signal exponent changed
340  0.1, // background 1 exponent
341  -0.5, // background 2 exponent
342  0.2, // energy resolution a term
343  0.01, // energy resolution b term
344  };
345 
346  intLumi=0.0;
347  while(intLumi<lumiMC) {
348  Int_t iTypeGen;
349  Bool_t isTriggered;
350  Double_t ptGen;
351  Double_t ptObs=GenerateEvent(parm_MC_SYS,triggerParm_MC,&intLumi,
352  &isTriggered,&ptGen,&iTypeGen);
353  if(!isTriggered) ptObs=0.0;
354 
355  // (1) for unfolding
356  if(iTypeGen==0) {
357  histUnfoldMatrixSys->Fill(ptGen,ptObs);
358  }
359 
360  // (2) for bin-by-bin
361  if(iTypeGen==0) {
362  if((ptGen>=xminGen)&&(ptGen<xmaxGen)) {
363  histBbbSignalRecSys->Fill(ptObs);
364  } else {
365  histBbbBgrPtSys->Fill(ptObs);
366  }
367  histBbbSignalGenSys->Fill(ptGen);
368  }
369  }
370 
371  // this method is new in version 16 of TUnfold
372  cout<<"TUnfold version is "<<TUnfold::GetTUnfoldVersion()<<"\n";
373 
374  //========================
375  // Step 5: unfolding
376  //
377 
378  // here one has to decide about the regularisation scheme
379 
380  // regularize curvature
381  TUnfold::ERegMode regMode =
383 
384  // preserve the area
385  TUnfold::EConstraint constraintMode=
387 
388  // bin content is divided by the bin width
389  TUnfoldDensity::EDensityMode densityFlags=
391 
392  // set up matrix of migrations
393  TUnfoldDensity unfold(histUnfoldMatrix,TUnfold::kHistMapOutputHoriz,
394  regMode,constraintMode,densityFlags);
395 
396  // define the input vector (the measured data distribution)
397  unfold.SetInput(histUnfoldInput);
398 
399  // subtract background, normalized to data luminosity
400  // with 10% scale error each
401  Double_t scale_bgr=1.0;
402  Double_t dscale_bgr=0.1;
403  unfold.SubtractBackground(histUnfoldBgr1,"background1",scale_bgr,dscale_bgr);
404  unfold.SubtractBackground(histUnfoldBgr2,"background2",scale_bgr,dscale_bgr);
405 
406  // add systematic error
407  unfold.AddSysError(histUnfoldMatrixSys,"signalshape_SYS",
410 
411  // run the unfolding
412  Int_t nScan=30;
413  TSpline *logTauX,*logTauY;
414  TGraph *lCurve;
415  // this method scans the parameter tau and finds the kink in the L curve
416  // finally, the unfolding is done for the best choice of tau
417  Int_t iBest=unfold.ScanLcurve(nScan,0.,0.,&lCurve,&logTauX,&logTauY);
418 
419  cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
420  <<" / "<<unfold.GetNdf()<<"\n";
421 
422  // save graphs with one point to visualize best choice of tau
423  Double_t t[1],x[1],y[1];
424  logTauX->GetKnot(iBest,t[0],x[0]);
425  logTauY->GetKnot(iBest,t[0],y[0]);
426  TGraph *bestLcurve=new TGraph(1,x,y);
427  TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
428 
429 
430  //===========================
431  // Step 6: retreive unfolding results
432 
433  // get unfolding output
434  // includes the statistical and background errors
435  // but not the other systematic uncertainties
436  TH1 *histUnfoldOutput=unfold.GetOutput("PT(unfold,stat+bgrerr)");
437 
438  // retreive error matrix of statistical errors
439  TH2 *histEmatStat=unfold.GetEmatrixInput("unfolding stat error matrix");
440 
441  // retreive full error matrix
442  // This includes all systematic errors
443  TH2 *histEmatTotal=unfold.GetEmatrixTotal("unfolding total error matrix");
444 
445  // create two copies of the unfolded data, one with statistical errors
446  // the other with total errors
447  TH1 *histUnfoldStat=new TH1D("PT(unfold,staterr)",";Pt(gen)",
448  nGen,xminGen,xmaxGen);
449  TH1 *histUnfoldTotal=new TH1D("PT(unfold,totalerr)",";Pt(gen)",
450  nGen,xminGen,xmaxGen);
451  for(Int_t i=0;i<nGen+2;i++) {
452  Double_t c=histUnfoldOutput->GetBinContent(i);
453 
454  // histogram with unfolded data and stat errors
455  histUnfoldStat->SetBinContent(i,c);
456  histUnfoldStat->SetBinError
457  (i,TMath::Sqrt(histEmatStat->GetBinContent(i,i)));
458 
459  // histogram with unfolded data and total errors
460  histUnfoldTotal->SetBinContent(i,c);
461  histUnfoldTotal->SetBinError
462  (i,TMath::Sqrt(histEmatTotal->GetBinContent(i,i)));
463  }
464 
465  // create histogram with correlation matrix
466  TH2D *histCorr=new TH2D("Corr(total)",";Pt(gen);Pt(gen)",
467  nGen,xminGen,xmaxGen,nGen,xminGen,xmaxGen);
468  for(Int_t i=0;i<nGen+2;i++) {
469  Double_t ei,ej;
470  ei=TMath::Sqrt(histEmatTotal->GetBinContent(i,i));
471  if(ei<=0.0) continue;
472  for(Int_t j=0;j<nGen+2;j++) {
473  ej=TMath::Sqrt(histEmatTotal->GetBinContent(j,j));
474  if(ej<=0.0) continue;
475  histCorr->SetBinContent(i,j,histEmatTotal->GetBinContent(i,j)/ei/ej);
476  }
477  }
478 
479  // retreive bgr source 1
480  TH1 *histDetNormBgr1=unfold.GetBackground("bgr1 normalized",
481  "background1");
482  TH1 *histDetNormBgrTotal=unfold.GetBackground("bgr total normalized");
483 
484  //========================
485  // Step 7: plots
486 
487  TCanvas *output = new TCanvas();
488  output->Divide(3,2);
489  output->cd(1);
490  // data, MC prediction, background
491  histUnfoldInput->SetMinimum(0.0);
492  histUnfoldInput->Draw("E");
493  histDetMC->SetMinimum(0.0);
494  histDetMC->SetLineColor(kBlue);
495  histDetNormBgrTotal->SetLineColor(kRed);
496  histDetNormBgr1->SetLineColor(kCyan);
497  histDetMC->Draw("SAME HIST");
498  histDetNormBgr1->Draw("SAME HIST");
499  histDetNormBgrTotal->Draw("SAME HIST");
500 
501  output->cd(2);
502  // unfolded data, data truth, MC truth
503  histUnfoldTotal->SetMinimum(0.0);
504  histUnfoldTotal->SetMaximum(histUnfoldTotal->GetMaximum()*1.5);
505  // outer error: total error
506  histUnfoldTotal->Draw("E");
507  // middle error: stat+bgr
508  histUnfoldOutput->Draw("SAME E1");
509  // inner error: stat only
510  histUnfoldStat->Draw("SAME E1");
511  histDataTruth->Draw("SAME HIST");
512  histBbbSignalGen->SetLineColor(kBlue);
513  histBbbSignalGen->Draw("SAME HIST");
514 
515  output->cd(3);
516  // unfolding matrix
517  histUnfoldMatrix->SetLineColor(kBlue);
518  histUnfoldMatrix->Draw("BOX");
519 
520  // show tau as a function of chi**2
521  output->cd(4);
522  logTauX->Draw();
523  bestLogTauLogChi2->SetMarkerColor(kRed);
524  bestLogTauLogChi2->Draw("*");
525 
526  // show the L curve
527  output->cd(5);
528  lCurve->Draw("AL");
529  bestLcurve->SetMarkerColor(kRed);
530  bestLcurve->Draw("*");
531 
532  // show correlation matrix
533  output->cd(6);
534  histCorr->Draw("BOX");
535 
536  //============================================================
537  // step 8: compare results to the so-called bin-by-bin "correction"
538 
539  std::cout<<"bin truth result (stat) (bgr) (sys)\n";
540  std::cout<<"===+=====+=========+========+========+=======\n";
541  for(Int_t i=1;i<=nGen;i++) {
542  // data contribution in this bin
543  Double_t data=histBbbInput->GetBinContent(i);
544  Double_t errData=histBbbInput->GetBinError(i);
545 
546  // subtract background contributions
547  Double_t data_bgr=data
548  - scale_bgr*(histBbbBgr1->GetBinContent(i)
549  + histBbbBgr2->GetBinContent(i)
550  + histBbbBgrPt->GetBinContent(i));
551  Double_t errData_bgr=TMath::Sqrt
552  (errData*errData+
553  TMath::Power(dscale_bgr*histBbbBgr1->GetBinContent(i),2)+
554  TMath::Power(scale_bgr*histBbbBgr1->GetBinError(i),2)+
555  TMath::Power(dscale_bgr*histBbbBgr2->GetBinContent(i),2)+
556  TMath::Power(scale_bgr*histBbbBgr2->GetBinError(i),2)+
557  TMath::Power(dscale_bgr*histBbbBgrPt->GetBinContent(i),2)+
558  TMath::Power(scale_bgr*histBbbBgrPt->GetBinError(i),2));
559  // "correct" the data, using the Monte Carlo and neglecting off-diagonals
560  Double_t fCorr=(histBbbSignalGen->GetBinContent(i)/
561  histBbbSignalRec->GetBinContent(i));
562  Double_t data_bbb= data_bgr *fCorr;
563  // stat only error
564  Double_t errData_stat_bbb = errData*fCorr;
565  // stat plus background subtraction
566  Double_t errData_statbgr_bbb = errData_bgr*fCorr;
567 
568  // estimate systematic error by repeating the exercise
569  // using the MC with systematic shifts applied
570  Double_t fCorr_sys=(histBbbSignalGenSys->GetBinContent(i)/
571  histBbbSignalRecSys->GetBinContent(i));
572  Double_t shift_sys= data_bgr*fCorr_sys - data_bbb;
573 
574  // add systematic shift quadratically and get total error
575  Double_t errData_total_bbb=
576  TMath::Sqrt(errData_statbgr_bbb*errData_statbgr_bbb
577  +shift_sys*shift_sys);
578 
579  // get results from real unfolding
580  Double_t data_unfold= histUnfoldStat->GetBinContent(i);
581  Double_t errData_stat_unfold=histUnfoldStat->GetBinError(i);
582  Double_t errData_statbgr_unfold=histUnfoldOutput->GetBinError(i);
583  Double_t errData_total_unfold=histUnfoldTotal->GetBinError(i);
584 
585  // compare
586 
587  std::cout<<TString::Format
588  ("%3d %5.0f %8.1f +/-%5.1f +/-%5.1f +/-%5.1f (unfolding)",
589  i,histDataTruth->GetBinContent(i),data_unfold,
590  errData_stat_unfold,TMath::Sqrt(errData_statbgr_unfold*
591  errData_statbgr_unfold-
592  errData_stat_unfold*
593  errData_stat_unfold),
594  TMath::Sqrt(errData_total_unfold*
595  errData_total_unfold-
596  errData_statbgr_unfold*
597  errData_statbgr_unfold))<<"\n";
598  std::cout<<TString::Format
599  (" %8.1f +/-%5.1f +/-%5.1f +/-%5.1f (bin-by-bin)",
600  data_bbb,errData_stat_bbb,TMath::Sqrt(errData_statbgr_bbb*
601  errData_statbgr_bbb-
602  errData_stat_bbb*
603  errData_stat_bbb),
604  TMath::Sqrt(errData_total_bbb*
605  errData_total_bbb-
606  errData_statbgr_bbb*
607  errData_statbgr_bbb))
608  <<"\n\n";
609  }
610 }
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3125
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition: TH1.cxx:7664
Random number generator class based on M.
Definition: TRandom3.h:29
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:399
Double_t Log(Double_t x)
Definition: TMath.h:526
return c
Definition: Rtypes.h:61
R__EXTERN TStyle * gStyle
Definition: TStyle.h:418
Base class for spline implementation containing the Draw/Paint methods //.
Definition: TSpline.h:22
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4638
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:400
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:747
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
static void SetDefaultSumw2(Bool_t sumw2=kTRUE)
When this static function is called with sumw2=kTRUE, all new histograms will automatically activate ...
Definition: TH1.cxx:6012
Double_t x[n]
Definition: legend1.C:17
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2335
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:31
EConstraint
Definition: TUnfold.h:103
double pow(double, double)
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:43
const Double_t sigma
virtual void SetBinError(Int_t bin, Double_t error)
See convention for numbering bins in TH1::GetBin.
Definition: TH1.cxx:8309
ERegMode
Definition: TUnfold.h:107
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:512
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
Service class for 2-Dim histogram classes.
Definition: TH2.h:36
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2851
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8323
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1209
static const char * GetTUnfoldVersion(void)
Definition: TUnfold.cxx:289
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:618
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0
The Canvas class.
Definition: TCanvas.h:41
Double_t Exp(Double_t x)
Definition: TMath.h:495
double f(double x)
double Double_t
Definition: RtypesCore.h:55
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TSpline.cxx:100
Double_t y[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:80
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:307
Definition: Rtypes.h:61
virtual Double_t BreitWigner(Double_t mean=0, Double_t gamma=1)
Return a number distributed following a BreitWigner function with mean and gamma. ...
Definition: TRandom.cxx:186
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1089
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:301
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:90
Definition: Rtypes.h:61
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2475
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition: TH1.cxx:8173
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:296