TSVDUnfoldExample.C

Go to the documentation of this file.
00001 //  Data unfolding using Singular Value Decomposition 
00002 // 
00003 //////////////////////////////////////////////////////////////////////////
00004 //                                                                      //
00005 // TSVDUnfold example                                                   //
00006 //                                                                      //
00007 // Data unfolding using Singular Value Decomposition (hep-ph/9509307)   //
00008 // Authors: Kerstin Tackmann, Andreas Hoecker, Heiko Lacker             //
00009 //                                                                      //
00010 // Example distribution and smearing model from Tim Adye (RAL)          //
00011 //                                                                      //
00012 // Nov 23, 2010                                                         //
00013 //////////////////////////////////////////////////////////////////////////
00014 
00015 #include <iostream>
00016 
00017 #include "TROOT.h"
00018 #include "TSystem.h"
00019 #include "TStyle.h"
00020 #include "TRandom3.h"
00021 #include "TString.h"
00022 #include "TMath.h"
00023 #include "TH1D.h"
00024 #include "TH2D.h"
00025 #include "TLegend.h"
00026 #include "TCanvas.h"
00027 #include "TColor.h"
00028 #include "TLine.h"
00029 
00030 #if not defined(__CINT__) || defined(__MAKECINT__)
00031 #include "TSVDUnfold.h"
00032 #endif
00033 
00034 Double_t Reconstruct( Double_t xt, TRandom3& R )
00035 {
00036    // apply some Gaussian smearing + bias and efficiency corrections to fake reconstruction
00037    const Double_t cutdummy = -99999.0;
00038    Double_t xeff = 0.3 + (1.0 - 0.3)/20.0*(xt + 10.0);  // efficiency
00039    Double_t x    = R.Rndm();
00040    if (x > xeff) return cutdummy;
00041    else {
00042       Double_t xsmear= R.Gaus(-2.5,0.2); // bias and smear
00043       return xt+xsmear;
00044    }
00045 }
00046 
00047 void TSVDUnfoldExample() 
00048 {
00049    gROOT->Reset();
00050    gROOT->SetStyle("Plain");
00051    gStyle->SetOptStat(0);
00052 
00053    TRandom3 R;
00054 
00055    const Double_t cutdummy= -99999.0;
00056    
00057    // --- Data/MC toy generation -----------------------------------
00058 
00059    // The MC input
00060    Int_t nbins = 40;
00061    TH1D *xini = new TH1D("xini", "MC truth", nbins, -10.0, 10.0);
00062    TH1D *bini = new TH1D("bini", "MC reco", nbins, -10.0, 10.0);
00063    TH2D *Adet = new TH2D("Adet", "detector response", nbins, -10.0, 10.0, nbins, -10.0, 10.0);
00064 
00065    // Data
00066    TH1D *data = new TH1D("data", "data", nbins, -10.0, 10.0);
00067    // Data "truth" distribution to test the unfolding
00068    TH1D *datatrue = new TH1D("datatrue", "data truth", nbins, -10.0, 10.0);
00069    // Statistical covariance matrix
00070    TH2D *statcov = new TH2D("statcov", "covariance matrix", nbins, -10.0, 10.0, nbins, -10.0, 10.0);
00071 
00072    // Fill the MC using a Breit-Wigner, mean 0.3 and width 2.5.
00073    for (Int_t i= 0; i<100000; i++) {
00074       Double_t xt = R.BreitWigner(0.3, 2.5);
00075       xini->Fill(xt);
00076       Double_t x = Reconstruct( xt, R );
00077       if (x != cutdummy) {
00078          Adet->Fill(x, xt);
00079          bini->Fill(x);
00080       }
00081    }
00082 
00083    // Fill the "data" with a Gaussian, mean 0 and width 2.
00084    for (Int_t i=0; i<10000; i++) {
00085       Double_t xt = R.Gaus(0.0, 2.0);
00086       datatrue->Fill(xt);
00087       Double_t x = Reconstruct( xt, R );
00088       if (x != cutdummy) data->Fill(x);
00089    }
00090 
00091    cout << "Created toy distributions and errors for: " << endl;
00092    cout << "... \"true MC\"   and \"reconstructed (smeared) MC\"" << endl;
00093    cout << "... \"true data\" and \"reconstructed (smeared) data\"" << endl;
00094    cout << "... the \"detector response matrix\"" << endl;
00095 
00096    // Fill the data covariance matrix
00097    for (int i=1; i<=data->GetNbinsX(); i++) {
00098       statcov->SetBinContent(i,i,data->GetBinError(i)*data->GetBinError(i)); 
00099    }
00100 
00101    // --- Here starts the actual unfolding -------------------------
00102 
00103    // Create TSVDUnfold object and initialise
00104    TSVDUnfold *tsvdunf = new TSVDUnfold( data, bini, xini, Adet );
00105 
00106    // It is possible to normalise unfolded spectrum to unit area
00107    tsvdunf->SetNormalize( kFALSE ); // no normalisation here
00108 
00109    // Perform the unfolding with regularisation parameter kreg = 13
00110    // - the larger kreg, the finer grained the unfolding, but the more fluctuations occur
00111    // - the smaller kreg, the stronger is the regularisation and the bias
00112    TH1D* unfres = tsvdunf->Unfold( 13 );
00113 
00114    // Get the distribution of the d to cross check the regularization
00115    // - choose kreg to be the point where |d_i| stop being statistically significantly >>1
00116    TH1D* ddist = tsvdunf->GetD();
00117 
00118    // Get the distribution of the singular values
00119    TH1D* svdist = tsvdunf->GetSV();
00120 
00121    // Compute the error matrix for the unfolded spectrum using toy MC
00122    // using the measured covariance matrix as input to generate the toys
00123    // 1000 toys is usually reasonable
00124    TH2D* ustatcov = tsvdunf->GetUnfoldCovMatrix( statcov, 1000 );   
00125 
00126    // Now compute the error matrix on the unfolded distribution originating
00127    // from the finite detector matrix statistics
00128    TH2D* uadetcov = tsvdunf->GetAdetCovMatrix( 1000 );   
00129 
00130    // Sum up the two (they are uncorrelated)
00131    ustatcov->Add( uadetcov );
00132    
00133    // --- Only plotting stuff below ------------------------------
00134 
00135    for (int i=1; i<=unfres->GetNbinsX(); i++) {
00136       unfres->SetBinError(i, TMath::Sqrt(ustatcov->GetBinContent(i,i)));
00137    }
00138 
00139    xini->Scale(0.7*datatrue->Integral()/xini->Integral());
00140 
00141    TLegend *leg = new TLegend(0.58,0.68,0.99,0.88);
00142    leg->SetBorderSize(0);
00143    leg->SetFillColor(0);
00144    leg->SetFillStyle(0);
00145    leg->AddEntry(unfres,"Unfolded Data","p");
00146    leg->AddEntry(datatrue,"True Data","l");
00147    leg->AddEntry(data,"Reconstructed Data","l");
00148    leg->AddEntry(xini,"True MC","l");
00149 
00150    TCanvas *c1 = new TCanvas( "c1", "Unfolding toy example with TSVDUnfold", 900, 800 );
00151 
00152    // --- Style settings -----------------------------------------
00153    Int_t c_Canvas    = TColor::GetColor( "#f0f0f0" );
00154    Int_t c_FrameFill = TColor::GetColor( "#fffffd" );
00155    Int_t c_TitleBox  = TColor::GetColor( "#6D7B8D" );
00156    Int_t c_TitleText = TColor::GetColor( "#FFFFFF" );
00157 
00158    c1->SetFrameFillColor( c_FrameFill );
00159    c1->SetFillColor     ( c_Canvas    );
00160    c1->Divide(1,2);
00161    c1->cd(1);
00162    
00163 
00164    gStyle->SetTitleFillColor( c_TitleBox  );
00165    gStyle->SetTitleTextColor( c_TitleText );
00166    gStyle->SetTitleBorderSize( 1 );
00167    gStyle->SetTitleH( 0.052 );
00168    gStyle->SetTitleX( c1->GetLeftMargin() );
00169    gStyle->SetTitleY( 1 - c1->GetTopMargin() + gStyle->GetTitleH() );
00170    gStyle->SetTitleW( 1 - c1->GetLeftMargin() - c1->GetRightMargin() );
00171 
00172    TH1D* frame = new TH1D( *unfres );
00173    frame->SetTitle( "Unfolding toy example with TSVDUnfold" );
00174    frame->GetXaxis()->SetTitle( "x variable" );
00175    frame->GetYaxis()->SetTitle( "Events" );
00176    frame->GetXaxis()->SetTitleOffset( 1.25 );
00177    frame->GetYaxis()->SetTitleOffset( 1.29 );
00178    frame->Draw();
00179 
00180    data->SetLineStyle(2);
00181    data->SetLineColor(4);
00182    data->SetLineWidth(2);
00183    unfres->SetMarkerStyle(20);
00184    datatrue->SetLineColor(2);
00185    datatrue->SetLineWidth(2);
00186    xini->SetLineStyle(2);
00187    xini->SetLineColor(8);
00188    xini->SetLineWidth(2);
00189    // ------------------------------------------------------------
00190 
00191    // add histograms
00192    unfres->Draw("same");
00193    datatrue->Draw("same");
00194    data->Draw("same");
00195    xini->Draw("same");
00196 
00197    leg->Draw();
00198 
00199 
00200    // distribution of the d quantity
00201    gStyle->SetOptLogy();
00202    TVirtualPad * c12 = c1->cd(2);
00203    c12->Divide(2,1);
00204    TVirtualPad * c2 = c12->cd(1);
00205    c2->SetFrameFillColor( c_FrameFill );
00206    c2->SetFillColor     ( c_Canvas    );
00207 
00208    TLine *line = new TLine( 0.,1.,40.,1. );
00209    line->SetLineStyle(2);
00210 
00211    TH1D* dframe = new TH1D( *ddist );
00212    dframe->SetTitle( "TSVDUnfold |d_{i}| for toy example" );
00213    dframe->GetXaxis()->SetTitle( "i" );
00214    dframe->GetYaxis()->SetTitle( "|d_{i}|" );
00215    dframe->GetXaxis()->SetTitleOffset( 1.25 );
00216    dframe->GetYaxis()->SetTitleOffset( 1.29 );
00217    dframe->SetMinimum( 0.001 );
00218    dframe->Draw();
00219 
00220    ddist->SetLineWidth( 2 );
00221    ddist->Draw( "same" );
00222    line->Draw();
00223 
00224    // distribution of the singular values
00225    gStyle->SetOptLogy();
00226    TVirtualPad * c3 = c12->cd(2);
00227    c3->SetFrameFillColor( c_FrameFill );
00228    c3->SetFillColor     ( c_Canvas    );
00229 
00230    TH1D* svframe = new TH1D( *svdist );
00231    svframe->SetTitle( "TSVDUnfold singular values for toy example" );
00232    svframe->GetXaxis()->SetTitle( "i" );
00233    svframe->GetYaxis()->SetTitle( "|s_{i}|" );
00234    svframe->GetXaxis()->SetTitleOffset( 1.25 );
00235    svframe->GetYaxis()->SetTitleOffset( 1.29 );
00236    svframe->SetMinimum( 0.001 );
00237    svframe->Draw();
00238 
00239    svdist->SetLineWidth( 2 );
00240    svdist->Draw( "same" );
00241    //c1->Print("TSVDUnfoldExample.pdf");
00242 }

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