ROOT logo
#include "TSPlot.h"
#include "TTree.h"
#include "TH1.h"
#include "TCanvas.h"
#include "TFile.h"
#include "TPaveLabel.h"
#include "TPad.h"
#include "TPaveText.h"
#include "Riostream.h"

void TestSPlot()
{
//This tutorial illustrates the use of class TSPlot and of the sPlots method
//
//It is an example of analysis of charmless B decays, performed for BABAR.
//One is dealing with a data sample in which two species are present:
//the first is termed signal and the second background. 
//A maximum Likelihood fit is performed to obtain the two yields N1 and N2
//The fit relies on two discriminating variables collectively denoted y, 
//which are chosen within three possible variables denoted Mes, dE and F.
//The variable which is not incorporated in y, is used as the control variable x.
//The distributions of discriminating variables and more details about the method 
//can be found in the TSPlot class description
//
// NOTE: This script requires a data file "TestSPlot_toyMC.dat".
//       This data file is not distributed in the standard ROOT binary tar file.
//       You can download it from ftp://root.cern.ch/root/TestSPlot_toyMC.dat
//
//Authors: Anna Kreshuk, Muriel Pivc

   TString dir = gSystem->UnixPathName(gInterpreter->GetCurrentMacroName());
   dir.ReplaceAll("TestSPlot.C","");
   dir.ReplaceAll("/./","/");
   TString dataFile = Form("%sTestSPlot_toyMC.dat",dir.Data());
   
   //Read the data and initialize a TSPlot object
   TTree *datatree = new TTree("datatree", "datatree");
   datatree->ReadFile(dataFile, "Mes/D:dE/D:F/D:MesSignal/D:MesBackground/D:\
      dESignal/D:dEBackground/D:FSignal/D:FBackground/D");

   TSPlot *splot = new TSPlot(0, 3, 5420, 2, datatree);

   //Set the selection for data tree
   //Note the order of the variables: 
   //first the control variables (not presented in this example),
   //then the 3 discriminating variables, then their probability distribution 
   //functions for the first species(signal) and then their pdfs for the 
   //second species(background)
   splot->SetTreeSelection("Mes:dE:F:MesSignal:dESignal:FSignal:MesBackground:\
      dEBackground:FBackground");

   //Set the initial estimates of the number of events in each species 
   //- used as initial parameter values for the Minuit likelihood fit
   Int_t ne[2];
   ne[0]=500; ne[1]=5000;
   splot->SetInitialNumbersOfSpecies(ne);

   //Compute the weights
   splot->MakeSPlot();

   //Fill the sPlots
   splot->FillSWeightsHists(25);

   //Now let's look at the sPlots
   //The first two histograms are sPlots for the Mes variable signal and 
   //background. dE and F were chosen as discriminating variables to determine 
   //N1 and N2, through a maximum Likelihood fit, and thus the sPlots for the 
   //control variable Mes, unknown to the fit, was contructed.
   //One can see that the sPlot for signal reproduces the PDF correctly, 
   //even when the latter vanishes.
   //
   //The lower two histograms are sPlots for the F variables signal and 
   //background. dE and Mes were chosen as discriminating variables to 
   //determine N1 and N2, through a maximum Likelihood fit, and thus the 
   //sPlots for the control variable F, unknown to the fit, was contructed.

   TCanvas *myc = new TCanvas("myc", 
   "sPlots of Mes and F signal and background", 800, 600);
   myc->SetFillColor(40);

   TPaveText *pt = new TPaveText(0.02,0.85,0.98,0.98);
   pt->SetFillColor(18);
   pt->SetTextFont(20);
   pt->SetTextColor(4);
   pt->AddText("sPlots of Mes and F signal and background,");
   pt->AddText("obtained by the tutorial TestSPlot.C on BABAR MC \
      data (sPlot_toyMC.fit)");
   TText *t3=pt->AddText(
      "M. Pivk and F. R. Le Diberder, Nucl.Inst.Meth.A, physics/0402083");
   t3->SetTextColor(1);
   t3->SetTextFont(30);
   pt->Draw();

   TPad* pad1 = new TPad("pad1","Mes signal",0.02,0.43,0.48,0.83,33);
   TPad* pad2 = new TPad("pad2","Mes background",0.5,0.43,0.98,0.83,33);
   TPad* pad3 = new TPad("pad3", "F signal", 0.02, 0.02, 0.48, 0.41,33);
   TPad* pad4 = new TPad("pad4", "F background", 0.5, 0.02, 0.98, 0.41,33);
   pad1->Draw();
   pad2->Draw();
   pad3->Draw();
   pad4->Draw();

   pad1->cd();
   pad1->SetGrid();
   TH1D *sweight00 = splot->GetSWeightsHist(-1, 0, 0);
   sweight00->SetTitle("Mes signal");
   sweight00->SetStats(kFALSE);
   sweight00->Draw("e");
   sweight00->SetMarkerStyle(21);
   sweight00->SetMarkerSize(0.7);
   sweight00->SetMarkerColor(2);
   sweight00->SetLineColor(2);
   sweight00->GetXaxis()->SetLabelSize(0.05);
   sweight00->GetYaxis()->SetLabelSize(0.06);
   sweight00->GetXaxis()->SetLabelOffset(0.02);

   pad2->cd();
   pad2->SetGrid();
   TH1D *sweight10 = splot->GetSWeightsHist(-1, 1, 0);
   sweight10->SetTitle("Mes background");
   sweight10->SetStats(kFALSE);
   sweight10->Draw("e");
   sweight10->SetMarkerStyle(21);
   sweight10->SetMarkerSize(0.7);
   sweight10->SetMarkerColor(2);
   sweight10->SetLineColor(2);
   sweight10->GetXaxis()->SetLabelSize(0.05);
   sweight10->GetYaxis()->SetLabelSize(0.06);
   sweight10->GetXaxis()->SetLabelOffset(0.02);
 
   pad3->cd();
   pad3->SetGrid();
   TH1D *sweight02 = splot->GetSWeightsHist(-1, 0, 2);
   sweight02->SetTitle("F signal");
   sweight02->SetStats(kFALSE);
   sweight02->Draw("e");
   sweight02->SetMarkerStyle(21);
   sweight02->SetMarkerSize(0.7);
   sweight02->SetMarkerColor(2);
   sweight02->SetLineColor(2);
   sweight02->GetXaxis()->SetLabelSize(0.06);
   sweight02->GetYaxis()->SetLabelSize(0.06);
   sweight02->GetXaxis()->SetLabelOffset(0.01);

   pad4->cd();
   pad4->SetGrid();
   TH1D *sweight12 = splot->GetSWeightsHist(-1, 1, 2);
   sweight12->SetTitle("F background");
   sweight12->SetStats(kFALSE);
   sweight12->Draw("e");
   sweight12->SetMarkerStyle(21);
   sweight12->SetMarkerSize(0.7);
   sweight12->SetMarkerColor(2);
   sweight12->SetLineColor(2);
   sweight12->GetXaxis()->SetLabelSize(0.06);
   sweight12->GetYaxis()->SetLabelSize(0.06);
   sweight02->GetXaxis()->SetLabelOffset(0.01);
   myc->cd();
}
 TestSPlot.C:1
 TestSPlot.C:2
 TestSPlot.C:3
 TestSPlot.C:4
 TestSPlot.C:5
 TestSPlot.C:6
 TestSPlot.C:7
 TestSPlot.C:8
 TestSPlot.C:9
 TestSPlot.C:10
 TestSPlot.C:11
 TestSPlot.C:12
 TestSPlot.C:13
 TestSPlot.C:14
 TestSPlot.C:15
 TestSPlot.C:16
 TestSPlot.C:17
 TestSPlot.C:18
 TestSPlot.C:19
 TestSPlot.C:20
 TestSPlot.C:21
 TestSPlot.C:22
 TestSPlot.C:23
 TestSPlot.C:24
 TestSPlot.C:25
 TestSPlot.C:26
 TestSPlot.C:27
 TestSPlot.C:28
 TestSPlot.C:29
 TestSPlot.C:30
 TestSPlot.C:31
 TestSPlot.C:32
 TestSPlot.C:33
 TestSPlot.C:34
 TestSPlot.C:35
 TestSPlot.C:36
 TestSPlot.C:37
 TestSPlot.C:38
 TestSPlot.C:39
 TestSPlot.C:40
 TestSPlot.C:41
 TestSPlot.C:42
 TestSPlot.C:43
 TestSPlot.C:44
 TestSPlot.C:45
 TestSPlot.C:46
 TestSPlot.C:47
 TestSPlot.C:48
 TestSPlot.C:49
 TestSPlot.C:50
 TestSPlot.C:51
 TestSPlot.C:52
 TestSPlot.C:53
 TestSPlot.C:54
 TestSPlot.C:55
 TestSPlot.C:56
 TestSPlot.C:57
 TestSPlot.C:58
 TestSPlot.C:59
 TestSPlot.C:60
 TestSPlot.C:61
 TestSPlot.C:62
 TestSPlot.C:63
 TestSPlot.C:64
 TestSPlot.C:65
 TestSPlot.C:66
 TestSPlot.C:67
 TestSPlot.C:68
 TestSPlot.C:69
 TestSPlot.C:70
 TestSPlot.C:71
 TestSPlot.C:72
 TestSPlot.C:73
 TestSPlot.C:74
 TestSPlot.C:75
 TestSPlot.C:76
 TestSPlot.C:77
 TestSPlot.C:78
 TestSPlot.C:79
 TestSPlot.C:80
 TestSPlot.C:81
 TestSPlot.C:82
 TestSPlot.C:83
 TestSPlot.C:84
 TestSPlot.C:85
 TestSPlot.C:86
 TestSPlot.C:87
 TestSPlot.C:88
 TestSPlot.C:89
 TestSPlot.C:90
 TestSPlot.C:91
 TestSPlot.C:92
 TestSPlot.C:93
 TestSPlot.C:94
 TestSPlot.C:95
 TestSPlot.C:96
 TestSPlot.C:97
 TestSPlot.C:98
 TestSPlot.C:99
 TestSPlot.C:100
 TestSPlot.C:101
 TestSPlot.C:102
 TestSPlot.C:103
 TestSPlot.C:104
 TestSPlot.C:105
 TestSPlot.C:106
 TestSPlot.C:107
 TestSPlot.C:108
 TestSPlot.C:109
 TestSPlot.C:110
 TestSPlot.C:111
 TestSPlot.C:112
 TestSPlot.C:113
 TestSPlot.C:114
 TestSPlot.C:115
 TestSPlot.C:116
 TestSPlot.C:117
 TestSPlot.C:118
 TestSPlot.C:119
 TestSPlot.C:120
 TestSPlot.C:121
 TestSPlot.C:122
 TestSPlot.C:123
 TestSPlot.C:124
 TestSPlot.C:125
 TestSPlot.C:126
 TestSPlot.C:127
 TestSPlot.C:128
 TestSPlot.C:129
 TestSPlot.C:130
 TestSPlot.C:131
 TestSPlot.C:132
 TestSPlot.C:133
 TestSPlot.C:134
 TestSPlot.C:135
 TestSPlot.C:136
 TestSPlot.C:137
 TestSPlot.C:138
 TestSPlot.C:139
 TestSPlot.C:140
 TestSPlot.C:141
 TestSPlot.C:142
 TestSPlot.C:143
 TestSPlot.C:144
 TestSPlot.C:145
 TestSPlot.C:146
 TestSPlot.C:147
 TestSPlot.C:148
 TestSPlot.C:149
 TestSPlot.C:150
 TestSPlot.C:151
 TestSPlot.C:152
 TestSPlot.C:153
 TestSPlot.C:154
 TestSPlot.C:155
 TestSPlot.C:156
 TestSPlot.C:157
 TestSPlot.C:158
 TestSPlot.C:159
 TestSPlot.C:160
thumb