ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
testUnfold5d.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_unfold5
3 ///
4 /// Version 17.0 example for multi-dimensional unfolding
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \author Stefan Schmitt, DESY
11 #include <iostream>
12 #include <cmath>
13 #include <map>
14 #include <TMath.h>
15 #include <TCanvas.h>
16 #include <TStyle.h>
17 #include <TGraph.h>
18 #include <TFile.h>
19 #include <TH1.h>
20 #include "TUnfoldDensity.h"
21 
22 using namespace std;
23 
24 // #define PRINT_MATRIX_L
25 
26 void testUnfold5d()
27 {
28  // switch on histogram errors
30 
31  //==============================================
32  // step 1 : open output file
33  TFile *outputFile=new TFile("testUnfold5_results.root","recreate");
34 
35  //==============================================
36  // step 2 : read binning schemes and input histograms
37  TFile *inputFile=new TFile("testUnfold5_histograms.root");
38 
39  outputFile->cd();
40 
41  TUnfoldBinning *detectorBinning,*generatorBinning;
42 
43  inputFile->GetObject("detector",detectorBinning);
44  inputFile->GetObject("generator",generatorBinning);
45 
46  if((!detectorBinning)||(!generatorBinning)) {
47  cout<<"problem to read binning schemes\n";
48  }
49 
50  // save binning schemes to output file
51  detectorBinning->Write();
52  generatorBinning->Write();
53 
54  // read histograms
55  TH1 *histDataReco,*histDataTruth;
56  TH2 *histMCGenRec;
57 
58  inputFile->GetObject("histDataReco",histDataReco);
59  inputFile->GetObject("histDataTruth",histDataTruth);
60  inputFile->GetObject("histMCGenRec",histMCGenRec);
61 
62  histDataReco->Write();
63  histDataTruth->Write();
64  histMCGenRec->Write();
65 
66  if((!histDataReco)||(!histDataTruth)||(!histMCGenRec)) {
67  cout<<"problem to read input histograms\n";
68  }
69 
70  //========================
71  // Step 3: unfolding
72 
73  // preserve the area
75 
76  // basic choice of regularisation scheme:
77  // curvature (second derivative)
79 
80  // density flags
81  TUnfoldDensity::EDensityMode densityFlags=
83 
84  // detailed steering for regularisation
85  const char *REGULARISATION_DISTRIBUTION=0;
86  const char *REGULARISATION_AXISSTEERING="*[B]";
87 
88  // set up matrix of migrations
89  TUnfoldDensity unfold(histMCGenRec,TUnfold::kHistMapOutputHoriz,
90  regMode,constraintMode,densityFlags,
91  generatorBinning,detectorBinning,
92  REGULARISATION_DISTRIBUTION,
93  REGULARISATION_AXISSTEERING);
94 
95  // define the input vector (the measured data distribution)
96  unfold.SetInput(histDataReco /* ,0.0,1.0 */);
97 
98  // print matrix of regularisation conditions
99 #ifdef PRINT_MATRIX_L
100  TH2 *histL= unfold.GetL("L");
101  for(Int_t j=1;j<=histL->GetNbinsY();j++) {
102  cout<<"L["<<unfold.GetLBinning()->GetBinName(j)<<"]";
103  for(Int_t i=1;i<=histL->GetNbinsX();i++) {
104  Double_t c=histL->GetBinContent(i,j);
105  if(c!=0.0) cout<<" ["<<i<<"]="<<c;
106  }
107  cout<<"\n";
108  }
109 #endif
110  // run the unfolding
111  //
112  // here, tau is determined by scanning the global correlation coefficients
113 
114  Int_t nScan=30;
115  TSpline *rhoLogTau=0;
116  TGraph *lCurve=0;
117 
118  // for determining tau, scan the correlation coefficients
119  // correlation coefficients may be probed for all distributions
120  // or only for selected distributions
121  // underflow/overflow bins may be included/excluded
122  //
123  const char *SCAN_DISTRIBUTION="signal";
124  const char *SCAN_AXISSTEERING=0;
125 
126  Int_t iBest=unfold.ScanTau(nScan,0.,0.,&rhoLogTau,
128  SCAN_DISTRIBUTION,SCAN_AXISSTEERING,
129  &lCurve);
130 
131  // create graphs with one point to visualize best choice of tau
132  Double_t t[1],rho[1],x[1],y[1];
133  rhoLogTau->GetKnot(iBest,t[0],rho[0]);
134  lCurve->GetPoint(iBest,x[0],y[0]);
135  TGraph *bestRhoLogTau=new TGraph(1,t,rho);
136  TGraph *bestLCurve=new TGraph(1,x,y);
137  Double_t *tAll=new Double_t[nScan],*rhoAll=new Double_t[nScan];
138  for(Int_t i=0;i<nScan;i++) {
139  rhoLogTau->GetKnot(i,tAll[i],rhoAll[i]);
140  }
141  TGraph *knots=new TGraph(nScan,tAll,rhoAll);
142 
143  cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
144  <<" / "<<unfold.GetNdf()<<"\n";
145 
146 
147  //===========================
148  // Step 4: retreive and plot unfolding results
149 
150  // get unfolding output
151  TH1 *histDataUnfold=unfold.GetOutput("unfolded signal",0,0,0,kFALSE);
152  // get MOnte Carlo reconstructed data
153  TH1 *histMCReco=histMCGenRec->ProjectionY("histMCReco",0,-1,"e");
154  TH1 *histMCTruth=histMCGenRec->ProjectionX("histMCTruth",0,-1,"e");
155  Double_t scaleFactor=histDataTruth->GetSumOfWeights()/
156  histMCTruth->GetSumOfWeights();
157  histMCReco->Scale(scaleFactor);
158  histMCTruth->Scale(scaleFactor);
159  // get matrix of probabilities
160  TH2 *histProbability=unfold.GetProbabilityMatrix("histProbability");
161  // get global correlation coefficients
162  TH1 *histGlobalCorr=unfold.GetRhoItotal("histGlobalCorr",0,0,0,kFALSE);
163  TH1 *histGlobalCorrScan=unfold.GetRhoItotal
164  ("histGlobalCorrScan",0,SCAN_DISTRIBUTION,SCAN_AXISSTEERING,kFALSE);
165  TH2 *histCorrCoeff=unfold.GetRhoIJtotal("histCorrCoeff",0,0,0,kFALSE);
166 
167  TCanvas *canvas = new TCanvas();
168  canvas->Print("testUnfold5.ps[");
169 
170  //========== page 1 ============
171  // unfolding control plots
172  // input, matrix, output
173  // tau-scan, global correlations, correlation coefficients
174  canvas->Clear();
175  canvas->Divide(3,2);
176 
177  // (1) all bins, compare to original MC distribution
178  canvas->cd(1);
179  histDataReco->SetMinimum(0.0);
180  histDataReco->Draw("E");
181  histMCReco->SetLineColor(kBlue);
182  histMCReco->Draw("SAME HIST");
183  // (2) matrix of probabilities
184  canvas->cd(2);
185  histProbability->Draw("BOX");
186  // (3) unfolded data, data truth, MC truth
187  canvas->cd(3);
188  gPad->SetLogy();
189  histDataUnfold->Draw("E");
190  histDataTruth->SetLineColor(kBlue);
191  histDataTruth->Draw("SAME HIST");
192  histMCTruth->SetLineColor(kRed);
193  histMCTruth->Draw("SAME HIST");
194  // (4) scan of correlation vs tau
195  canvas->cd(4);
196  rhoLogTau->Draw();
197  knots->Draw("*");
198  bestRhoLogTau->SetMarkerColor(kRed);
199  bestRhoLogTau->Draw("*");
200  // (5) global correlation coefficients for the distributions
201  // used during the scan
202  canvas->cd(5);
203  //histCorrCoeff->Draw("BOX");
204  histGlobalCorrScan->Draw("HIST");
205  // (6) L-curve
206  canvas->cd(6);
207  lCurve->Draw("AL");
208  bestLCurve->SetMarkerColor(kRed);
209  bestLCurve->Draw("*");
210 
211 
212  canvas->Print("testUnfold5.ps");
213 
214  canvas->Print("testUnfold5.ps]");
215 
216 }
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:823
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
return c
Definition: Rtypes.h:61
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
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
static void SetDefaultSumw2(Bool_t sumw2=kTRUE)
static function.
Definition: TH1.cxx:6253
EConstraint
Definition: TUnfold.h:103
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:2523
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:2563
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 Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition: TH1.cxx:7354
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0
This class serves as a container of analysis bins analysis bins are specified by defining the axes of...
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
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
#define gPad
Definition: TVirtualPad.h:288
Definition: Rtypes.h:61