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