Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold2.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program as an example for a user specific regularisation scheme.
5///
6/// 1. Generate Monte Carlo and Data events
7/// The events consist of:
8/// - signal
9/// - background
10///
11/// The signal is a resonance. It is generated with a Breit-Wigner,
12/// smeared by a Gaussian
13///
14/// 2. Unfold the data. The result is:
15/// - The background level
16/// - The shape of the resonance, corrected for detector effects
17///
18/// The regularisation is done on the curvature, excluding the bins
19/// near the peak.
20///
21/// 3. produce some plots
22///
23/// \macro_output
24/// \macro_code
25///
26/// **Version 17.6, in parallel to changes in TUnfold**
27///
28/// #### History:
29/// - Version 17.5, in parallel to changes in TUnfold
30/// - Version 17.4, in parallel to changes in TUnfold
31/// - Version 17.3, in parallel to changes in TUnfold
32/// - Version 17.2, in parallel to changes in TUnfold
33/// - Version 17.1, in parallel to changes in TUnfold
34/// - Version 17.0, updated for changed methods in TUnfold
35/// - Version 16.1, parallel to changes in TUnfold
36/// - Version 16.0, parallel to changes in TUnfold
37/// - Version 15, with automatic L-curve scan, simplified example
38/// - Version 14, with changes in TUnfoldSys.cxx
39/// - Version 13, with changes to TUnfold.C
40/// - Version 12, with improvements to TUnfold.cxx
41/// - Version 11, print chi**2 and number of degrees of freedom
42/// - Version 10, with bug-fix in TUnfold.cxx
43/// - Version 9, with bug-fix in TUnfold.cxx, TUnfold.h
44/// - Version 8, with bug-fix in TUnfold.cxx, TUnfold.h
45/// - Version 7, with bug-fix in TUnfold.cxx, TUnfold.h
46/// - Version 6a, fix problem with dynamic array allocation under windows
47/// - Version 6, re-include class MyUnfold in the example
48/// - Version 5, move class MyUnfold to separate files
49/// - Version 4, with bug-fix in TUnfold.C
50/// - Version 3, with bug-fix in TUnfold.C
51/// - Version 2, with changed ScanLcurve() arguments
52/// - Version 1, remove L curve analysis, use ScanLcurve() method instead
53/// - Version 0, L curve analysis included here
54///
55/// This file is part of TUnfold.
56///
57/// TUnfold is free software: you can redistribute it and/or modify
58/// it under the terms of the GNU General Public License as published by
59/// the Free Software Foundation, either version 3 of the License, or
60/// (at your option) any later version.
61///
62/// TUnfold is distributed in the hope that it will be useful,
63/// but WITHOUT ANY WARRANTY; without even the implied warranty of
64/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
65/// GNU General Public License for more details.
66///
67/// You should have received a copy of the GNU General Public License
68/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
69///
70/// \author Stefan Schmitt DESY, 14.10.2008
71
72#include <TMath.h>
73#include <TCanvas.h>
74#include <TRandom3.h>
75#include <TF1.h>
76#include <TStyle.h>
77#include <TVector.h>
78#include <TGraph.h>
79#include "TUnfold.h"
80
81
82TRandom *rnd=nullptr;
83
84// generate an event
85// output:
86// negative mass: background event
87// positive mass: signal event
88Double_t GenerateEvent(Double_t bgr, // relative fraction of background
89 Double_t mass, // peak position
90 Double_t gamma) // peak width
91{
92 Double_t t;
93 if(rnd->Rndm()>bgr) {
94 // generate signal event
95 // with positive mass
96 do {
97 do {
98 t=rnd->Rndm();
99 } while(t>=1.0);
100 t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
101 } while(t<=0.0);
102 return t;
103 } else {
104 // generate background event
105 // generate events following a power-law distribution
106 // f(E) = K * TMath::power((E0+E),N0)
107 static Double_t const E0=2.4;
108 static Double_t const N0=2.9;
109 do {
110 do {
111 t=rnd->Rndm();
112 } while(t>=1.0);
113 // the mass is returned negative
114 // In our example a convenient way to indicate it is a background event.
115 t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
116 } while(t>=0.0);
117 return t;
118 }
119}
120
121// smear the event to detector level
122// input:
123// mass on generator level (mTrue>0 !)
124// output:
125// mass on detector level
127 // smear by double-gaussian
128 static Double_t frac=0.1;
129 static Double_t wideBias=0.03;
130 static Double_t wideSigma=0.5;
131 static Double_t smallBias=0.0;
132 static Double_t smallSigma=0.1;
133 if(rnd->Rndm()>frac) {
134 return rnd->Gaus(mTrue+smallBias,smallSigma);
135 } else {
136 return rnd->Gaus(mTrue+wideBias,wideSigma);
137 }
138}
139
140int testUnfold2()
141{
142 // switch on histogram errors
144
145 // random generator
146 rnd=new TRandom3();
147
148 // data and MC luminosity, cross-section
149 Double_t const luminosityData=100000;
150 Double_t const luminosityMC=1000000;
151 Double_t const crossSection=1.0;
152
153 Int_t const nDet=250;
154 Int_t const nGen=100;
155 Double_t const xminDet=0.0;
156 Double_t const xmaxDet=10.0;
157 Double_t const xminGen=0.0;
158 Double_t const xmaxGen=10.0;
159
160 //============================================
161 // generate MC distribution
162 //
163 TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
164 TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
165 TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",nDet,xminDet,xmaxDet,
168 for(Int_t i=0;i<neventMC;i++) {
169 Double_t mGen=GenerateEvent(0.3, // relative fraction of background
170 4.0, // peak position in MC
171 0.2); // peak width in MC
173 // the generated mass is negative for background
174 // and positive for signal
175 // so it will be filled in the underflow bin
176 // this is very convenient for the unfolding:
177 // the unfolded result will contain the number of background
178 // events in the underflow bin
179
180 // generated MC distribution (for comparison only)
182 // reconstructed MC distribution (for comparison only)
184
185 // matrix describing how the generator input migrates to the
186 // reconstructed level. Unfolding input.
187 // NOTE on underflow/overflow bins:
188 // (1) the detector level under/overflow bins are used for
189 // normalisation ("efficiency" correction)
190 // in our toy example, these bins are populated from tails
191 // of the initial MC distribution.
192 // (2) the generator level underflow/overflow bins are
193 // unfolded. In this example:
194 // underflow bin: background events reconstructed in the detector
195 // overflow bin: signal events generated at masses > xmaxDet
196 // for the unfolded result these bins will be filled
197 // -> the background normalisation will be contained in the underflow bin
199 }
200
201 //============================================
202 // generate data distribution
203 //
204 TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
205 TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
207 for(Int_t i=0;i<neventData;i++) {
208 Double_t mGen=GenerateEvent(0.4, // relative fraction of background
209 3.8, // peak position
210 0.15); // peak width
212 // generated data mass for comparison plots
213 // for real data, we do not have this histogram
214 histMgenData->Fill(mGen);
215
216 // reconstructed mass, unfolding input
217 histMdetData->Fill(mDet);
218 }
219
220 //=========================================================================
221 // set up the unfolding
224 // regularisation
225 //----------------
226 // the regularisation is done on the curvature (2nd derivative) of
227 // the output distribution
228 //
229 // One has to exclude the bins near the peak of the Breit-Wigner,
230 // because there the curvature is high
231 // (and the regularisation eventually could enforce a small
232 // curvature, thus biasing result)
233 //
234 // in real life, the parameters below would have to be optimized,
235 // depending on the data peak position and width
236 // Or maybe one finds a different regularisation scheme... this is
237 // just an example...
239 Int_t nPeek=3;
241 // calculate bin number corresponding to estimated peak position
243 // offset 1.5
244 // accounts for start bin 1
245 // and rounding errors +0.5
246 +1.5);
247 // regularize output bins 1..iPeek-nPeek
248 unfold.RegularizeBins(1,1,iPeek-nPeek,regMode);
249 // regularize output bins iPeek+nPeek..nGen
250 unfold.RegularizeBins(iPeek+nPeek,1,nGen-(iPeek+nPeek),regMode);
251
252 // unfolding
253 //-----------
254
255 // set input distribution and bias scale (=0)
256 if(unfold.SetInput(histMdetData,0.0)>=10000) {
257 std::cout<<"Unfolding result may be wrong\n";
258 }
259
260 // do the unfolding here
261 Double_t tauMin=0.0;
262 Double_t tauMax=0.0;
263 Int_t nScan=30;
264 Int_t iBest;
266 TGraph *lCurve;
267 // this method scans the parameter tau and finds the kink in the L curve
268 // finally, the unfolding is done for the "best" choice of tau
270 std::cout<<"tau="<<unfold.GetTau()<<"\n";
271 std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
272 <<" / "<<unfold.GetNdf()<<"\n";
273
274 // save point corresponding to the kink in the L curve as TGraph
275 Double_t t[1],x[1],y[1];
276 logTauX->GetKnot(iBest,t[0],x[0]);
277 logTauY->GetKnot(iBest,t[0],y[0]);
278 TGraph *bestLcurve=new TGraph(1,x,y);
279 TGraph *bestLogTauX=new TGraph(1,t,x);
280
281 //============================================================
282 // extract unfolding results into histograms
283
284 // set up a bin map, excluding underflow and overflow bins
285 // the binMap relates the output of the unfolding to the final
286 // histogram bins
287 Int_t *binMap=new Int_t[nGen+2];
288 for(Int_t i=1;i<=nGen;i++) binMap[i]=i;
289 binMap[0]=-1;
290 binMap[nGen+1]=-1;
291
292 TH1D *histMunfold=new TH1D("Unfolded",";mass(gen)",nGen,xminGen,xmaxGen);
293 unfold.GetOutput(histMunfold,binMap);
294 TH1D *histMdetFold=new TH1D("FoldedBack","mass(det)",nDet,xminDet,xmaxDet);
295 unfold.GetFoldedOutput(histMdetFold);
296
297 // store global correlation coefficients
298 TH1D *histRhoi=new TH1D("rho_I","mass",nGen,xminGen,xmaxGen);
299 unfold.GetRhoI(histRhoi,binMap);
300
301 delete[] binMap;
302 binMap=nullptr;
303
304 //=====================================================================
305 // plot some histograms
307
308 // produce some plots
309 output.Divide(3,2);
310
311 // Show the matrix which connects input and output
312 // There are overflow bins at the bottom, not shown in the plot
313 // These contain the background shape.
314 // The overflow bins to the left and right contain
315 // events which are not reconstructed. These are necessary for proper MC
316 // normalisation
317 output.cd(1);
318 histMdetGenMC->Draw("BOX");
319
320 // draw generator-level distribution:
321 // data (red) [for real data this is not available]
322 // MC input (black) [with completely wrong peak position and shape]
323 // unfolded data (blue)
324 output.cd(2);
325 histMunfold->SetLineColor(kBlue);
326 histMunfold->Draw();
327 histMgenData->SetLineColor(kRed);
328 histMgenData->Draw("SAME");
329 histMgenMC->Draw("SAME HIST");
330
331 // show detector level distributions
332 // data (red)
333 // MC (black)
334 // unfolded data (blue)
335 output.cd(3);
336 histMdetFold->SetLineColor(kBlue);
337 histMdetFold->Draw();
338 histMdetData->SetLineColor(kRed);
339 histMdetData->Draw("SAME");
340 histMdetMC->Draw("SAME HIST");
341
342 // show correlation coefficients
343 // all bins outside the peak are found to be highly correlated
344 // But they are compatible with zero anyway
345 // If the peak shape is fitted,
346 // these correlations have to be taken into account, see example
347 output.cd(4);
348 histRhoi->Draw();
349
350 // show rhoi_max(tau) distribution
351 output.cd(5);
352 logTauX->Draw();
353 bestLogTauX->SetMarkerColor(kRed);
354 bestLogTauX->Draw("*");
355
356 output.cd(6);
357 lCurve->Draw("AL");
358 bestLcurve->SetMarkerColor(kRed);
359 bestLcurve->Draw("*");
360
361 output.SaveAs("testUnfold2.ps");
362 return 0;
363}
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
@ kRed
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:926
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
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
An algorithm to unfold distributions from detector to truth level.
Definition TUnfold.h:107
ERegMode
choice of regularisation scheme
Definition TUnfold.h:123
@ kRegModeNone
no regularisation, or defined later by RegularizeXXX() methods
Definition TUnfold.h:126
@ kRegModeCurvature
regularize the 2nd derivative of the output distribution
Definition TUnfold.h:135
@ kHistMapOutputVert
truth level on y-axis of the response matrix
Definition TUnfold.h:149
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
double gamma(double x)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:732
constexpr Double_t Pi()
Definition TMath.h:40
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:611
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
static void output()