testUnfold2.C

Go to the documentation of this file.
00001 // TUnfold test program as an example for a more complex regularisation scheme
00002 // Author: Stefan Schmitt
00003 // DESY, 14.10.2008
00004 
00005 //  Version 16, parallel to changes in TUnfold
00006 //
00007 //  History:
00008 //    Version 15, with automatic L-curve scan, simplified example
00009 //    Version 14, with changes in TUnfoldSys.cxx
00010 //    Version 13,  with changes to TUnfold.C
00011 //    Version 12,  with improvements to TUnfold.cxx
00012 //    Version 11,  print chi**2 and number of degrees of freedom
00013 //    Version 10, with bug-fix in TUnfold.cxx
00014 //    Version 9, with bug-fix in TUnfold.cxx, TUnfold.h
00015 //    Version 8, with bug-fix in TUnfold.cxx, TUnfold.h
00016 //    Version 7, with bug-fix in TUnfold.cxx, TUnfold.h
00017 //    Version 6a, fix problem with dynamic array allocation under windows
00018 //    Version 6, re-include class MyUnfold in the example
00019 //    Version 5, move class MyUnfold to seperate files
00020 //    Version 4, with bug-fix in TUnfold.C
00021 //    Version 3, with bug-fix in TUnfold.C
00022 //    Version 2, with changed ScanLcurve() arguments
00023 //    Version 1, remove L curve analysis, use ScanLcurve() method instead
00024 //    Version 0, L curve analysis included here
00025 
00026 #include <TMath.h>
00027 #include <TCanvas.h>
00028 #include <TRandom3.h>
00029 #include <TFitter.h>
00030 #include <TF1.h>
00031 #include <TStyle.h>
00032 #include <TVector.h>
00033 #include <TGraph.h>
00034 #include <TUnfold.h>
00035 
00036 using namespace std;
00037 
00038 ///////////////////////////////////////////////////////////////////////
00039 // 
00040 //  Test program as an example for a more complex regularisation scheme
00041 //
00042 //  (1) Generate Monte Carlo and Data events
00043 //      The events consist of
00044 //        signal
00045 //        background
00046 //
00047 //      The signal is a resonance. It is generated with a Breit-Wigner,
00048 //      smeared by a Gaussian
00049 //
00050 //  (2) Unfold the data. The result is:
00051 //      The background level
00052 //      The shape of the resonance, corrected for detector effects
00053 //
00054 //      The regularisation is done on the curvature, excluding the bins
00055 //      near the peak.
00056 //
00057 //  (3) produce some plots
00058 //
00059 ///////////////////////////////////////////////////////////////////////
00060 
00061 TRandom *rnd=0;
00062 
00063 // generate an event
00064 // output:
00065 //  negative mass: background event
00066 //  positive mass: signal event
00067 Double_t GenerateEvent(Double_t bgr, // relative fraction of background
00068                        Double_t mass, // peak position
00069                        Double_t gamma) // peak width
00070 {
00071   Double_t t;
00072   if(rnd->Rndm()>bgr) {
00073     // generate signal event
00074     // with positive mass
00075     do {
00076       do {
00077         t=rnd->Rndm();
00078       } while(t>=1.0); 
00079       t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
00080     } while(t<=0.0);
00081     return t;
00082   } else {
00083     // generate background event
00084     // generate events following a power-law distribution
00085     //   f(E) = K * TMath::power((E0+E),N0)
00086     static Double_t const E0=2.4;
00087     static Double_t const N0=2.9;
00088     do {
00089       do {
00090         t=rnd->Rndm();
00091       } while(t>=1.0);
00092       // the mass is returned negative
00093       // In our example a convenient way to indicate it is a background event.
00094       t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
00095     } while(t>=0.0);
00096     return t;
00097   }
00098 }
00099 
00100 // smear the event to detector level
00101 // input:
00102 //   mass on generator level (mTrue>0 !)
00103 // output:
00104 //   mass on detector level
00105 Double_t DetectorEvent(Double_t mTrue) {
00106   // smear by double-gaussian
00107   static Double_t frac=0.1;
00108   static Double_t wideBias=0.03;
00109   static Double_t wideSigma=0.5;
00110   static Double_t smallBias=0.0;
00111   static Double_t smallSigma=0.1;
00112   if(rnd->Rndm()>frac) {
00113     return rnd->Gaus(mTrue+smallBias,smallSigma);
00114   } else {
00115     return rnd->Gaus(mTrue+wideBias,wideSigma);
00116   }
00117 }
00118 
00119 int testUnfold2() 
00120 {
00121   // switch on histogram errors
00122   TH1::SetDefaultSumw2();
00123 
00124   // random generator
00125   rnd=new TRandom3();
00126 
00127   // data and MC luminosity, cross-section
00128   Double_t const luminosityData=100000;
00129   Double_t const luminosityMC=1000000;
00130   Double_t const crossSection=1.0;
00131 
00132   Int_t const nDet=250;
00133   Int_t const nGen=100;
00134   Double_t const xminDet=0.0;
00135   Double_t const xmaxDet=10.0;
00136   Double_t const xminGen=0.0;
00137   Double_t const xmaxGen=10.0;
00138 
00139   //============================================
00140   // generate MC distribution
00141   //
00142   TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
00143   TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
00144   TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",nDet,xminDet,xmaxDet,
00145                               nGen,xminGen,xmaxGen);
00146   Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
00147   for(Int_t i=0;i<neventMC;i++) {
00148     Double_t mGen=GenerateEvent(0.3, // relative fraction of background
00149                                 4.0, // peak position in MC
00150                                 0.2); // peak width in MC
00151     Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00152     // the generated mass is negative for background
00153     // and positive for signal
00154     // so it will be filled in the underflow bin
00155     // this is very convenient for the unfolding:
00156     // the unfolded result will contain the number of background
00157     // events in the underflow bin
00158 
00159     // generated MC distribution (for comparison only)
00160     histMgenMC->Fill(mGen,luminosityData/luminosityMC);
00161     // reconstructed MC distribution (for comparison only)
00162     histMdetMC->Fill(mDet,luminosityData/luminosityMC);
00163 
00164     // matrix describing how the generator input migrates to the
00165     // reconstructed level. Unfolding input.
00166     // NOTE on underflow/overflow bins:
00167     //  (1) the detector level under/overflow bins are used for
00168     //       normalisation ("efficiency" correction)
00169     //       in our toy example, these bins are populated from tails
00170     //       of the initial MC distribution.
00171     //  (2) the generator level underflow/overflow bins are
00172     //       unfolded. In this example:
00173     //       underflow bin: background events reconstructed in the detector
00174     //       overflow bin: signal events generated at masses > xmaxDet
00175     // for the unfolded result these bins will be filled
00176     //  -> the background normalisation will be contained in the underflow bin
00177     histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
00178   }
00179 
00180   //============================================
00181   // generate data distribution
00182   //
00183   TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
00184   TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
00185   Int_t neventData=rnd->Poisson(luminosityData*crossSection);
00186   for(Int_t i=0;i<neventData;i++) {
00187     Double_t mGen=GenerateEvent(0.4, // relative fraction of background
00188                                 3.8, // peak position
00189                                 0.15); // peak width
00190     Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00191     // generated data mass for comparison plots
00192     // for real data, we do not have this histogram
00193     histMgenData->Fill(mGen);
00194 
00195     // reconstructed mass, unfolding input
00196     histMdetData->Fill(mDet);
00197   }
00198 
00199   //=========================================================================
00200   // set up the unfolding
00201   TUnfold unfold(histMdetGenMC,TUnfold::kHistMapOutputVert,
00202                  TUnfold::kRegModeNone);
00203   // regularisation
00204   //----------------
00205   // the regularisation is done on the curvature (2nd derivative) of
00206   // the output distribution
00207   //
00208   // One has to exclude the bins near the peak of the Breit-Wigner,
00209   // because there the curvature is high
00210   // (and the regularisation eventually could enforce a small
00211   //  curvature, thus biasing result)
00212   //
00213   // in real life, the parameters below would have to be optimized,
00214   // depending on the data peak position and width
00215   // Or maybe one finds a different regularisation scheme... this is
00216   // just an example...
00217   Double_t estimatedPeakPosition=3.8;
00218   Int_t nPeek=3;
00219   TUnfold::ERegMode regMode=TUnfold::kRegModeCurvature;
00220   // calculate bin number correspoinding to estimated peak position
00221   Int_t iPeek=(Int_t)(nGen*(estimatedPeakPosition-xminGen)/(xmaxGen-xminGen)
00222                       // offset 1.5
00223                       // accounts for start bin 1
00224                       // and rounding errors +0.5
00225                       +1.5);
00226   // regularize output bins 1..iPeek-nPeek
00227   unfold.RegularizeBins(1,1,iPeek-nPeek,regMode);
00228   // regularize output bins iPeek+nPeek..nGen
00229   unfold.RegularizeBins(iPeek+nPeek,1,nGen-(iPeek+nPeek),regMode);
00230 
00231   // unfolding
00232   //-----------
00233 
00234   // set input distribution and bias scale (=0)
00235   if(unfold.SetInput(histMdetData,0.0)>=10000) {
00236     std::cout<<"Unfolding result may be wrong\n";
00237   }
00238 
00239   // do the unfolding here
00240   Double_t tauMin=0.0;
00241   Double_t tauMax=0.0;
00242   Int_t nScan=30;
00243   Int_t iBest;
00244   TSpline *logTauX,*logTauY;
00245   TGraph *lCurve;
00246   // this method scans the parameter tau and finds the kink in the L curve
00247   // finally, the unfolding is done for the "best" choice of tau
00248   iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
00249   std::cout<<"tau="<<unfold.GetTau()<<"\n";  
00250   std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
00251            <<" / "<<unfold.GetNdf()<<"\n";
00252 
00253   // save point corresponding to the kink in the L curve as TGraph
00254   Double_t t[1],x[1],y[1];
00255   logTauX->GetKnot(iBest,t[0],x[0]);
00256   logTauY->GetKnot(iBest,t[0],y[0]);
00257   TGraph *bestLcurve=new TGraph(1,x,y);
00258   TGraph *bestLogTauX=new TGraph(1,t,x);
00259 
00260   //============================================================
00261   // extract unfolding results into histograms
00262 
00263   // set up a bin map, excluding underflow and overflow bins
00264   // the binMap relates the the output of the unfolding to the final
00265   // histogram bins
00266   Int_t *binMap=new Int_t[nGen+2];
00267   for(Int_t i=1;i<=nGen;i++) binMap[i]=i;
00268   binMap[0]=-1;
00269   binMap[nGen+1]=-1;
00270 
00271   TH1D *histMunfold=new TH1D("Unfolded",";mass(gen)",nGen,xminGen,xmaxGen);
00272   unfold.GetOutput(histMunfold,binMap);
00273   TH1D *histMdetFold=unfold.GetFoldedOutput("FoldedBack","mass(det)",
00274                                               xminDet,xmaxDet);
00275 
00276   // store global correlation coefficients
00277   TH1D *histRhoi=new TH1D("rho_I","mass",nGen,xminGen,xmaxGen);  
00278   unfold.GetRhoI(histRhoi,0,binMap);
00279 
00280   delete[] binMap;
00281   binMap=0;
00282 
00283   //=====================================================================
00284   // plot some histograms
00285   TCanvas output;
00286 
00287   // produce some plots
00288   output.Divide(3,2);
00289 
00290   // Show the matrix which connects input and output
00291   // There are overflow bins at the bottom, not shown in the plot
00292   // These contain the background shape.
00293   // The overflow bins to the left and right contain
00294   // events which are not reconstructed. These are necessary for proper MC
00295   // normalisation
00296   output.cd(1);
00297   histMdetGenMC->Draw("BOX");
00298 
00299   // draw generator-level distribution:
00300   //   data (red) [for real data this is not available]
00301   //   MC input (black) [with completely wrong peak position and shape]
00302   //   unfolded data (blue)
00303   output.cd(2);
00304   histMunfold->SetLineColor(kBlue);
00305   histMunfold->Draw();
00306   histMgenData->SetLineColor(kRed);
00307   histMgenData->Draw("SAME");
00308   histMgenMC->Draw("SAME HIST");
00309 
00310   // show detector level distributions
00311   //    data (red)
00312   //    MC (black)
00313   //    unfolded data (blue)
00314   output.cd(3);
00315   histMdetFold->SetLineColor(kBlue);
00316   histMdetFold->Draw();
00317   histMdetData->SetLineColor(kRed);
00318   histMdetData->Draw("SAME");
00319   histMdetMC->Draw("SAME HIST");
00320 
00321   // show correlation coefficients
00322   //     all bins outside the peak are found to be highly correlated
00323   //     But they are compatible with zero anyway
00324   //     If the peak shape is fitted,
00325   //     these correlations have to be taken into account, see example
00326   output.cd(4);
00327   histRhoi->Draw();
00328 
00329   // show rhoi_max(tau) distribution
00330   output.cd(5);
00331   logTauX->Draw();
00332   bestLogTauX->SetMarkerColor(kRed);
00333   bestLogTauX->Draw("*");
00334 
00335   output.cd(6);
00336   lCurve->Draw("AL");
00337   bestLcurve->SetMarkerColor(kRed);
00338   bestLcurve->Draw("*");
00339  
00340   output.SaveAs("testUnfold2.ps");
00341   return 0;
00342 }

Generated on Tue Jul 5 15:44:51 2011 for ROOT_528-00b_version by  doxygen 1.5.1