Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold4.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program for the class TUnfoldSys.
5///
6/// Simple toy tests of the TUnfold package
7///
8/// Pseudo data (5000 events) are unfolded into three components
9/// The unfolding is performed once without and once with area constraint
10///
11/// Ideally, the pulls may show that the result is biased if no constraint
12/// is applied. This is expected because the true data errors are not known,
13/// and instead the sqrt(data) errors are used.
14///
15/// \macro_output
16/// \macro_code
17///
18/// **Version 17.6, in parallel to changes in TUnfold**
19///
20/// #### History:
21/// - Version 17.5, in parallel to changes in TUnfold
22/// - Version 17.4, in parallel to changes in TUnfold
23/// - Version 17.3, in parallel to changes in TUnfold
24/// - Version 17.2, in parallel to changes in TUnfold
25/// - Version 17.1, in parallel to changes in TUnfold
26/// - Version 16.1, parallel to changes in TUnfold
27/// - Version 16.0, parallel to changes in TUnfold
28/// - Version 15, use L-curve scan to scan the average correlation
29///
30/// This file is part of TUnfold.
31///
32/// TUnfold is free software: you can redistribute it and/or modify
33/// it under the terms of the GNU General Public License as published by
34/// the Free Software Foundation, either version 3 of the License, or
35/// (at your option) any later version.
36///
37/// TUnfold is distributed in the hope that it will be useful,
38/// but WITHOUT ANY WARRANTY; without even the implied warranty of
39/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
40/// GNU General Public License for more details.
41///
42/// You should have received a copy of the GNU General Public License
43/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
44///
45/// \author Stefan Schmitt DESY, 14.10.2008
46
47#include <TMath.h>
48#include <TCanvas.h>
49#include <TRandom3.h>
50#include <TF1.h>
51#include <TStyle.h>
52#include <TVector.h>
53#include <TGraph.h>
54#include <TError.h>
55
56#include "TUnfoldDensity.h"
57
58using std::cout;
59
60TRandom *rnd=nullptr;
61
63 // choose an integer random number in the range [0,nmax]
64 // (the generator level bin)
65 // depending on the probabilities
66 // probability[0],probability[1],...probability[nmax-1]
67 Double_t f=rnd->Rndm();
68 Int_t r=0;
69 while((r<nmax)&&(f>=probability[r])) {
70 f -= probability[r];
71 r++;
72 }
73 return r;
74}
75
77 // return a coordinate (the reconstructed variable)
78 // depending on shapeParm[]
79 // shapeParm[0]: fraction of events with Gaussian distribution
80 // shapeParm[1]: mean of Gaussian
81 // shapeParm[2]: width of Gaussian
82 // (1-shapeParm[0]): fraction of events with flat distribution
83 // shapeParm[3]: minimum of flat component
84 // shapeParm[4]: maximum of flat component
85 Double_t f=rnd->Rndm();
86 Double_t r;
87 if(f<shapeParm[0]) {
88 r=rnd->Gaus(shapeParm[1],shapeParm[2]);
89 } else {
90 r=rnd->Rndm()*(shapeParm[4]-shapeParm[3])+shapeParm[3];
91 }
92 return r;
93}
94
95void testUnfold4(bool printInfo = false)
96{
97
98 // switch off printing Info messages
100
101 // switch on histogram errors
103
104 // random generator
105 rnd=new TRandom3();
106
107 // data and MC number of events
108 Double_t const nData0= 500.0;
109 Double_t const nMC0 = 50000.0;
110
111 // Binning
112 // reconstructed variable (0-10)
113 Int_t const nDet=15;
114 Double_t const xminDet=0.0;
115 Double_t const xmaxDet=15.0;
116
117 // signal binning (three shapes: 0,1,2)
118 Int_t const nGen=3;
119 Double_t const xminGen=-0.5;
120 Double_t const xmaxGen= 2.5;
121
122 // parameters
123 // fraction of events per signal shape
124 static const Double_t genFrac[]={0.3,0.6,0.1};
125
126 // signal shapes
127 static const Double_t genShape[][5]=
128 {{1.0,2.0,1.5,0.,15.},
129 {1.0,7.0,2.5,0.,15.},
130 {0.0,0.0,0.0,0.,15.}};
131
132 // define DATA histograms
133 // observed data distribution
134 TH1D *histDetDATA=new TH1D("Yrec",";DATA(Yrec)",nDet,xminDet,xmaxDet);
135
136 // define MC histograms
137 // matrix of migrations
138 TH2D *histGenDetMC=new TH2D("Yrec%Xgen","MC(Xgen,Yrec)",
140
141 TH1D *histUnfold=new TH1D("Xgen",";DATA(Xgen)",nGen,xminGen,xmaxGen);
142
143 TH1D **histPullNC=new TH1D* [nGen];
144 TH1D **histPullArea=new TH1D* [nGen];
145 for(int i=0;i<nGen;i++) {
146 histPullNC[i]=new TH1D(TString::Format("PullNC%d",i),"pull",15,-3.,3.);
147 histPullArea[i]=new TH1D(TString::Format("PullArea%d",i),"pull",15,-3.,3.);
148 }
149
150 // this method is new in version 16 of TUnfold
151 cout<<"TUnfold version is "<<TUnfold::GetTUnfoldVersion()<<"\n";
152
153 for(int itoy=0;itoy<1000;itoy++) {
154 if(!(itoy %10)) cout<<"toy iteration: "<<itoy<<"\n";
155 histDetDATA->Reset();
156 histGenDetMC->Reset();
157
158 Int_t nData=rnd->Poisson(nData0);
159 for(Int_t i=0;i<nData;i++) {
162 histDetDATA->Fill(yObs);
163 }
164
165 Int_t nMC=rnd->Poisson(nMC0);
166 for(Int_t i=0;i<nMC;i++) {
169 histGenDetMC->Fill(iGen,yObs);
170 }
171 /* for(Int_t ix=0;ix<=histGenDetMC->GetNbinsX()+1;ix++) {
172 for(Int_t iy=0;iy<=histGenDetMC->GetNbinsY()+1;iy++) {
173 cout<<ix<<iy<<" : "<<histGenDetMC->GetBinContent(ix,iy)<<"\n";
174 }
175 } */
176 //========================
177 // unfolding
178
181 // define the input vector (the measured data distribution)
182 unfold.SetInput(histDetDATA,0.0,1.0);
183
184 // run the unfolding
185 unfold.ScanLcurve(50,0.,0.,nullptr,nullptr,nullptr);
186
187 // fill pull distributions without constraint
188 unfold.GetOutput(histUnfold);
189
190 for(int i=0;i<nGen;i++) {
191 histPullNC[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
192 histUnfold->GetBinError(i+1));
193 }
194
195 // repeat unfolding on the same data, now with Area constraint
196 unfold.SetConstraint(TUnfold::kEConstraintArea);
197
198 // run the unfolding
199 unfold.ScanLcurve(50,0.,0.,nullptr,nullptr,nullptr);
200
201 // fill pull distributions with constraint
202 unfold.GetOutput(histUnfold);
203
204 for(int i=0;i<nGen;i++) {
205 histPullArea[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
206 histUnfold->GetBinError(i+1));
207 }
208
209 }
211 output.Divide(3,2);
212
213 gStyle->SetOptFit(1111);
214
215 for(int i=0;i<nGen;i++) {
216 output.cd(i+1);
217 histPullNC[i]->Fit("gaus");
218 histPullNC[i]->Draw();
219 }
220 for(int i=0;i<nGen;i++) {
221 output.cd(i+4);
222 histPullArea[i]->Fit("gaus");
223 histPullArea[i]->Draw();
224 }
225 output.SaveAs("testUnfold4.ps");
226}
#define f(i)
Definition RSha256.hxx:104
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
constexpr Int_t kWarning
Definition TError.h:46
Int_t gErrorIgnoreLevel
errors with level below this value will be ignored. Default is kUnset.
Definition TError.cxx:33
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:927
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:6751
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:400
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
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:2384
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:1594
An algorithm to unfold distributions from detector to truth level, with background subtraction and pr...
Definition TUnfoldSys.h:59
static const char * GetTUnfoldVersion(void)
return a string describing the TUnfold version
Definition TUnfold.cxx:3717
@ kEConstraintArea
enforce preservation of the area
Definition TUnfold.h:119
@ kEConstraintNone
use no extra constraint
Definition TUnfold.h:116
@ kRegModeSize
regularise the amplitude of the output distribution
Definition TUnfold.h:129
@ kHistMapOutputHoriz
truth level on x-axis of the response matrix
Definition TUnfold.h:146
static void output()