Logo ROOT  
Reference Guide
testUnfold5d.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program for the classes TUnfoldDensity and TUnfoldBinning.
5///
6/// A toy test of the TUnfold package
7///
8/// This is an example of unfolding a two-dimensional distribution
9/// also using an auxiliary measurement to constrain some background
10///
11/// The example comprises several macros
12/// testUnfold5a.C create root files with TTree objects for
13/// signal, background and data
14/// -> write files testUnfold5_signal.root
15/// testUnfold5_background.root
16/// testUnfold5_data.root
17///
18/// testUnfold5b.C create a root file with the TUnfoldBinning objects
19/// -> write file testUnfold5_binning.root
20///
21/// testUnfold5c.C loop over trees and fill histograms based on the
22/// TUnfoldBinning objects
23/// -> read testUnfold5_binning.root
24/// testUnfold5_signal.root
25/// testUnfold5_background.root
26/// testUnfold5_data.root
27///
28/// -> write testUnfold5_histograms.root
29///
30/// testUnfold5d.C run the unfolding
31/// -> read testUnfold5_histograms.root
32/// -> write testUnfold5_result.root
33/// testUnfold5_result.ps
34///
35/// \macro_output
36/// \macro_code
37///
38/// **Version 17.6, in parallel to changes in TUnfold**
39///
40/// #### History:
41/// - Version 17.5, in parallel to changes in TUnfold
42/// - Version 17.4, in parallel to changes in TUnfold
43/// - Version 17.3, in parallel to changes in TUnfold
44/// - Version 17.2, in parallel to changes in TUnfold
45/// - Version 17.1, in parallel to changes in TUnfold
46/// - Version 17.0 example for multi-dimensional unfolding
47///
48/// This file is part of TUnfold.
49///
50/// TUnfold is free software: you can redistribute it and/or modify
51/// it under the terms of the GNU General Public License as published by
52/// the Free Software Foundation, either version 3 of the License, or
53/// (at your option) any later version.
54///
55/// TUnfold is distributed in the hope that it will be useful,
56/// but WITHOUT ANY WARRANTY; without even the implied warranty of
57/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
58/// GNU General Public License for more details.
59///
60/// You should have received a copy of the GNU General Public License
61/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
62///
63/// \author Stefan Schmitt DESY, 14.10.2008
64
65#include <iostream>
66#include <cmath>
67#include <map>
68#include <TMath.h>
69#include <TCanvas.h>
70#include <TStyle.h>
71#include <TGraph.h>
72#include <TFile.h>
73#include <TH1.h>
74#include "TUnfoldDensity.h"
75
76using namespace std;
77
78// #define PRINT_MATRIX_L
79
80#define TEST_INPUT_COVARIANCE
81
82void testUnfold5d()
83{
84 // switch on histogram errors
86
87 //==============================================
88 // step 1 : open output file
89 TFile *outputFile=new TFile("testUnfold5_results.root","recreate");
90
91 //==============================================
92 // step 2 : read binning schemes and input histograms
93 TFile *inputFile=new TFile("testUnfold5_histograms.root");
94
95 outputFile->cd();
96
97 TUnfoldBinning *detectorBinning,*generatorBinning;
98
99 inputFile->GetObject("detector",detectorBinning);
100 inputFile->GetObject("generator",generatorBinning);
101
102 if((!detectorBinning)||(!generatorBinning)) {
103 cout<<"problem to read binning schemes\n";
104 }
105
106 // save binning schemes to output file
107 detectorBinning->Write();
108 generatorBinning->Write();
109
110 // read histograms
111 TH1 *histDataReco,*histDataTruth;
112 TH2 *histMCGenRec;
113
114 inputFile->GetObject("histDataReco",histDataReco);
115 inputFile->GetObject("histDataTruth",histDataTruth);
116 inputFile->GetObject("histMCGenRec",histMCGenRec);
117
118#ifdef TEST_ZERO_UNCORR_ERROR
119 // special test (bug in version 17.2 and below)
120 // set all errors in hisMCGenRec to zero
121 // -> program will crash
122 for(int i=0;i<=histMCGenRec->GetNbinsX()+1;i++) {
123 for(int j=0;j<=histMCGenRec->GetNbinsY()+1;j++) {
124 histMCGenRec->SetBinError(i,j,0.0);
125 }
126 }
127#endif
128
129 histDataReco->Write();
130 histDataTruth->Write();
131 histMCGenRec->Write();
132
133 if((!histDataReco)||(!histDataTruth)||(!histMCGenRec)) {
134 cout<<"problem to read input histograms\n";
135 }
136
137 //========================
138 // Step 3: unfolding
139
140 // preserve the area
142
143 // basic choice of regularisation scheme:
144 // curvature (second derivative)
146
147 // density flags
148 TUnfoldDensity::EDensityMode densityFlags=
150
151 // detailed steering for regularisation
152 const char *REGULARISATION_DISTRIBUTION=0;
153 const char *REGULARISATION_AXISSTEERING="*[B]";
154
155 // set up matrix of migrations
157 regMode,constraintMode,densityFlags,
158 generatorBinning,detectorBinning,
159 REGULARISATION_DISTRIBUTION,
160 REGULARISATION_AXISSTEERING);
161
162 // define the input vector (the measured data distribution)
163
164#ifdef TEST_INPUT_COVARIANCE
165 // special test to use input covariance matrix
166 TH2D *inputEmatrix=
167 detectorBinning->CreateErrorMatrixHistogram("input_covar",true);
168 for(int i=1;i<=inputEmatrix->GetNbinsX();i++) {
169 Double_t e=histDataReco->GetBinError(i);
170 inputEmatrix->SetBinContent(i,i,e*e);
171 // test: non-zero covariance where variance is zero
172 // if(e<=0.) inputEmatrix->SetBinContent(i,i+1,1.0);
173 }
174 unfold.SetInput(histDataReco,0.0,0.0,inputEmatrix);
175#else
176 unfold.SetInput(histDataReco /* ,0.0,1.0 */);
177#endif
178 // print matrix of regularisation conditions
179#ifdef PRINT_MATRIX_L
180 TH2 *histL= unfold.GetL("L");
181 for(Int_t j=1;j<=histL->GetNbinsY();j++) {
182 cout<<"L["<<unfold.GetLBinning()->GetBinName(j)<<"]";
183 for(Int_t i=1;i<=histL->GetNbinsX();i++) {
184 Double_t c=histL->GetBinContent(i,j);
185 if(c!=0.0) cout<<" ["<<i<<"]="<<c;
186 }
187 cout<<"\n";
188 }
189#endif
190 // run the unfolding
191 //
192 // here, tau is determined by scanning the global correlation coefficients
193
194 Int_t nScan=30;
195 TSpline *rhoLogTau=0;
196 TGraph *lCurve=0;
197
198 // for determining tau, scan the correlation coefficients
199 // correlation coefficients may be probed for all distributions
200 // or only for selected distributions
201 // underflow/overflow bins may be included/excluded
202 //
203 const char *SCAN_DISTRIBUTION="signal";
204 const char *SCAN_AXISSTEERING=0;
205
206 Int_t iBest=unfold.ScanTau(nScan,0.,0.,&rhoLogTau,
208 SCAN_DISTRIBUTION,SCAN_AXISSTEERING,
209 &lCurve);
210
211 // create graphs with one point to visualize best choice of tau
212 Double_t t[1],rho[1],x[1],y[1];
213 rhoLogTau->GetKnot(iBest,t[0],rho[0]);
214 lCurve->GetPoint(iBest,x[0],y[0]);
215 TGraph *bestRhoLogTau=new TGraph(1,t,rho);
216 TGraph *bestLCurve=new TGraph(1,x,y);
217 Double_t *tAll=new Double_t[nScan],*rhoAll=new Double_t[nScan];
218 for(Int_t i=0;i<nScan;i++) {
219 rhoLogTau->GetKnot(i,tAll[i],rhoAll[i]);
220 }
221 TGraph *knots=new TGraph(nScan,tAll,rhoAll);
222
223 cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
224 <<" / "<<unfold.GetNdf()<<"\n";
225
226
227 //===========================
228 // Step 4: retrieve and plot unfolding results
229
230 // get unfolding output
231 TH1 *histDataUnfold=unfold.GetOutput("unfolded signal",0,0,0,kFALSE);
232 // get Monte Carlo reconstructed data
233 TH1 *histMCReco=histMCGenRec->ProjectionY("histMCReco",0,-1,"e");
234 TH1 *histMCTruth=histMCGenRec->ProjectionX("histMCTruth",0,-1,"e");
235 Double_t scaleFactor=histDataTruth->GetSumOfWeights()/
236 histMCTruth->GetSumOfWeights();
237 histMCReco->Scale(scaleFactor);
238 histMCTruth->Scale(scaleFactor);
239 // get matrix of probabilities
240 TH2 *histProbability=unfold.GetProbabilityMatrix("histProbability");
241 // get global correlation coefficients
242 /* TH1 *histGlobalCorr=*/ unfold.GetRhoItotal("histGlobalCorr",0,0,0,kFALSE);
243 TH1 *histGlobalCorrScan=unfold.GetRhoItotal
244 ("histGlobalCorrScan",0,SCAN_DISTRIBUTION,SCAN_AXISSTEERING,kFALSE);
245 /* TH2 *histCorrCoeff=*/ unfold.GetRhoIJtotal("histCorrCoeff",0,0,0,kFALSE);
246
247 TCanvas canvas;
248 canvas.Print("testUnfold5.ps[");
249
250 //========== page 1 ============
251 // unfolding control plots
252 // input, matrix, output
253 // tau-scan, global correlations, correlation coefficients
254 canvas.Clear();
255 canvas.Divide(3,2);
256
257 // (1) all bins, compare to original MC distribution
258 canvas.cd(1);
259 histDataReco->SetMinimum(0.0);
260 histDataReco->Draw("E");
261 histMCReco->SetLineColor(kBlue);
262 histMCReco->Draw("SAME HIST");
263 // (2) matrix of probabilities
264 canvas.cd(2);
265 histProbability->Draw("BOX");
266 // (3) unfolded data, data truth, MC truth
267 canvas.cd(3);
268 gPad->SetLogy();
269 histDataUnfold->Draw("E");
270 histDataTruth->SetLineColor(kBlue);
271 histDataTruth->Draw("SAME HIST");
272 histMCTruth->SetLineColor(kRed);
273 histMCTruth->Draw("SAME HIST");
274 // (4) scan of correlation vs tau
275 canvas.cd(4);
276 rhoLogTau->Draw();
277 knots->Draw("*");
278 bestRhoLogTau->SetMarkerColor(kRed);
279 bestRhoLogTau->Draw("*");
280 // (5) global correlation coefficients for the distributions
281 // used during the scan
282 canvas.cd(5);
283 //histCorrCoeff->Draw("BOX");
284 histGlobalCorrScan->Draw("HIST");
285 // (6) L-curve
286 canvas.cd(6);
287 lCurve->Draw("AL");
288 bestLCurve->SetMarkerColor(kRed);
289 bestLCurve->Draw("*");
290
291
292 canvas.Print("testUnfold5.ps");
293
294 canvas.Print("testUnfold5.ps]");
295
296}
#define c(i)
Definition: RSha256.hxx:101
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
@ kRed
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
#define gPad
Definition: TVirtualPad.h:286
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:31
void Clear(Option_t *option="")
Remove all primitives from the canvas.
Definition: TCanvas.cxx:715
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:696
Bool_t cd(const char *path=nullptr) override
Change current directory to "this" directory.
void GetObject(const char *namecycle, T *&ptr)
Definition: TDirectory.h:149
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:753
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition: TGraph.cxx:1586
The TH1 histogram class.
Definition: TH1.h:56
virtual Int_t GetNbinsY() const
Definition: TH1.h:293
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition: TH1.cxx:8507
virtual Int_t GetNbinsX() const
Definition: TH1.h:292
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:8650
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:395
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:6318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6234
virtual Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition: TH1.cxx:7414
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:292
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
TH1D * ProjectionY(const char *name="_py", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
Project a 2-D histogram into a 1-D histogram along Y.
Definition: TH2.cxx:2301
TH1D * ProjectionX(const char *name="_px", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
Project a 2-D histogram into a 1-D histogram along X.
Definition: TH2.cxx:2261
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:88
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2441
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition: TObject.cxx:785
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:1163
virtual void Print(const char *filename="") const
Save Pad contents in a file in one of various formats.
Definition: TPad.cxx:4657
Base class for spline implementation containing the Draw/Paint methods.
Definition: TSpline.h:22
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TSpline.cxx:97
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0
Binning schemes for use with the unfolding algorithm TUnfoldDensity.
TH2D * CreateErrorMatrixHistogram(const char *histogramName, Bool_t originalAxisBinning, Int_t **binMap=0, const char *histogramTitle=0, const char *axisSteering=0) const
Create a TH2D histogram capable to hold a covariance matrix.
An algorithm to unfold distributions from detector to truth level.
@ kEScanTauRhoMax
maximum global correlation coefficient (from TUnfold::GetRhoI())
EDensityMode
choice of regularisation scale factors to cinstruct the matrix L
@ kDensityModeBinWidth
scale factors from multidimensional bin width
EConstraint
type of extra constraint
Definition: TUnfold.h:110
@ kEConstraintArea
enforce preservation of the area
Definition: TUnfold.h:116
ERegMode
choice of regularisation scheme
Definition: TUnfold.h:120
@ kRegModeCurvature
regularize the 2nd derivative of the output distribution
Definition: TUnfold.h:132
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
Definition: TUnfold.h:143
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17