TestSPlot.C

Go to the documentation of this file.
00001 #include "TSPlot.h"
00002 #include "TTree.h"
00003 #include "TH1.h"
00004 #include "TCanvas.h"
00005 #include "TFile.h"
00006 #include "TPaveLabel.h"
00007 #include "TPad.h"
00008 #include "TPaveText.h"
00009 #include "Riostream.h"
00010 
00011 void TestSPlot()
00012 {
00013 //This tutorial illustrates the use of class TSPlot and of the sPlots method
00014 //
00015 //It is an example of analysis of charmless B decays, performed for BABAR.
00016 //One is dealing with a data sample in which two species are present:
00017 //the first is termed signal and the second background. 
00018 //A maximum Likelihood fit is performed to obtain the two yields N1 and N2
00019 //The fit relies on two discriminating variables collectively denoted y, 
00020 //which are chosen within three possible variables denoted Mes, dE and F.
00021 //The variable which is not incorporated in y, is used as the control variable x.
00022 //The distributions of discriminating variables and more details about the method 
00023 //can be found in the TSPlot class description
00024 //
00025 // NOTE: This script requires a data file "TestSPlot_toyMC.dat".
00026 //       This data file is not distributed in the standard ROOT binary tar file.
00027 //       You can download it from ftp://root.cern.ch/root/TestSPlot_toyMC.dat
00028 //
00029 //Authors: Anna Kreshuk, Muriel Pivc
00030 
00031    TString dir = gSystem->UnixPathName(gInterpreter->GetCurrentMacroName());
00032    dir.ReplaceAll("TestSPlot.C","");
00033    dir.ReplaceAll("/./","/");
00034    TString dataFile = Form("%sTestSPlot_toyMC.dat",dir.Data());
00035    
00036    //Read the data and initialize a TSPlot object
00037    TTree *datatree = new TTree("datatree", "datatree");
00038    datatree->ReadFile(dataFile, 
00039       "Mes/D:dE/D:F/D:MesSignal/D:MesBackground/D:"
00040       "dESignal/D:dEBackground/D:FSignal/D:FBackground/D");
00041 
00042    TSPlot *splot = new TSPlot(0, 3, 5420, 2, datatree);
00043 
00044    //Set the selection for data tree
00045    //Note the order of the variables: 
00046    //first the control variables (not presented in this example),
00047    //then the 3 discriminating variables, then their probability distribution 
00048    //functions for the first species(signal) and then their pdfs for the 
00049    //second species(background)
00050    splot->SetTreeSelection(
00051       "Mes:dE:F:MesSignal:dESignal:FSignal:MesBackground:"
00052       "dEBackground:FBackground");
00053 
00054    //Set the initial estimates of the number of events in each species 
00055    //- used as initial parameter values for the Minuit likelihood fit
00056    Int_t ne[2];
00057    ne[0]=500; ne[1]=5000;
00058    splot->SetInitialNumbersOfSpecies(ne);
00059 
00060    //Compute the weights
00061    splot->MakeSPlot();
00062 
00063    //Fill the sPlots
00064    splot->FillSWeightsHists(25);
00065 
00066    //Now let's look at the sPlots
00067    //The first two histograms are sPlots for the Mes variable signal and 
00068    //background. dE and F were chosen as discriminating variables to determine 
00069    //N1 and N2, through a maximum Likelihood fit, and thus the sPlots for the 
00070    //control variable Mes, unknown to the fit, was contructed.
00071    //One can see that the sPlot for signal reproduces the PDF correctly, 
00072    //even when the latter vanishes.
00073    //
00074    //The lower two histograms are sPlots for the F variables signal and 
00075    //background. dE and Mes were chosen as discriminating variables to 
00076    //determine N1 and N2, through a maximum Likelihood fit, and thus the 
00077    //sPlots for the control variable F, unknown to the fit, was contructed.
00078 
00079    TCanvas *myc = new TCanvas("myc", 
00080    "sPlots of Mes and F signal and background", 800, 600);
00081    myc->SetFillColor(40);
00082 
00083    TPaveText *pt = new TPaveText(0.02,0.85,0.98,0.98);
00084    pt->SetFillColor(18);
00085    pt->SetTextFont(20);
00086    pt->SetTextColor(4);
00087    pt->AddText("sPlots of Mes and F signal and background,");
00088    pt->AddText("obtained by the tutorial TestSPlot.C on BABAR MC "
00089                "data (sPlot_toyMC.fit)");
00090    TText *t3=pt->AddText(
00091       "M. Pivk and F. R. Le Diberder, Nucl.Inst.Meth.A, physics/0402083");
00092    t3->SetTextColor(1);
00093    t3->SetTextFont(30);
00094    pt->Draw();
00095 
00096    TPad* pad1 = new TPad("pad1","Mes signal",0.02,0.43,0.48,0.83,33);
00097    TPad* pad2 = new TPad("pad2","Mes background",0.5,0.43,0.98,0.83,33);
00098    TPad* pad3 = new TPad("pad3", "F signal", 0.02, 0.02, 0.48, 0.41,33);
00099    TPad* pad4 = new TPad("pad4", "F background", 0.5, 0.02, 0.98, 0.41,33);
00100    pad1->Draw();
00101    pad2->Draw();
00102    pad3->Draw();
00103    pad4->Draw();
00104 
00105    pad1->cd();
00106    pad1->SetGrid();
00107    TH1D *sweight00 = splot->GetSWeightsHist(-1, 0, 0);
00108    sweight00->SetTitle("Mes signal");
00109    sweight00->SetStats(kFALSE);
00110    sweight00->Draw("e");
00111    sweight00->SetMarkerStyle(21);
00112    sweight00->SetMarkerSize(0.7);
00113    sweight00->SetMarkerColor(2);
00114    sweight00->SetLineColor(2);
00115    sweight00->GetXaxis()->SetLabelSize(0.05);
00116    sweight00->GetYaxis()->SetLabelSize(0.06);
00117    sweight00->GetXaxis()->SetLabelOffset(0.02);
00118 
00119    pad2->cd();
00120    pad2->SetGrid();
00121    TH1D *sweight10 = splot->GetSWeightsHist(-1, 1, 0);
00122    sweight10->SetTitle("Mes background");
00123    sweight10->SetStats(kFALSE);
00124    sweight10->Draw("e");
00125    sweight10->SetMarkerStyle(21);
00126    sweight10->SetMarkerSize(0.7);
00127    sweight10->SetMarkerColor(2);
00128    sweight10->SetLineColor(2);
00129    sweight10->GetXaxis()->SetLabelSize(0.05);
00130    sweight10->GetYaxis()->SetLabelSize(0.06);
00131    sweight10->GetXaxis()->SetLabelOffset(0.02);
00132  
00133    pad3->cd();
00134    pad3->SetGrid();
00135    TH1D *sweight02 = splot->GetSWeightsHist(-1, 0, 2);
00136    sweight02->SetTitle("F signal");
00137    sweight02->SetStats(kFALSE);
00138    sweight02->Draw("e");
00139    sweight02->SetMarkerStyle(21);
00140    sweight02->SetMarkerSize(0.7);
00141    sweight02->SetMarkerColor(2);
00142    sweight02->SetLineColor(2);
00143    sweight02->GetXaxis()->SetLabelSize(0.06);
00144    sweight02->GetYaxis()->SetLabelSize(0.06);
00145    sweight02->GetXaxis()->SetLabelOffset(0.01);
00146 
00147    pad4->cd();
00148    pad4->SetGrid();
00149    TH1D *sweight12 = splot->GetSWeightsHist(-1, 1, 2);
00150    sweight12->SetTitle("F background");
00151    sweight12->SetStats(kFALSE);
00152    sweight12->Draw("e");
00153    sweight12->SetMarkerStyle(21);
00154    sweight12->SetMarkerSize(0.7);
00155    sweight12->SetMarkerColor(2);
00156    sweight12->SetLineColor(2);
00157    sweight12->GetXaxis()->SetLabelSize(0.06);
00158    sweight12->GetYaxis()->SetLabelSize(0.06);
00159    sweight02->GetXaxis()->SetLabelOffset(0.01);
00160    myc->cd();
00161 }

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