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