Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
testUnfold1.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_unfold
3/// \notebook
4/// Test program for the classes TUnfold and related.
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/// Systematic errors from the MC shape variation are included
19/// and propagated to the result
20///
21/// 3. fit the unfolded distribution, including the correlation matrix
22///
23/// 4. save six plots to a file `testUnfold1.ps`
24/// ~~~
25/// 1 2 3
26/// 4 5 6
27/// ~~~
28/// 1. 2d-plot of the matrix describing the migrations
29/// 2. generator-level distributions
30/// - blue: unfolded data, total errors
31/// - green: unfolded data, statistical errors
32/// - red: generated data
33/// - black: fit to green data points
34/// 3. detector level distributions
35/// - blue: unfolded data, folded back through the matrix
36/// - black: Monte Carlo (with wrong peal position)
37/// - blue: data
38/// 4. global correlation coefficients
39/// 5. \f$ \chi^2 \f$ as a function of \f$ log(\tau) \f$
40/// the star indicates the final choice of \f$ \tau \f$
41/// 6. the L curve
42///
43/// \macro_output
44/// \macro_image
45/// \macro_code
46///
47/// **Version 17.6, in parallel to changes in TUnfold**
48///
49/// #### History:
50/// - Version 17.5, in parallel to changes in TUnfold
51/// - Version 17.4, in parallel to changes in TUnfold
52/// - Version 17.3, in parallel to changes in TUnfold
53/// - Version 17.2, in parallel to changes in TUnfold
54/// - Version 17.1, in parallel to changes in TUnfold
55/// - Version 17.0, updated for using the classes TUnfoldDensity, TUnfoldBinning
56/// - Version 16.1, parallel to changes in TUnfold
57/// - Version 16.0, parallel to changes in TUnfold
58/// - Version 15, with automated L-curve scan
59/// - Version 14, with changes in TUnfoldSys.cxx
60/// - Version 13, include test of systematic errors
61/// - Version 12, catch error when defining the input
62/// - Version 11, print chi**2 and number of degrees of freedom
63/// - Version 10, with bug-fix in TUnfold.cxx
64/// - Version 9, with bug-fix in TUnfold.cxx and TUnfold.h
65/// - Version 8, with bug-fix in TUnfold.cxx and TUnfold.h
66/// - Version 7, with bug-fix in TUnfold.cxx and TUnfold.h
67/// - Version 6a, fix problem with dynamic array allocation under windows
68/// - Version 6, bug-fixes in TUnfold.C
69/// - Version 5, replace main() by testUnfold1()
70/// - Version 4, with bug-fix in TUnfold.C
71/// - Version 3, with bug-fix in TUnfold.C
72/// - Version 2, with changed ScanLcurve() arguments
73/// - Version 1, remove L curve analysis, use ScanLcurve() method instead
74/// - Version 0, L curve analysis included here
75///
76/// This file is part of TUnfold.
77///
78/// TUnfold is free software: you can redistribute it and/or modify
79/// it under the terms of the GNU General Public License as published by
80/// the Free Software Foundation, either version 3 of the License, or
81/// (at your option) any later version.
82///
83/// TUnfold is distributed in the hope that it will be useful,
84/// but WITHOUT ANY WARRANTY; without even the implied warranty of
85/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
86/// GNU General Public License for more details.
87///
88/// You should have received a copy of the GNU General Public License
89/// along with TUnfold. If not, see <http://www.gnu.org/licenses/>.
90///
91/// \author Stefan Schmitt DESY, 14.10.2008
92
93
94#include <TError.h>
95#include <TMath.h>
96#include <TCanvas.h>
97#include <TRandom3.h>
98#include <TF1.h>
99#include <TStyle.h>
100#include <TVector.h>
101#include <TGraph.h>
102
103#include "TUnfoldDensity.h"
104
105// #define VERBOSE_LCURVE_SCAN
106
107TRandom *rnd=nullptr;
108
110
111TVirtualFitter *gFitter=nullptr;
112
113void chisquare_corr(Int_t &npar, Double_t * /*gin */, Double_t &f, Double_t *u, Int_t /*flag */) {
114 // Minimization function for H1s using a Chisquare method
115 // only one-dimensional histograms are supported
116 // Correlated errors are taken from an external inverse covariance matrix
117 // stored in a 2-dimensional histogram
118 Double_t x;
119 TH1 *hfit = (TH1*)gFitter->GetObjectFit();
120 TF1 *f1 = (TF1*)gFitter->GetUserFunc();
121
122 f1->InitArgs(&x,u);
123 npar = f1->GetNpar();
124 f = 0;
125
126 Int_t npfit = 0;
127 Int_t nPoints=hfit->GetNbinsX();
129 for (Int_t i=0;i<nPoints;i++) {
130 x = hfit->GetBinCenter(i+1);
133 if (TF1::RejectedPoint()) df[i]=0.0;
134 else npfit++;
135 }
136 for (Int_t i=0;i<nPoints;i++) {
137 for (Int_t j=0;j<nPoints;j++) {
138 f += df[i]*df[j]*gHistInvEMatrix->GetBinContent(i+1,j+1);
139 }
140 }
141 delete[] df;
143}
144
146 Double_t dm=x[0]-par[1];
147 return par[0]/(dm*dm+par[2]*par[2]);
148}
149
150
151// generate an event
152// output:
153// negative mass: background event
154// positive mass: signal event
155Double_t GenerateEvent(Double_t bgr, // relative fraction of background
156 Double_t mass, // peak position
157 Double_t gamma) // peak width
158{
159 Double_t t;
160 if(rnd->Rndm()>bgr) {
161 // generate signal event
162 // with positive mass
163 do {
164 do {
165 t=rnd->Rndm();
166 } while(t>=1.0);
167 t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
168 } while(t<=0.0);
169 return t;
170 } else {
171 // generate background event
172 // generate events following a power-law distribution
173 // f(E) = K * TMath::power((E0+E),N0)
174 static Double_t const E0=2.4;
175 static Double_t const N0=2.9;
176 do {
177 do {
178 t=rnd->Rndm();
179 } while(t>=1.0);
180 // the mass is returned negative
181 // In our example a convenient way to indicate it is a background event.
182 t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
183 } while(t>=0.0);
184 return t;
185 }
186}
187
188// smear the event to detector level
189// input:
190// mass on generator level (mTrue>0 !)
191// output:
192// mass on detector level
194 // smear by double-gaussian
195 static Double_t frac=0.1;
196 static Double_t wideBias=0.03;
197 static Double_t wideSigma=0.5;
198 static Double_t smallBias=0.0;
199 static Double_t smallSigma=0.1;
200 if(rnd->Rndm()>frac) {
201 return rnd->Gaus(mTrue+smallBias,smallSigma);
202 } else {
203 return rnd->Gaus(mTrue+wideBias,wideSigma);
204 }
205}
206
207int testUnfold1()
208{
209 // switch on histogram errors
211
212 // show fit result
213 gStyle->SetOptFit(1111);
214
215 // random generator
216 rnd=new TRandom3();
217
218 // data and MC luminosity, cross-section
219 Double_t const luminosityData=100000;
220 Double_t const luminosityMC=1000000;
221 Double_t const crossSection=1.0;
222
223 Int_t const nDet=250;
224 Int_t const nGen=100;
225 Double_t const xminDet=0.0;
226 Double_t const xmaxDet=10.0;
227 Double_t const xminGen=0.0;
228 Double_t const xmaxGen=10.0;
229
230 //============================================
231 // generate MC distribution
232 //
233 TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
234 TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
235 TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",
238 for(Int_t i=0;i<neventMC;i++) {
239 Double_t mGen=GenerateEvent(0.3, // relative fraction of background
240 4.0, // peak position in MC
241 0.2); // peak width in MC
243 // the generated mass is negative for background
244 // and positive for signal
245 // so it will be filled in the underflow bin
246 // this is very convenient for the unfolding:
247 // the unfolded result will contain the number of background
248 // events in the underflow bin
249
250 // generated MC distribution (for comparison only)
252 // reconstructed MC distribution (for comparison only)
254
255 // matrix describing how the generator input migrates to the
256 // reconstructed level. Unfolding input.
257 // NOTE on underflow/overflow bins:
258 // (1) the detector level under/overflow bins are used for
259 // normalisation ("efficiency" correction)
260 // in our toy example, these bins are populated from tails
261 // of the initial MC distribution.
262 // (2) the generator level underflow/overflow bins are
263 // unfolded. In this example:
264 // underflow bin: background events reconstructed in the detector
265 // overflow bin: signal events generated at masses > xmaxDet
266 // for the unfolded result these bins will be filled
267 // -> the background normalisation will be contained in the underflow bin
269 }
270
271 //============================================
272 // generate alternative MC
273 // this will be used to derive a systematic error due to MC
274 // parameter uncertainties
275 TH2D *histMdetGenSysMC=new TH2D("MdetgenSysMC",";mass(det);mass(gen)",
278 for(Int_t i=0;i<neventMC;i++) {
279 Double_t mGen=GenerateEvent
280 (0.5, // relative fraction of background
281 3.6, // peak position in MC with systematic shift
282 0.15); // peak width in MC
285 }
286
287 //============================================
288 // generate data distribution
289 //
290 TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
291 TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
293 for(Int_t i=0;i<neventData;i++) {
294 Double_t mGen=GenerateEvent(0.4, // relative fraction of background
295 3.8, // peak position in data
296 0.15); // peak width in data
298 // generated data mass for comparison plots
299 // for real data, we do not have this histogram
300 histMgenData->Fill(mGen);
301
302 // reconstructed mass, unfolding input
303 histMdetData->Fill(mDet);
304 }
305
306 //=========================================================================
307 // divide by bin width to get density distributions
308 TH1D *histDensityGenData=new TH1D("DensityGenData",";mass(gen)",
310 TH1D *histDensityGenMC=new TH1D("DensityGenMC",";mass(gen)",
312 for(Int_t i=1;i<=nGen;i++) {
313 histDensityGenData->SetBinContent(i,histMgenData->GetBinContent(i)/
314 histMgenData->GetBinWidth(i));
315 histDensityGenMC->SetBinContent(i,histMgenMC->GetBinContent(i)/
316 histMgenMC->GetBinWidth(i));
317 }
318
319 //=========================================================================
320 // set up the unfolding
321 // define migration matrix
323
324 // define input and bias scheme
325 // do not use the bias, because MC peak may be at the wrong place
326 // watch out for error codes returned by the SetInput method
327 // errors larger or equal 10000 are fatal:
328 // the data points specified as input are not sufficient to constrain the
329 // unfolding process
330 if(unfold.SetInput(histMdetData)>=10000) {
331 std::cout<<"Unfolding result may be wrong\n";
332 }
333
334 //========================================================================
335 // the unfolding is done here
336 //
337 // scan L curve and find best point
338 Int_t nScan=30;
339 // use automatic L-curve scan: start with taumin=taumax=0.0
340 Double_t tauMin=0.0;
341 Double_t tauMax=0.0;
342 Int_t iBest;
344 TGraph *lCurve;
345
346 // if required, report Info messages (for debugging the L-curve scan)
347#ifdef VERBOSE_LCURVE_SCAN
350#endif
351 // this method scans the parameter tau and finds the kink in the L curve
352 // finally, the unfolding is done for the best choice of tau
354
355 // if required, switch to previous log-level
356#ifdef VERBOSE_LCURVE_SCAN
358#endif
359
360 //==========================================================================
361 // define a correlated systematic error
362 // for example, assume there is a 10% correlated error for all reconstructed
363 // masses larger than 7
366 TH2D *histMdetGenSys1=new TH2D("Mdetgensys1",";mass(det);mass(gen)",
368 for(Int_t i=0;i<=nDet+1;i++) {
369 if(histMdetData->GetBinCenter(i)>=SYS_ERROR1_MSTART) {
370 for(Int_t j=0;j<=nGen+1;j++) {
371 histMdetGenSys1->SetBinContent(i,j,SYS_ERROR1_SIZE);
372 }
373 }
374 }
375 unfold.AddSysError(histMdetGenSysMC,"SYSERROR_MC",TUnfold::kHistMapOutputVert,
377 unfold.AddSysError(histMdetGenSys1,"SYSERROR1",TUnfold::kHistMapOutputVert,
379
380 //==========================================================================
381 // print some results
382 //
383 std::cout<<"tau="<<unfold.GetTau()<<"\n";
384 std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
385 <<" / "<<unfold.GetNdf()<<"\n";
386 std::cout<<"chi**2(sys)="<<unfold.GetChi2Sys()<<"\n";
387
388
389 //==========================================================================
390 // create graphs with one point to visualize the best choice of tau
391 //
392 Double_t t[1],x[1],y[1];
393 logTauX->GetKnot(iBest,t[0],x[0]);
394 logTauY->GetKnot(iBest,t[0],y[0]);
395 TGraph *bestLcurve=new TGraph(1,x,y);
397
398 //==========================================================================
399 // retrieve results into histograms
400
401 // get unfolded distribution
402 TH1 *histMunfold=unfold.GetOutput("Unfolded");
403
404 // get unfolding result, folded back
405 TH1 *histMdetFold=unfold.GetFoldedOutput("FoldedBack");
406
407 // get error matrix (input distribution [stat] errors only)
408 // TH2D *histEmatData=unfold.GetEmatrix("EmatData");
409
410 // get total error matrix:
411 // migration matrix uncorrelated and correlated systematic errors
412 // added in quadrature to the data statistical errors
413 TH2 *histEmatTotal=unfold.GetEmatrixTotal("EmatTotal");
414
415 // create data histogram with the total errors
417 new TH1D("TotalError",";mass(gen)",nGen,xminGen,xmaxGen);
418 for(Int_t bin=1;bin<=nGen;bin++) {
419 histTotalError->SetBinContent(bin,histMunfold->GetBinContent(bin));
420 histTotalError->SetBinError
421 (bin,TMath::Sqrt(histEmatTotal->GetBinContent(bin,bin)));
422 }
423
424 // get global correlation coefficients
425 // for this calculation one has to specify whether the
426 // underflow/overflow bins are included or not
427 // default: include all bins
428 // here: exclude underflow and overflow bins
429 TH1 *histRhoi=unfold.GetRhoItotal("rho_I",
430 nullptr, // use default title
431 nullptr, // all distributions
432 "*[UO]", // discard underflow and overflow bins on all axes
433 kTRUE, // use original binning
434 &gHistInvEMatrix // store inverse of error matrix
435 );
436
437 //======================================================================
438 // fit Breit-Wigner shape to unfolded data, using the full error matrix
439 // here we use a "user" chi**2 function to take into account
440 // the full covariance matrix
441
443 gFitter->SetFCN(chisquare_corr);
444
445 TF1 *bw=new TF1("bw",bw_func,xminGen,xmaxGen,3);
446 bw->SetParameter(0,1000.);
447 bw->SetParameter(1,3.8);
448 bw->SetParameter(2,0.2);
449
450 // for (wrong!) fitting without correlations, drop the option "U"
451 // here.
452 histMunfold->Fit(bw,"UE");
453
454 //=====================================================================
455 // plot some histograms
457 output.Divide(3,2);
458
459 // Show the matrix which connects input and output
460 // There are overflow bins at the bottom, not shown in the plot
461 // These contain the background shape.
462 // The overflow bins to the left and right contain
463 // events which are not reconstructed. These are necessary for proper MC
464 // normalisation
465 output.cd(1);
466 histMdetGenMC->Draw("BOX");
467
468 // draw generator-level distribution:
469 // data (red) [for real data this is not available]
470 // MC input (black) [with completely wrong peak position and shape]
471 // unfolded data (blue)
472 output.cd(2);
473 histTotalError->SetLineColor(kBlue);
474 histTotalError->Draw("E");
475 histMunfold->SetLineColor(kGreen);
476 histMunfold->Draw("SAME E1");
477 histDensityGenData->SetLineColor(kRed);
478 histDensityGenData->Draw("SAME");
479 histDensityGenMC->Draw("SAME HIST");
480
481 // show detector level distributions
482 // data (red)
483 // MC (black) [with completely wrong peak position and shape]
484 // unfolded data (blue)
485 output.cd(3);
486 histMdetFold->SetLineColor(kBlue);
487 histMdetFold->Draw();
488 histMdetMC->Draw("SAME HIST");
489
490 TH1 *histInput=unfold.GetInput("Minput",";mass(det)");
491
492 histInput->SetLineColor(kRed);
493 histInput->Draw("SAME");
494
495 // show correlation coefficients
496 output.cd(4);
497 histRhoi->Draw();
498
499 // show tau as a function of chi**2
500 output.cd(5);
501 logTauX->Draw();
502 bestLogTauLogChi2->SetMarkerColor(kRed);
503 bestLogTauLogChi2->Draw("*");
504
505 // show the L curve
506 output.cd(6);
507 lCurve->Draw("AL");
508 bestLcurve->SetMarkerColor(kRed);
509 bestLcurve->Draw("*");
510
511 output.SaveAs("testUnfold1.ps");
512
513 return 0;
514}
#define f(i)
Definition RSha256.hxx:104
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
@ kRed
Definition Rtypes.h:67
@ kGreen
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.
constexpr Int_t kInfo
Definition TError.h:45
Int_t gErrorIgnoreLevel
errors with level below this value will be ignored. Default is kUnset.
Definition TError.cxx:33
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:182
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition TF1.cxx:3699
virtual Int_t GetNpar() const
Definition TF1.h:461
virtual void SetNumberFitPoints(Int_t npfits)
Definition TF1.h:608
virtual void InitArgs(const Double_t *x, const Double_t *params)
Initialize parameters addresses.
Definition TF1.cxx:2507
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
Evaluate function with given coordinates and parameters.
Definition TF1.cxx:1475
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition TF1.cxx:3708
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
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
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
Service class for 2-D histogram classes.
Definition TH2.h:30
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
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.
@ kSysErrModeRelative
matrix gives the relative shifts
Definition TUnfoldSys.h:112
@ kSysErrModeMatrix
matrix is an alternative to the default matrix, the errors are the difference to the original matrix
Definition TUnfoldSys.h:108
@ kHistMapOutputVert
truth level on y-axis of the response matrix
Definition TUnfold.h:149
Abstract Base Class for Fitting.
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
TF1 * f1
Definition legend1.C:11
double gamma(double x)
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
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()