Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TestSPlot.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_splot
3/// This tutorial illustrates the use of class TSPlot and of the sPlots method
4///
5/// It is an example of analysis of charmless B decays, performed for BABAR.
6/// One is dealing with a data sample in which two species are present:
7/// the first is termed signal and the second background.
8/// A maximum Likelihood fit is performed to obtain the two yields N1 and N2
9/// The fit relies on two discriminating variables collectively denoted y,
10/// which are chosen within three possible variables denoted Mes, dE and F.
11/// The variable which is not incorporated in y, is used as the control variable x.
12/// The distributions of discriminating variables and more details about the method
13/// can be found in the TSPlot class description
14///
15/// NOTE: This script requires a data file `$ROOTSYS/tutorials/splot/TestSPlot_toyMC.dat`.
16///
17/// \macro_image
18/// \macro_output
19/// \macro_code
20///
21/// \authors Anna Kreshuk, Muriel Pivc
22
23#include "TSPlot.h"
24#include "TTree.h"
25#include "TH1.h"
26#include "TCanvas.h"
27#include "TFile.h"
28#include "TPaveLabel.h"
29#include "TPad.h"
30#include "TPaveText.h"
31#include "Riostream.h"
32
33void TestSPlot()
34{
36 dir.ReplaceAll("TestSPlot.C", "");
37 dir.ReplaceAll("/./", "/");
38 TString dataFile = Form("%sTestSPlot_toyMC.dat", dir.Data());
39
40 // Read the data and initialize a TSPlot object
41 TTree *datatree = new TTree("datatree", "datatree");
42 datatree->ReadFile(
43 dataFile, "Mes/D:dE/D:F/D:MesSignal/D:MesBackground/D:dESignal/D:dEBackground/D:FSignal/D:FBackground/D", ' ');
44
45 TSPlot *splot = new TSPlot(0, 3, 5420, 2, datatree);
46
47 // Set the selection for data tree
48 // Note the order of the variables:
49 // first the control variables (not presented in this example),
50 // then the 3 discriminating variables, then their probability distribution
51 // functions for the first species(signal) and then their pdfs for the
52 // second species(background)
53 splot->SetTreeSelection("Mes:dE:F:MesSignal:dESignal:FSignal:MesBackground:"
54 "dEBackground:FBackground");
55
56 // Set the initial estimates of the number of events in each species
57 //- used as initial parameter values for the Minuit likelihood fit
58 Int_t ne[2];
59 ne[0] = 500;
60 ne[1] = 5000;
61 splot->SetInitialNumbersOfSpecies(ne);
62
63 // Compute the weights
64 splot->MakeSPlot();
65
66 // Fill the sPlots
67 splot->FillSWeightsHists(25);
68
69 // Now let's look at the sPlots
70 // The first two histograms are sPlots for the Mes variable signal and
71 // background. dE and F were chosen as discriminating variables to determine
72 // N1 and N2, through a maximum Likelihood fit, and thus the sPlots for the
73 // control variable Mes, unknown to the fit, was constructed.
74 // One can see that the sPlot for signal reproduces the PDF correctly,
75 // even when the latter vanishes.
76 //
77 // The lower two histograms are sPlots for the F variables signal and
78 // background. dE and Mes were chosen as discriminating variables to
79 // determine N1 and N2, through a maximum Likelihood fit, and thus the
80 // sPlots for the control variable F, unknown to the fit, was constructed.
81
82 TCanvas *myc = new TCanvas("myc", "sPlots of Mes and F signal and background", 800, 600);
83 myc->SetFillColor(40);
84
85 TPaveText *pt = new TPaveText(0.02, 0.85, 0.98, 0.98);
86 pt->SetFillColor(18);
87 pt->SetTextFont(20);
88 pt->SetTextColor(4);
89 pt->AddText("sPlots of Mes and F signal and background,");
90 pt->AddText("obtained by the tutorial TestSPlot.C on BABAR MC "
91 "data (sPlot_toyMC.fit)");
92 TText *t3 = pt->AddText("M. Pivk and F. R. Le Diberder, Nucl.Inst.Meth.A, physics/0402083");
93 t3->SetTextColor(1);
94 t3->SetTextFont(30);
95 pt->Draw();
96
97 TPad *pad1 = new TPad("pad1", "Mes signal", 0.02, 0.43, 0.48, 0.83, 33);
98 TPad *pad2 = new TPad("pad2", "Mes background", 0.5, 0.43, 0.98, 0.83, 33);
99 TPad *pad3 = new TPad("pad3", "F signal", 0.02, 0.02, 0.48, 0.41, 33);
100 TPad *pad4 = new TPad("pad4", "F background", 0.5, 0.02, 0.98, 0.41, 33);
101 pad1->Draw();
102 pad2->Draw();
103 pad3->Draw();
104 pad4->Draw();
105
106 pad1->cd();
107 pad1->SetGrid();
108 TH1D *sweight00 = splot->GetSWeightsHist(-1, 0, 0);
109 sweight00->SetTitle("Mes signal");
110 sweight00->SetStats(kFALSE);
111 sweight00->Draw("e");
112 sweight00->SetMarkerStyle(21);
113 sweight00->SetMarkerSize(0.7);
114 sweight00->SetMarkerColor(2);
115 sweight00->SetLineColor(2);
116 sweight00->GetXaxis()->SetLabelSize(0.05);
117 sweight00->GetYaxis()->SetLabelSize(0.06);
118 sweight00->GetXaxis()->SetLabelOffset(0.02);
119
120 pad2->cd();
121 pad2->SetGrid();
122 TH1D *sweight10 = splot->GetSWeightsHist(-1, 1, 0);
123 sweight10->SetTitle("Mes background");
124 sweight10->SetStats(kFALSE);
125 sweight10->Draw("e");
126 sweight10->SetMarkerStyle(21);
127 sweight10->SetMarkerSize(0.7);
128 sweight10->SetMarkerColor(2);
129 sweight10->SetLineColor(2);
130 sweight10->GetXaxis()->SetLabelSize(0.05);
131 sweight10->GetYaxis()->SetLabelSize(0.06);
132 sweight10->GetXaxis()->SetLabelOffset(0.02);
133
134 pad3->cd();
135 pad3->SetGrid();
136 TH1D *sweight02 = splot->GetSWeightsHist(-1, 0, 2);
137 sweight02->SetTitle("F signal");
138 sweight02->SetStats(kFALSE);
139 sweight02->Draw("e");
140 sweight02->SetMarkerStyle(21);
141 sweight02->SetMarkerSize(0.7);
142 sweight02->SetMarkerColor(2);
143 sweight02->SetLineColor(2);
144 sweight02->GetXaxis()->SetLabelSize(0.06);
145 sweight02->GetYaxis()->SetLabelSize(0.06);
146 sweight02->GetXaxis()->SetLabelOffset(0.01);
147
148 pad4->cd();
149 pad4->SetGrid();
150 TH1D *sweight12 = splot->GetSWeightsHist(-1, 1, 2);
151 sweight12->SetTitle("F background");
152 sweight12->SetStats(kFALSE);
153 sweight12->Draw("e");
154 sweight12->SetMarkerStyle(21);
155 sweight12->SetMarkerSize(0.7);
156 sweight12->SetMarkerColor(2);
157 sweight12->SetLineColor(2);
158 sweight12->GetXaxis()->SetLabelSize(0.06);
159 sweight12->GetYaxis()->SetLabelSize(0.06);
160 sweight02->GetXaxis()->SetLabelOffset(0.01);
161 myc->cd();
162}
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition TAttText.h:46
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:48
The Canvas class.
Definition TCanvas.h:23
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:925
The most important graphics class in the ROOT system.
Definition TPad.h:28
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
virtual TText * AddText(Double_t x1, Double_t y1, const char *label)
Add a new Text line to this pavetext at given coordinates.
void Draw(Option_t *option="") override
Draw this pavetext with its current attributes.
<div class="legacybox"><h2>Legacy Code</h2> TSPlot is a legacy interface: there will be no bug fixes ...
Definition TSPlot.h:21
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
virtual const char * UnixPathName(const char *unixpathname)
Convert from a local pathname to a Unix pathname.
Definition TSystem.cxx:1075
Base class for several text objects.
Definition TText.h:22
A TTree represents a columnar dataset.
Definition TTree.h:84
TPaveText * pt