CMS 3D CMS Logo

LEPBandPlot.cc
Go to the documentation of this file.
1 // Author: Danilo.Piparo@cern.ch 01/06/2008
2 
3 #include <cassert>
4 #include <cmath>
5 
6 #if (defined (STANDALONE) or defined (__CINT__) )
7  #include "LEPBandPlot.h"
8 #else
10 #endif
11 #include "TGraphAsymmErrors.h"
12 #include "TStyle.h"
13 #include "TAxis.h"
14 
15 
16 //For Cint
17 #if (defined (STANDALONE) or defined (__CINT__) )
19 #endif
20 /*----------------------------------------------------------------------------*/
21 
23 LEPBandPlot::LEPBandPlot(const char* name,
24  const char* title,
25  const int n_points,
26  double* x_vals,
27  double* sb_vals,
28  double* b_vals,
29  double* b_rms,
30  double* exp_vals):
31  StatisticalPlot(name,title,false){
32 
33  // bkg hypothesis line
34  m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
35  m_b_line_graph->SetLineWidth(2);
36  m_b_line_graph->SetLineStyle(2);
37  m_b_line_graph->SetFillColor(kWhite);
38 
39 
40  // bkg hypothesis band 1 sigma
41  m_b_band_graph_1sigma = new TGraphErrors(n_points, x_vals, b_vals ,0, b_rms);
42  m_b_band_graph_1sigma->SetFillColor(kGreen);
43  m_b_band_graph_1sigma->SetLineColor(kGreen);
44  m_b_band_graph_1sigma->SetMarkerColor(kGreen);
45 
46  // Make the band 2 times wider:
47  double* b_2rms = new double[n_points];
48  for (int i=0;i<n_points;++i)
49  b_2rms[i]=2*b_rms[i];
50 
51  // bkg hypothesis band 2 sigma
52  m_b_band_graph_2sigma = new TGraphErrors(n_points, x_vals, b_vals ,0, b_2rms);
53  m_b_band_graph_2sigma->SetFillColor(kYellow);
54  m_b_band_graph_2sigma->SetFillColor(kYellow);
55  m_b_band_graph_2sigma->SetLineColor(kYellow);
56  m_b_band_graph_2sigma->SetMarkerColor(kYellow);
57  m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");
58 
59  // sig+bkg hypothesis line
60  m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
61  m_sb_line_graph->SetLineWidth(2);
62  m_sb_line_graph->SetLineStyle(4);
63  m_sb_line_graph->SetLineColor(kRed);
64  m_sb_line_graph->SetFillColor(kWhite);
65 
66  // The points of the data
67  if (exp_vals!=0){
68  m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
69  m_data_line_graph->SetLineWidth(2);
70  m_data_line_graph->SetFillColor(kWhite);
71  }
72  else
73  m_data_line_graph =NULL;
74 
75 
76  // Line for 0
77  m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(),0,
78  m_b_line_graph->GetXaxis()->GetXmax(),0);
79 
80  // The legend
81 
82  m_legend = new TLegend(0.75,0.78,0.98,0.98);
83  m_legend->SetName("Confidence Levels");
84  m_legend->AddEntry(m_b_band_graph_1sigma,"-2lnQ #pm 1#sigma");
85  m_legend->AddEntry(m_b_band_graph_2sigma,"-2lnQ #pm 2#sigma");
86  m_legend->AddEntry(m_b_line_graph,"-2lnQ_{B}");
87  m_legend->AddEntry(m_sb_line_graph,"-2lnQ_{SB}");
88  if (m_data_line_graph!=NULL)
89  m_legend->AddEntry(m_data_line_graph,"-2lnQ_{Obs}");
90 
91  m_legend->SetFillColor(0);
92 
93  delete[] b_2rms;
94  }
95 
96 /*----------------------------------------------------------------------------*/
97 
99  const char* title,
100  const int n_points,
101  double* x_vals,
102  double* sb_vals,
103  double* b_vals,
104  double* b_up_bars1,
105  double* b_down_bars1,
106  double* b_up_bars2,
107  double* b_down_bars2,
108  double* exp_vals):
109  StatisticalPlot(name,title,false){
110 
111  // bkg hypothesis line
112  m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
113  m_b_line_graph->SetLineWidth(2);
114  m_b_line_graph->SetLineStyle(2);
115  m_b_line_graph->SetFillColor(kWhite);
116 
117 
118 
119  // bkg hypothesis band 1 sigma
120  m_b_band_graph_1sigma = new TGraphAsymmErrors(n_points,
121  x_vals,
122  b_vals,
123  0,
124  0,
125  b_down_bars1,
126  b_up_bars1);
127  m_b_band_graph_1sigma->SetFillColor(kGreen);
128  m_b_band_graph_1sigma->SetLineColor(kGreen);
129  m_b_band_graph_1sigma->SetMarkerColor(kGreen);
130 
131  // bkg hypothesis band 2 sigma
132  m_b_band_graph_2sigma = new TGraphAsymmErrors(n_points,
133  x_vals,
134  b_vals,
135  0,
136  0,
137  b_down_bars2,
138  b_up_bars2);
139  m_b_band_graph_2sigma->SetFillColor(kYellow);
140  m_b_band_graph_2sigma->SetFillColor(kYellow);
141  m_b_band_graph_2sigma->SetLineColor(kYellow);
142  m_b_band_graph_2sigma->SetMarkerColor(kYellow);
143  m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");
144 
145  // sig+bkg hypothesis line
146  m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
147  m_sb_line_graph->SetLineWidth(2);
148  m_sb_line_graph->SetLineStyle(4);
149  m_sb_line_graph->SetLineColor(kRed);
150  m_sb_line_graph->SetFillColor(kWhite);
151 
152  // The points of the data
153  if (exp_vals!=0){
154  m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
155  m_data_line_graph->SetLineWidth(2);
156  m_data_line_graph->SetFillColor(kWhite);
157  }
158  else
160 
161  // Line for 0
162  m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(),0,
163  m_b_line_graph->GetXaxis()->GetXmax(),0);
164 
165  // The legend
166 
167  m_legend = new TLegend(0.75,0.78,0.98,0.98);
168  m_legend->SetName("Confidence Levels");
169  m_legend->AddEntry(m_b_band_graph_1sigma,"-2lnQ #pm 1#sigma");
170  m_legend->AddEntry(m_b_band_graph_2sigma,"-2lnQ #pm 2#sigma");
171  m_legend->AddEntry(m_b_line_graph,"-2lnQ_{B}");
172  m_legend->AddEntry(m_sb_line_graph,"-2lnQ_{SB}");
173  if (m_data_line_graph!=NULL)
174  m_legend->AddEntry(m_data_line_graph,"-2lnQ_{Obs}");
175 
176  m_legend->SetFillColor(0);
177  }
178 
179 /*----------------------------------------------------------------------------*/
180 
182 
183  delete m_b_line_graph;
184  delete m_sb_line_graph;
185 
186  delete m_b_band_graph_1sigma;
187  delete m_b_band_graph_2sigma;
188 
189  if (m_data_line_graph!=NULL)
190  delete m_data_line_graph;
191 
192  delete m_zero_line;
193 
194  delete m_legend;
195 
196  }
197 
198 /*----------------------------------------------------------------------------*/
199 
205  m_b_band_graph_2sigma->GetXaxis()->SetTitle(title);
206  }
207 
208 /*----------------------------------------------------------------------------*/
209 
214 void LEPBandPlot::setTitle(const char* title){
215  m_b_band_graph_2sigma->SetTitle(title);
216  }
217 
218 /*----------------------------------------------------------------------------*/
219 
220 void LEPBandPlot::draw (const char* options){
221 
222  setCanvas(new TCanvas(GetName(),GetTitle()));
223  getCanvas()->cd();
224 
225  getCanvas()->SetGridx();
226  getCanvas()->SetGridy();
227 
228  TString opt(options);
229  // Bands
230  if (opt.Contains("4")==0){
231  m_b_band_graph_2sigma->Draw("A3");
232  m_b_band_graph_1sigma->Draw("3");
233  }
234  else{
235  m_b_band_graph_2sigma->Draw("A4");
236  m_b_band_graph_1sigma->Draw("4");
237  }
238 
239  // Lines
240  if (opt.Contains("4")==0){
241  m_b_line_graph->Draw("L");
242  m_sb_line_graph->Draw("L");
243  }
244  else{
245  m_b_line_graph->Draw("C");
246  m_sb_line_graph->Draw("C");
247  }
248 
249  if (m_data_line_graph!=NULL)
250  m_data_line_graph->Draw("L");
251 
252  m_zero_line->Draw("Same");
253 
254  // Legend
255  m_legend->Draw("Same");
256 
257  }
258 
259 /*----------------------------------------------------------------------------*/
260 
261 void LEPBandPlot::dumpToFile (const char* RootFileName, const char* options){
262 
263  TFile ofile(RootFileName,options);
264  ofile.cd();
265 
266  // Bands
267  m_b_band_graph_2sigma->Write("bkg_band_2sigma");
268  m_b_band_graph_1sigma->Draw("bkg_band_1sigma");
269 
270  // Lines
271  m_b_line_graph->Draw("bkg_line");
272  m_sb_line_graph->Draw("sigbkg_line");
273  if (m_data_line_graph)
274  m_data_line_graph->Draw("observed_line");
275 
276  // Legend
277  m_legend->Draw("IamTheLegend");
278 
279  ofile.Close();
280 
281  }
282 
283 /*----------------------------------------------------------------------------*/
284 
285 void LEPBandPlot::print (const char* options){
286  std::cout << "\nLEPBandPlot object " << GetName() << ":\n";
287  }
288 
289 /*----------------------------------------------------------------------------*/
290 // Automatically converted from the standalone version Wed Apr 15 11:36:34 2009
LEPBandPlot(const char *name, const char *title, const int n_points, double *x_vals, double *sb_vals, double *b_vals, double *b_rms, double *exp_vals=0)
Constructor.
void draw(const char *options="")
Draw on canvas.
Definition: LEPBandPlot.cc:220
void setCanvas(TCanvas *new_canvas)
Set the canvas.
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:95
#define NULL
Definition: scimark2.h:8
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:101
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:107
void print(const char *options="")
Print the relevant information.
Definition: LEPBandPlot.cc:285
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
void dumpToFile(const char *RootFileName, const char *options)
All the objects are written to rootfile.
Definition: LEPBandPlot.cc:261
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:89
void setTitle(const char *title)
Set the title of the plot.
Definition: LEPBandPlot.cc:214
StatisticalPlot: the base class for the statistical plots.
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:104
~LEPBandPlot()
Destructor.
Definition: LEPBandPlot.cc:181
void setXaxisTitle(const char *title)
Set the title of the x axis.
Definition: LEPBandPlot.cc:204
ClassImp(LEPBandPlot) LEPBandPlot
Constructor.
Definition: LEPBandPlot.cc:18
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:92
TCanvas * getCanvas()
Get the canvas.
LEPBandPlot: The plot for the CL bands "a la LEP".
Definition: LEPBandPlot.h:41