CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
LEPBandPlot Class Reference

LEPBandPlot: The plot for the CL bands "a la LEP". More...

#include <LEPBandPlot.h>

Inheritance diagram for LEPBandPlot:
StatisticalPlot

Public Member Functions

void draw (const char *options="")
 Draw on canvas. More...
 
void dumpToFile (const char *RootFileName, const char *options)
 All the objects are written to rootfile. More...
 
 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. More...
 
 LEPBandPlot (const char *name, const char *title, const int n_points, double *x_vals, double *sb_vals, double *b_vals, double *b_up_bars1, double *b_down_bars1, double *b_up_bars2, double *b_down_bars2, double *exp_vals=0)
 Constructor. More...
 
void print (const char *options="")
 Print the relevant information. More...
 
void setTitle (const char *title)
 Set the title of the plot. More...
 
void setXaxisTitle (const char *title)
 Set the title of the x axis. More...
 
 ~LEPBandPlot ()
 Destructor. More...
 
- Public Member Functions inherited from StatisticalPlot
void dumpToImage (const char *filename)
 Write an image on disk. More...
 
TCanvas * getCanvas ()
 Get the canvas. More...
 
bool is_verbose ()
 get the verbosity More...
 
void setCanvas (TCanvas *new_canvas)
 Set the canvas. More...
 
void setVerbosity (bool verbosity)
 Set the verbosity. More...
 
 StatisticalPlot (const char *name, const char *title, bool verbosity=true)
 Constructor. More...
 
 ~StatisticalPlot ()
 Destructor. More...
 

Private Attributes

TGraph * m_b_band_graph_1sigma
 The b band 1 sigma. More...
 
TGraph * m_b_band_graph_2sigma
 The b band 2 sigma. More...
 
TGraph * m_b_line_graph
 The b line. More...
 
TGraph * m_data_line_graph
 The data line. More...
 
TLegend * m_legend
 The legend. More...
 
TGraph * m_sb_line_graph
 The sb line. More...
 
TLine * m_zero_line
 The line at 0. More...
 

Detailed Description

LEPBandPlot: The plot for the CL bands "a la LEP".

Revision:
1.4
Date:
2009/05/15 09:55:43
Author
D. Piparo (danilo.piparo<at>cern.ch), G. Schott - Universitaet Karlsruhe

This class allows to produce plots like the ones of the Tevatron HWG and originally introduced by LEP HWG. It is thought as a "enhanced" TGraph. The input to give to obtain the plot are:

LEP_band_plot.png

Definition at line 41 of file LEPBandPlot.h.

Constructor & Destructor Documentation

LEPBandPlot::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.

Definition at line 23 of file LEPBandPlot.cc.

References i, m_b_band_graph_1sigma, m_b_band_graph_2sigma, m_b_line_graph, m_data_line_graph, m_legend, m_sb_line_graph, m_zero_line, and NULL.

30  :
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
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}");
89  m_legend->AddEntry(m_data_line_graph,"-2lnQ_{Obs}");
90 
91  m_legend->SetFillColor(0);
92 
93  delete[] b_2rms;
94  }
int i
Definition: DBlmapReader.cc:9
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
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:89
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:104
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:92
StatisticalPlot(const char *name, const char *title, bool verbosity=true)
Constructor.
LEPBandPlot::LEPBandPlot ( const char *  name,
const char *  title,
const int  n_points,
double *  x_vals,
double *  sb_vals,
double *  b_vals,
double *  b_up_bars1,
double *  b_down_bars1,
double *  b_up_bars2,
double *  b_down_bars2,
double *  exp_vals = 0 
)

Constructor.

Definition at line 98 of file LEPBandPlot.cc.

References m_b_band_graph_1sigma, m_b_band_graph_2sigma, m_b_line_graph, m_data_line_graph, m_legend, m_sb_line_graph, m_zero_line, and NULL.

108  :
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  }
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
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:89
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:104
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:92
StatisticalPlot(const char *name, const char *title, bool verbosity=true)
Constructor.
LEPBandPlot::~LEPBandPlot ( )

Destructor.

Definition at line 181 of file LEPBandPlot.cc.

