ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
testUnfold1.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_unfold
3 /// Test program for the classes TUnfold and related
4 ///
5 /// 1. Generate Monte Carlo and Data events
6 /// The events consist of
7 /// - signal
8 /// - background
9 ///
10 /// The signal is a resonance. It is generated with a Breit-Wigner,
11 /// smeared by a Gaussian
12 ///
13 /// 2. Unfold the data. The result is:
14 /// The background level
15 /// The shape of the resonance, corrected for detector effects
16 ///
17 /// Systematic errors from the MC shape variation are included
18 /// and propagated to the result
19 ///
20 /// 3. fit the unfolded distribution, including the correlation matrix
21 ///
22 /// 4. save six plots to a file testUnfold1.ps
23 /// - 1 2 3
24 /// - 4 5 6
25 /// 1. 2d-plot of the matrix decsribing the migrations
26 /// 2. generator-level distributions
27 /// - blue: unfolded data, total errors
28 /// - green: unfolded data, statistical errors
29 /// - red: generated data
30 /// - black: fit to green data points
31 /// 3. detector level distributions
32 /// - blue: unfoldede data, folded back through the matrix
33 /// - black: Monte Carlo (with wrong peal position)
34 /// - blue: data
35 /// 4. global correlation coefficients
36 /// 5. \f$ \chi^2 \f$ as a function of log(tau)
37 /// - the star indicates the final choice of tau
38 /// 6. the L curve
39 ///
40 /// Version 17.0, updated for using the classes TUnfoldDensity, TUnfoldBinning
41 ///
42 /// History:
43 /// - Version 16.1, parallel to changes in TUnfold
44 /// - Version 16.0, parallel to changes in TUnfold
45 /// - Version 15, with automated L-curve scan
46 /// - Version 14, with changes in TUnfoldSys.cxx
47 /// - Version 13, include test of systematic errors
48 /// - Version 12, catch error when defining the input
49 /// - Version 11, print chi**2 and number of degrees of freedom
50 /// - Version 10, with bug-fix in TUnfold.cxx
51 /// - Version 9, with bug-fix in TUnfold.cxx and TUnfold.h
52 /// - Version 8, with bug-fix in TUnfold.cxx and TUnfold.h
53 /// - Version 7, with bug-fix in TUnfold.cxx and TUnfold.h
54 /// - Version 6a, fix problem with dynamic array allocation under windows
55 /// - Version 6, bug-fixes in TUnfold.C
56 /// - Version 5, replace main() by testUnfold1()
57 /// - Version 4, with bug-fix in TUnfold.C
58 /// - Version 3, with bug-fix in TUnfold.C
59 /// - Version 2, with changed ScanLcurve() arguments
60 /// - Version 1, remove L curve analysis, use ScanLcurve() method instead
61 /// - Version 0, L curve analysis included here
62 ///
63 /// \macro_image
64 /// \macro_output
65 /// \macro_code
66 ///
67 /// \author Stefan Schmitt, DESY
68 
69 #include <TError.h>
70 #include <TMath.h>
71 #include <TCanvas.h>
72 #include <TRandom3.h>
73 #include <TFitter.h>
74 #include <TF1.h>
75 #include <TStyle.h>
76 #include <TVector.h>
77 #include <TGraph.h>
78 
79 #include "TUnfoldDensity.h"
80 
81 // #define VERBOSE_LCURVE_SCAN
82 
83 using namespace std;
84 
85 TRandom *rnd=0;
86 
87 TH2 *gHistInvEMatrix;
88 
89 TVirtualFitter *gFitter=0;
90 
91 void chisquare_corr(Int_t &npar, Double_t * /*gin */, Double_t &f, Double_t *u, Int_t /*flag */) {
92  // Minimization function for H1s using a Chisquare method
93  // only one-dimensional histograms are supported
94  // Corelated errors are taken from an external inverse covariance matrix
95  // stored in a 2-dimensional histogram
96  Double_t x;
97  TH1 *hfit = (TH1*)gFitter->GetObjectFit();
98  TF1 *f1 = (TF1*)gFitter->GetUserFunc();
99 
100  f1->InitArgs(&x,u);
101  npar = f1->GetNpar();
102  f = 0;
103 
104  Int_t npfit = 0;
105  Int_t nPoints=hfit->GetNbinsX();
106  Double_t *df=new Double_t[nPoints];
107  for (Int_t i=0;i<nPoints;i++) {
108  x = hfit->GetBinCenter(i+1);
110  df[i] = f1->EvalPar(&x,u)-hfit->GetBinContent(i+1);
111  if (TF1::RejectedPoint()) df[i]=0.0;
112  else npfit++;
113  }
114  for (Int_t i=0;i<nPoints;i++) {
115  for (Int_t j=0;j<nPoints;j++) {
116  f += df[i]*df[j]*gHistInvEMatrix->GetBinContent(i+1,j+1);
117  }
118  }
119  delete[] df;
120  f1->SetNumberFitPoints(npfit);
121 }
122 
123 Double_t bw_func(Double_t *x,Double_t *par) {
124  Double_t dm=x[0]-par[1];
125  return par[0]/(dm*dm+par[2]*par[2]);
126 }
127 
128 
129 // generate an event
130 // output:
131 // negative mass: background event
132 // positive mass: signal event
133 Double_t GenerateEvent(Double_t bgr, // relative fraction of background
134  Double_t mass, // peak position
135  Double_t gamma) // peak width
136 {
137  Double_t t;
138  if(rnd->Rndm()>bgr) {
139  // generate signal event
140  // with positive mass
141  do {
142  do {
143  t=rnd->Rndm();
144  } while(t>=1.0);
145  t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
146  } while(t<=0.0);
147  return t;
148  } else {
149  // generate background event
150  // generate events following a power-law distribution
151  // f(E) = K * TMath::power((E0+E),N0)
152  static Double_t const E0=2.4;
153  static Double_t const N0=2.9;
154  do {
155  do {
156  t=rnd->Rndm();
157  } while(t>=1.0);
158  // the mass is returned negative
159  // In our example a convenient way to indicate it is a background event.
160  t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
161  } while(t>=0.0);
162  return t;
163  }
164 }
165 
166 // smear the event to detector level
167 // input:
168 // mass on generator level (mTrue>0 !)
169 // output:
170 // mass on detector level
171 Double_t DetectorEvent(Double_t mTrue) {
172  // smear by double-gaussian
173  static Double_t frac=0.1;
174  static Double_t wideBias=0.03;
175  static Double_t wideSigma=0.5;
176  static Double_t smallBias=0.0;
177  static Double_t smallSigma=0.1;
178  if(rnd->Rndm()>frac) {
179  return rnd->Gaus(mTrue+smallBias,smallSigma);
180  } else {
181  return rnd->Gaus(mTrue+wideBias,wideSigma);
182  }
183 }
184 
185 int testUnfold1()
186 {
187  // switch on histogram errors
189 
190  // show fit result
191  gStyle->SetOptFit(1111);
192 
193  // random generator
194  rnd=new TRandom3();
195 
196  // data and MC luminosity, cross-section
197  Double_t const luminosityData=100000;
198  Double_t const luminosityMC=1000000;
199  Double_t const crossSection=1.0;
200 
201  Int_t const nDet=250;
202  Int_t const nGen=100;
203  Double_t const xminDet=0.0;
204  Double_t const xmaxDet=10.0;
205  Double_t const xminGen=0.0;
206  Double_t const xmaxGen=10.0;
207 
208  //============================================
209  // generate MC distribution
210  //
211  TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
212  TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
213  TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",
214  nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
215  Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
216  for(Int_t i=0;i<neventMC;i++) {
217  Double_t mGen=GenerateEvent(0.3, // relative fraction of background
218  4.0, // peak position in MC
219  0.2); // peak width in MC
220  Double_t mDet=DetectorEvent(TMath::Abs(mGen));
221  // the generated mass is negative for background
222  // and positive for signal
223  // so it will be filled in the underflow bin
224  // this is very convenient for the unfolding:
225  // the unfolded result will contain the number of background
226  // events in the underflow bin
227 
228  // generated MC distribution (for comparison only)
229  histMgenMC->Fill(mGen,luminosityData/luminosityMC);
230  // reconstructed MC distribution (for comparison only)
231  histMdetMC->Fill(mDet,luminosityData/luminosityMC);
232 
233  // matrix describing how the generator input migrates to the
234  // reconstructed level. Unfolding input.
235  // NOTE on underflow/overflow bins:
236  // (1) the detector level under/overflow bins are used for
237  // normalisation ("efficiency" correction)
238  // in our toy example, these bins are populated from tails
239  // of the initial MC distribution.
240  // (2) the generator level underflow/overflow bins are
241  // unfolded. In this example:
242  // underflow bin: background events reconstructed in the detector
243  // overflow bin: signal events generated at masses > xmaxDet
244  // for the unfolded result these bins will be filled
245  // -> the background normalisation will be contained in the underflow bin
246  histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
247  }
248 
249  //============================================
250  // generate alternative MC
251  // this will be used to derive a systematic error due to MC
252  // parameter uncertainties
253  TH2D *histMdetGenSysMC=new TH2D("MdetgenSysMC",";mass(det);mass(gen)",
254  nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
255  neventMC=rnd->Poisson(luminosityMC*crossSection);
256  for(Int_t i=0;i<neventMC;i++) {
257  Double_t mGen=GenerateEvent
258  (0.5, // relative fraction of background
259  3.6, // peak position in MC with systematic shift
260  0.15); // peak width in MC
261  Double_t mDet=DetectorEvent(TMath::Abs(mGen));
262  histMdetGenSysMC->Fill(mDet,mGen,luminosityData/luminosityMC);
263  }
264 
265  //============================================
266  // generate data distribution
267  //
268  TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
269  TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
270  Int_t neventData=rnd->Poisson(luminosityData*crossSection);
271  for(Int_t i=0;i<neventData;i++) {
272  Double_t mGen=GenerateEvent(0.4, // relative fraction of background
273  3.8, // peak position in data
274  0.15); // peak width in data
275  Double_t mDet=DetectorEvent(TMath::Abs(mGen));
276  // generated data mass for comparison plots
277  // for real data, we do not have this histogram
278  histMgenData->Fill(mGen);
279 
280  // reconstructed mass, unfolding input
281  histMdetData->Fill(mDet);
282  }
283 
284  //=========================================================================
285  // divide by bin withd to get density distributions
286  TH1D *histDensityGenData=new TH1D("DensityGenData",";mass(gen)",
287  nGen,xminGen,xmaxGen);
288  TH1D *histDensityGenMC=new TH1D("DensityGenMC",";mass(gen)",
289  nGen,xminGen,xmaxGen);
290  for(Int_t i=1;i<=nGen;i++) {
291  histDensityGenData->SetBinContent(i,histMgenData->GetBinContent(i)/
292  histMgenData->GetBinWidth(i));
293  histDensityGenMC->SetBinContent(i,histMgenMC->GetBinContent(i)/
294  histMgenMC->GetBinWidth(i));
295  }
296 
297  //=========================================================================
298  // set up the unfolding
299  // define migration matrix
300  TUnfoldDensity unfold(histMdetGenMC,TUnfold::kHistMapOutputVert);
301 
302  // define input and bias scame
303  // do not use the bias, because MC peak may be at the wrong place
304  // watch out for error codes returned by the SetInput method
305  // errors larger or equal 10000 are fatal:
306  // the data points specified as input are not sufficient to constrain the
307  // unfolding process
308  if(unfold.SetInput(histMdetData)>=10000) {
309  std::cout<<"Unfolding result may be wrong\n";
310  }
311 
312  //========================================================================
313  // the unfolding is done here
314  //
315  // scan L curve and find best point
316  Int_t nScan=30;
317  // use automatic L-curve scan: start with taumin=taumax=0.0
318  Double_t tauMin=0.0;
319  Double_t tauMax=0.0;
320  Int_t iBest;
321  TSpline *logTauX,*logTauY;
322  TGraph *lCurve;
323 
324  // if required, report Info messages (for debugging the L-curve scan)
325 #ifdef VERBOSE_LCURVE_SCAN
326  Int_t oldinfo=gErrorIgnoreLevel;
328 #endif
329  // this method scans the parameter tau and finds the kink in the L curve
330  // finally, the unfolding is done for the best choice of tau
331  iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
332 
333  // if required, switch to previous log-level
334 #ifdef VERBOSE_LCURVE_SCAN
335  gErrorIgnoreLevel=oldinfo;
336 #endif
337 
338  //==========================================================================
339  // define a correlated systematic error
340  // for example, assume there is a 10% correlated error for all reconstructed
341  // masses larger than 7
342  Double_t SYS_ERROR1_MSTART=6;
343  Double_t SYS_ERROR1_SIZE=0.1;
344  TH2D *histMdetGenSys1=new TH2D("Mdetgensys1",";mass(det);mass(gen)",
345  nDet,xminDet,xmaxDet,nGen,xminGen,xmaxGen);
346  for(Int_t i=0;i<=nDet+1;i++) {
347  if(histMdetData->GetBinCenter(i)>=SYS_ERROR1_MSTART) {
348  for(Int_t j=0;j<=nGen+1;j++) {
349  histMdetGenSys1->SetBinContent(i,j,SYS_ERROR1_SIZE);
350  }
351  }
352  }
353  unfold.AddSysError(histMdetGenSysMC,"SYSERROR_MC",TUnfold::kHistMapOutputVert,
355  unfold.AddSysError(histMdetGenSys1,"SYSERROR1",TUnfold::kHistMapOutputVert,
357 
358  //==========================================================================
359  // print some results
360  //
361  std::cout<<"tau="<<unfold.GetTau()<<"\n";
362  std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
363  <<" / "<<unfold.GetNdf()<<"\n";
364  std::cout<<"chi**2(sys)="<<unfold.GetChi2Sys()<<"\n";
365 
366 
367  //==========================================================================
368  // create graphs with one point to visualize the best choice of tau
369  //
370  Double_t t[1],x[1],y[1];
371  logTauX->GetKnot(iBest,t[0],x[0]);
372  logTauY->GetKnot(iBest,t[0],y[0]);
373  TGraph *bestLcurve=new TGraph(1,x,y);
374  TGraph *bestLogTauLogChi2=new TGraph(1,t,x);
375 
376  //==========================================================================
377  // retreive results into histograms
378 
379  // get unfolded distribution
380  TH1 *histMunfold=unfold.GetOutput("Unfolded");
381 
382  // get unfolding result, folded back
383  TH1 *histMdetFold=unfold.GetFoldedOutput("FoldedBack");
384 
385  // get error matrix (input distribution [stat] errors only)
386  // TH2D *histEmatData=unfold.GetEmatrix("EmatData");
387 
388  // get total error matrix:
389  // migration matrix uncorrelated and correlated systematic errors
390  // added in quadrature to the data statistical errors
391  TH2 *histEmatTotal=unfold.GetEmatrixTotal("EmatTotal");
392 
393  // create data histogram with the total errors
394  TH1D *histTotalError=
395  new TH1D("TotalError",";mass(gen)",nGen,xminGen,xmaxGen);
396  for(Int_t bin=1;bin<=nGen;bin++) {
397  histTotalError->SetBinContent(bin,histMunfold->GetBinContent(bin));
398  histTotalError->SetBinError
399  (bin,TMath::Sqrt(histEmatTotal->GetBinContent(bin,bin)));
400  }
401 
402  // get global correlation coefficients
403  // for this calculation one has to specify whether the
404  // underflow/overflow bins are included or not
405  // default: include all bins
406  // here: exclude underflow and overflow bins
407  TH1 *histRhoi=unfold.GetRhoItotal("rho_I",
408  0, // use default title
409  0, // all distributions
410  "*[UO]", // discard underflow and overflow bins on all axes
411  kTRUE, // use original binning
412  &gHistInvEMatrix // store inverse of error matrix
413  );
414 
415  //======================================================================
416  // fit Breit-Wigner shape to unfolded data, using the full error matrix
417  // here we use a "user" chi**2 function to take into account
418  // the full covariance matrix
419 
420  gFitter=TVirtualFitter::Fitter(histMunfold);
421  gFitter->SetFCN(chisquare_corr);
422 
423  TF1 *bw=new TF1("bw",bw_func,xminGen,xmaxGen,3);
424  bw->SetParameter(0,1000.);
425  bw->SetParameter(1,3.8);
426  bw->SetParameter(2,0.2);
427 
428  // for (wrong!) fitting without correlations, drop the option "U"
429  // here.
430  histMunfold->Fit(bw,"UE");
431 
432  //=====================================================================
433  // plot some histograms
434  TCanvas output;
435  output.Divide(3,2);
436 
437  // Show the matrix which connects input and output
438  // There are overflow bins at the bottom, not shown in the plot
439  // These contain the background shape.
440  // The overflow bins to the left and right contain
441  // events which are not reconstructed. These are necessary for proper MC
442  // normalisation
443  output.cd(1);
444  histMdetGenMC->Draw("BOX");
445 
446  // draw generator-level distribution:
447  // data (red) [for real data this is not available]
448  // MC input (black) [with completely wrong peak position and shape]
449  // unfolded data (blue)
450  output.cd(2);
451  histTotalError->SetLineColor(kBlue);
452  histTotalError->Draw("E");
453  histMunfold->SetLineColor(kGreen);
454  histMunfold->Draw("SAME E1");
455  histDensityGenData->SetLineColor(kRed);
456  histDensityGenData->Draw("SAME");
457  histDensityGenMC->Draw("SAME HIST");
458 
459  // show detector level distributions
460  // data (red)
461  // MC (black) [with completely wrong peak position and shape]
462  // unfolded data (blue)
463  output.cd(3);
464  histMdetFold->SetLineColor(kBlue);
465  histMdetFold->Draw();
466  histMdetMC->Draw("SAME HIST");
467 
468  TH1 *histInput=unfold.GetInput("Minput",";mass(det)");
469 
470  histInput->SetLineColor(kRed);
471  histInput->Draw("SAME");
472 
473  // show correlation coefficients
474  output.cd(4);
475  histRhoi->Draw();
476 
477  // show tau as a function of chi**2
478  output.cd(5);
479  logTauX->Draw();
480  bestLogTauLogChi2->SetMarkerColor(kRed);
481  bestLogTauLogChi2->Draw("*");
482 
483  // show the L curve
484  output.cd(6);
485  lCurve->Draw("AL");
486  bestLcurve->SetMarkerColor(kRed);
487  bestLcurve->Draw("*");
488 
489  output.SaveAs("testUnfold1.ps");
490 
491  return 0;
492 }
double par[1]
Definition: unuranDistr.cxx:38
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
Random number generator class based on M.
Definition: TRandom3.h:29
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4629
R__EXTERN Int_t gErrorIgnoreLevel
Definition: TError.h:107
Definition: Rtypes.h:61
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
Base class for spline implementation containing the Draw/Paint methods //.
Definition: TSpline.h:22
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: Rtypes.h:92
Definition: Rtypes.h:61
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
static void SetDefaultSumw2(Bool_t sumw2=kTRUE)
static function.
Definition: TH1.cxx:6253
virtual Double_t GetBinWidth(Int_t bin) const
return bin width for 1D historam Better to use h1.GetXaxis().GetBinWidth(bin)
Definition: TH1.cxx:8492
Double_t x[n]
Definition: legend1.C:17
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:29
virtual void SetBinError(Int_t bin, Double_t error)
see convention for numbering bins in TH1::GetBin
Definition: TH1.cxx:8528
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
virtual Double_t GetBinCenter(Int_t bin) const
return bin center for 1D historam Better to use h1.GetXaxis().GetBinCenter(bin)
Definition: TH1.cxx:8470
TThread * t[5]
Definition: threadsh1.C:13
Service class for 2-Dim histogram classes.
Definition: TH2.h:36
const Int_t kInfo
Definition: TError.h:39
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
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:1204
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:3367
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0
Double_t Pi()
Definition: TMath.h:44
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
Abstract Base Class for Fitting.
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
static Bool_t RejectedPoint()
See TF1::RejectPoint above.
Definition: TF1.cxx:3376
Definition: Rtypes.h:61
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2699
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
static void output(int code)
Definition: gifencode.c:226
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:287
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3607
Double_t Tan(Double_t)
Definition: TMath.h:427
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:297