Logo ROOT   6.07/09
Reference Guide
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 /// Version 17.0, updated for changed methods in TUnfold
24 ///
25 /// History:
26 /// - Version 16.1, parallel to changes in TUnfold
27 /// - Version 16.0, parallel to changes in TUnfold
28 /// - Version 15, with automatic L-curve scan, simplified example
29 /// - Version 14, with changes in TUnfoldSys.cxx
30 /// - Version 13, with changes to TUnfold.C
31 /// - Version 12, with improvements to TUnfold.cxx
32 /// - Version 11, print chi**2 and number of degrees of freedom
33 /// - Version 10, with bug-fix in TUnfold.cxx
34 /// - Version 9, with bug-fix in TUnfold.cxx, TUnfold.h
35 /// - Version 8, with bug-fix in TUnfold.cxx, TUnfold.h
36 /// - Version 7, with bug-fix in TUnfold.cxx, TUnfold.h
37 /// - Version 6a, fix problem with dynamic array allocation under windows
38 /// - Version 6, re-include class MyUnfold in the example
39 /// - Version 5, move class MyUnfold to seperate files
40 /// - Version 4, with bug-fix in TUnfold.C
41 /// - Version 3, with bug-fix in TUnfold.C
42 /// - Version 2, with changed ScanLcurve() arguments
43 /// - Version 1, remove L curve analysis, use ScanLcurve() method instead
44 /// - Version 0, L curve analysis included here
45 ///
46 /// \macro_image
47 /// \macro_output
48 /// \macro_code
49 ///
50 /// \author Stefan Schmitt, DESY
51 
52 #include <TMath.h>
53 #include <TCanvas.h>
54 #include <TRandom3.h>
55 #include <TFitter.h>
56 #include <TF1.h>
57 #include <TStyle.h>
58 #include <TVector.h>
59 #include <TGraph.h>
60 #include "TUnfold.h"
61 
62 using namespace std;
63 
64 TRandom *rnd=0;
65 
66 // generate an event
67 // output:
68 // negative mass: background event
69 // positive mass: signal event
70 Double_t GenerateEvent(Double_t bgr, // relative fraction of background
71  Double_t mass, // peak position
72  Double_t gamma /* peak width*/ )
73 {
74  Double_t t;
75  if(rnd->Rndm()>bgr) {
76  // generate signal event
77  // with positive mass
78  do {
79  do {
80  t=rnd->Rndm();
81  } while(t>=1.0);
82  t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
83  } while(t<=0.0);
84  return t;
85  } else {
86  // generate background event
87  // generate events following a power-law distribution
88  // f(E) = K * TMath::power((E0+E),N0)
89  static Double_t const E0=2.4;
90  static Double_t const N0=2.9;
91  do {
92  do {
93  t=rnd->Rndm();
94  } while(t>=1.0);
95  // the mass is returned negative
96  // In our example a convenient way to indicate it is a background event.
97  t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
98  } while(t>=0.0);
99  return t;
100  }
101 }
102 
103 // smear the event to detector level
104 // input:
105 // mass on generator level (mTrue>0 !)
106 // output:
107 // mass on detector level
108 Double_t DetectorEvent(Double_t mTrue) {
109  // smear by double-gaussian
110  static Double_t frac=0.1;
111  static Double_t wideBias=0.03;
112  static Double_t wideSigma=0.5;
113  static Double_t smallBias=0.0;
114  static Double_t smallSigma=0.1;
115  if(rnd->Rndm()>frac) {
116  return rnd->Gaus(mTrue+smallBias,smallSigma);
117  } else {
118  return rnd->Gaus(mTrue+wideBias,wideSigma);
119  }
120 }
121 
122 int testUnfold2()
123 {
124  // switch on histogram errors
126 
127  // random generator
128  rnd=new TRandom3();
129 
130  // data and MC luminosity, cross-section
131  Double_t const luminosityData=100000;
132  Double_t const luminosityMC=1000000;
133  Double_t const crossSection=1.0;
134 
135  Int_t const nDet=250;
136  Int_t const nGen=100;
137  Double_t const xminDet=0.0;
138  Double_t const xmaxDet=10.0;
139  Double_t const xminGen=0.0;
140  Double_t const xmaxGen=10.0;
141 
142  //============================================
143  // generate MC distribution
144  //
145  TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
146  TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
147  TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",nDet,xminDet,xmaxDet,
148  nGen,xminGen,xmaxGen);
149  Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
150  for(Int_t i=0;i<neventMC;i++) {
151  Double_t mGen=GenerateEvent(0.3, // relative fraction of background
152  4.0, // peak position in MC
153  0.2); // peak width in MC
154  Double_t mDet=DetectorEvent(TMath::Abs(mGen));
155  // the generated mass is negative for background
156  // and positive for signal
157  // so it will be filled in the underflow bin
158  // this is very convenient for the unfolding:
159  // the unfolded result will contain the number of background
160  // events in the underflow bin
161 
162  // generated MC distribution (for comparison only)
163  histMgenMC->Fill(mGen,luminosityData/luminosityMC);
164  // reconstructed MC distribution (for comparison only)
165  histMdetMC->Fill(mDet,luminosityData/luminosityMC);
166 
167  // matrix describing how the generator input migrates to the
168  // reconstructed level. Unfolding input.
169  // NOTE on underflow/overflow bins:
170  // (1) the detector level under/overflow bins are used for
171  // normalisation ("efficiency" correction)
172  // in our toy example, these bins are populated from tails
173  // of the initial MC distribution.
174  // (2) the generator level underflow/overflow bins are
175  // unfolded. In this example:
176  // underflow bin: background events reconstructed in the detector
177  // overflow bin: signal events generated at masses > xmaxDet
178  // for the unfolded result these bins will be filled
179  // -> the background normalisation will be contained in the underflow bin
180  histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
181  }
182 
183  //============================================
184  // generate data distribution
185  //
186  TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
187  TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
188  Int_t neventData=rnd->Poisson(luminosityData*crossSection);
189  for(Int_t i=0;i<neventData;i++) {
190  Double_t mGen=GenerateEvent(0.4, // relative fraction of background
191  3.8, // peak position
192  0.15); // peak width
193  Double_t mDet=DetectorEvent(TMath::Abs(mGen));
194  // generated data mass for comparison plots
195  // for real data, we do not have this histogram
196  histMgenData->Fill(mGen);
197 
198  // reconstructed mass, unfolding input
199  histMdetData->Fill(mDet);
200  }
201 
202  //=========================================================================
203  // set up the unfolding
204  TUnfold unfold(histMdetGenMC,TUnfold::kHistMapOutputVert,
206  // regularisation
207  //----------------
208  // the regularisation is done on the curvature (2nd derivative) of
209  // the output distribution
210  //
211  // One has to exclude the bins near the peak of the Breit-Wigner,
212  // because there the curvature is high
213  // (and the regularisation eventually could enforce a small
214  // curvature, thus biasing result)
215  //
216  // in real life, the parameters below would have to be optimized,
217  // depending on the data peak position and width
218  // Or maybe one finds a different regularisation scheme... this is
219  // just an example...
220  Double_t estimatedPeakPosition=3.8;
221  Int_t nPeek=3;
223  // calculate bin number correspoinding to estimated peak position
224  Int_t iPeek=(Int_t)(nGen*(estimatedPeakPosition-xminGen)/(xmaxGen-xminGen)
225  // offset 1.5
226  // accounts for start bin 1
227  // and rounding errors +0.5
228  +1.5);
229  // regularize output bins 1..iPeek-nPeek
230  unfold.RegularizeBins(1,1,iPeek-nPeek,regMode);
231  // regularize output bins iPeek+nPeek..nGen
232  unfold.RegularizeBins(iPeek+nPeek,1,nGen-(iPeek+nPeek),regMode);
233 
234  // unfolding
235  //-----------
236 
237  // set input distribution and bias scale (=0)
238  if(unfold.SetInput(histMdetData,0.0)>=10000) {
239  std::cout<<"Unfolding result may be wrong\n";
240  }
241 
242  // do the unfolding here
243  Double_t tauMin=0.0;
244  Double_t tauMax=0.0;
245  Int_t nScan=30;
246  Int_t iBest;
247  TSpline *logTauX,*logTauY;
248  TGraph *lCurve;
249  // this method scans the parameter tau and finds the kink in the L curve
250  // finally, the unfolding is done for the "best" choice of tau
251  iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
252  std::cout<<"tau="<<unfold.GetTau()<<"\n";
253  std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
254  <<" / "<<unfold.GetNdf()<<"\n";
255 
256  // save point corresponding to the kink in the L curve as TGraph
257  Double_t t[1],x[1],y[1];
258  logTauX->GetKnot(iBest,t[0],x[0]);
259  logTauY->GetKnot(iBest,t[0],y[0]);
260  TGraph *bestLcurve=new TGraph(1,x,y);
261  TGraph *bestLogTauX=new TGraph(1,t,x);
262 
263  //============================================================
264  // extract unfolding results into histograms
265 
266  // set up a bin map, excluding underflow and overflow bins
267  // the binMap relates the the output of the unfolding to the final
268  // histogram bins
269  Int_t *binMap=new Int_t[nGen+2];
270  for(Int_t i=1;i<=nGen;i++) binMap[i]=i;
271  binMap[0]=-1;
272  binMap[nGen+1]=-1;
273 
274  TH1D *histMunfold=new TH1D("Unfolded",";mass(gen)",nGen,xminGen,xmaxGen);
275  unfold.GetOutput(histMunfold,binMap);
276  TH1D *histMdetFold=new TH1D("FoldedBack","mass(det)",nDet,xminDet,xmaxDet);
277  unfold.GetFoldedOutput(histMdetFold);
278 
279  // store global correlation coefficients
280  TH1D *histRhoi=new TH1D("rho_I","mass",nGen,xminGen,xmaxGen);
281  unfold.GetRhoI(histRhoi,binMap);
282 
283  delete[] binMap;
284  binMap=0;
285 
286  //=====================================================================
287  // plot some histograms
288  TCanvas *output = new TCanvas();
289 
290  // produce some plots
291  output->Divide(3,2);
292 
293  // Show the matrix which connects input and output
294  // There are overflow bins at the bottom, not shown in the plot
295  // These contain the background shape.
296  // The overflow bins to the left and right contain
297  // events which are not reconstructed. These are necessary for proper MC
298  // normalisation
299  output->cd(1);
300  histMdetGenMC->Draw("BOX");
301 
302  // draw generator-level distribution:
303  // data (red) [for real data this is not available]
304  // MC input (black) [with completely wrong peak position and shape]
305  // unfolded data (blue)
306  output->cd(2);
307  histMunfold->SetLineColor(kBlue);
308  histMunfold->Draw();
309  histMgenData->SetLineColor(kRed);
310  histMgenData->Draw("SAME");
311  histMgenMC->Draw("SAME HIST");
312 
313  // show detector level distributions
314  // data (red)
315  // MC (black)
316  // unfolded data (blue)
317  output->cd(3);
318  histMdetFold->SetLineColor(kBlue);
319  histMdetFold->Draw();
320  histMdetData->SetLineColor(kRed);
321  histMdetData->Draw("SAME");
322  histMdetMC->Draw("SAME HIST");
323 
324  // show correlation coefficients
325  // all bins outside the peak are found to be highly correlated
326  // But they are compatible with zero anyway
327  // If the peak shape is fitted,
328  // these correlations have to be taken into account, see example
329  output->cd(4);
330  histRhoi->Draw();
331 
332  // show rhoi_max(tau) distribution
333  output->cd(5);
334  logTauX->Draw();
335  bestLogTauX->SetMarkerColor(kRed);
336  bestLogTauX->Draw("*");
337 
338  output->cd(6);
339  lCurve->Draw("AL");
340  bestLcurve->SetMarkerColor(kRed);
341  bestLcurve->Draw("*");
342 
343  return 0;
344 }
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3127
Random number generator class based on M.
Definition: TRandom3.h:29
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
Definition: Rtypes.h:61
Base class for spline implementation containing the Draw/Paint methods //.
Definition: TSpline.h:22
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:747
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
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:5969
Double_t x[n]
Definition: legend1.C:17
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:31
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:43
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0
ERegMode
Definition: TUnfold.h:107
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:512
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
double gamma(double x)
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2853
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:618
Double_t Pi()
Definition: TMath.h:44
The Canvas class.
Definition: TCanvas.h:41
double Double_t
Definition: RtypesCore.h:55
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TSpline.cxx:100
Double_t y[n]
Definition: legend1.C:17
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:307
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1089
TUnfold is used to decompose a measurement y into several sources x given the measurement uncertainti...
Definition: TUnfold.h:99
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:301
Definition: Rtypes.h:61
virtual Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition: TRandom.cxx:362
static void output(int code)
Definition: gifencode.c:226
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
Double_t Tan(Double_t)
Definition: TMath.h:427
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:296