References m_b_band_graph_1sigma, m_b_band_graph_2sigma, m_b_line_graph, m_data_line_graph, m_legend, m_sb_line_graph, m_zero_line, and NULL.

181  {
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  }
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
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:89
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:104
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:92

Member Function Documentation

void LEPBandPlot::draw ( const char *  options = "")
virtual

Draw on canvas.

Implements StatisticalPlot.

Definition at line 220 of file LEPBandPlot.cc.

References StatisticalPlot::getCanvas(), m_b_band_graph_1sigma, m_b_band_graph_2sigma, m_b_line_graph, m_data_line_graph, m_legend, m_sb_line_graph, m_zero_line, NULL, and StatisticalPlot::setCanvas().

220  {
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  }
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
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:89
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:104
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:92
TCanvas * getCanvas()
Get the canvas.
void LEPBandPlot::dumpToFile ( const char *  RootFileName,
const char *  options 
)
virtual

All the objects are written to rootfile.

Implements StatisticalPlot.

Definition at line 261 of file LEPBandPlot.cc.

References m_b_band_graph_1sigma, m_b_band_graph_2sigma, m_b_line_graph, m_data_line_graph, m_legend, m_sb_line_graph, and checklumidiff::ofile.

261  {
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  }
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:95
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:101
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:107
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:89
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:92
void LEPBandPlot::print ( const char *  options = "")
virtual

Print the relevant information.

Implements StatisticalPlot.

Definition at line 285 of file LEPBandPlot.cc.

References gather_cfg::cout.

285  {
286  std::cout << "\nLEPBandPlot object " << GetName() << ":\n";
287  }
tuple cout
Definition: gather_cfg.py:121
void LEPBandPlot::setTitle ( const char *  title)

Set the title of the plot.

The title is here set only for the plot of the 2sigma band plot. Indeed its axes will be the only one to be drawn.

Definition at line 214 of file LEPBandPlot.cc.

References m_b_band_graph_2sigma.

214  {
215  m_b_band_graph_2sigma->SetTitle(title);
216  }
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98
void LEPBandPlot::setXaxisTitle ( const char *  title)

Set the title of the x axis.

The title of the x axis is here set only for the plot of the 2sigma band plot. Indeed its axes will be the only one to be drawn.

Definition at line 204 of file LEPBandPlot.cc.

References m_b_band_graph_2sigma.

204  {
205  m_b_band_graph_2sigma->GetXaxis()->SetTitle(title);
206  }
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:98

Member Data Documentation

TGraph* LEPBandPlot::m_b_band_graph_1sigma
private

The b band 1 sigma.

Definition at line 95 of file LEPBandPlot.h.

Referenced by draw(), dumpToFile(), LEPBandPlot(), and ~LEPBandPlot().

TGraph* LEPBandPlot::m_b_band_graph_2sigma
private

The b band 2 sigma.

Definition at line 98 of file LEPBandPlot.h.

Referenced by draw(), dumpToFile(), LEPBandPlot(), setTitle(), setXaxisTitle(), and ~LEPBandPlot().

TGraph* LEPBandPlot::m_b_line_graph
private

The b line.

Definition at line 92 of file LEPBandPlot.h.

Referenced by draw(), dumpToFile(), LEPBandPlot(), and ~LEPBandPlot().

TGraph* LEPBandPlot::m_data_line_graph
private

The data line.

Definition at line 89 of file LEPBandPlot.h.

Referenced by draw(), dumpToFile(), LEPBandPlot(), and ~LEPBandPlot().

TLegend* LEPBandPlot::m_legend
private

The legend.

Definition at line 107 of file LEPBandPlot.h.

Referenced by draw(), dumpToFile(), LEPBandPlot(), and ~LEPBandPlot().

TGraph* LEPBandPlot::m_sb_line_graph
private

The sb line.

Definition at line 101 of file LEPBandPlot.h.

Referenced by draw(), dumpToFile(), LEPBandPlot(), and ~LEPBandPlot().

TLine* LEPBandPlot::m_zero_line
private

The line at 0.

Definition at line 104 of file LEPBandPlot.h.

Referenced by draw(), LEPBandPlot(), and ~LEPBandPlot().