CMS 3D CMS Logo

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="") override
 Draw on canvas. More...
 
void dumpToFile (const char *RootFileName, const char *options) override
 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=nullptr)
 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=nullptr)
 Constructor. More...
 
void print (const char *options="") override
 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 () override
 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 () override
 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 = nullptr 
)

Constructor.

Definition at line 22 of file LEPBandPlot.cc.

References mps_fire::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, and m_zero_line.

30  : StatisticalPlot(name, title, false) {
31  // bkg hypothesis line
32  m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
33  m_b_line_graph->SetLineWidth(2);
34  m_b_line_graph->SetLineStyle(2);
35  m_b_line_graph->SetFillColor(kWhite);
36 
37  // bkg hypothesis band 1 sigma
38  m_b_band_graph_1sigma = new TGraphErrors(n_points, x_vals, b_vals, nullptr, b_rms);
39  m_b_band_graph_1sigma->SetFillColor(kGreen);
40  m_b_band_graph_1sigma->SetLineColor(kGreen);
41  m_b_band_graph_1sigma->SetMarkerColor(kGreen);
42 
43  // Make the band 2 times wider:
44  double* b_2rms = new double[n_points];
45  for (int i = 0; i < n_points; ++i)
46  b_2rms[i] = 2 * b_rms[i];
47 
48  // bkg hypothesis band 2 sigma
49  m_b_band_graph_2sigma = new TGraphErrors(n_points, x_vals, b_vals, nullptr, b_2rms);
50  m_b_band_graph_2sigma->SetFillColor(kYellow);
51  m_b_band_graph_2sigma->SetFillColor(kYellow);
52  m_b_band_graph_2sigma->SetLineColor(kYellow);
53  m_b_band_graph_2sigma->SetMarkerColor(kYellow);
54  m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");
55 
56  // sig+bkg hypothesis line
57  m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
58  m_sb_line_graph->SetLineWidth(2);
59  m_sb_line_graph->SetLineStyle(4);
60  m_sb_line_graph->SetLineColor(kRed);
61  m_sb_line_graph->SetFillColor(kWhite);
62 
63  // The points of the data
64  if (exp_vals != nullptr) {
65  m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
66  m_data_line_graph->SetLineWidth(2);
67  m_data_line_graph->SetFillColor(kWhite);
68  } else
69  m_data_line_graph = nullptr;
70 
71  // Line for 0
72  m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(), 0, m_b_line_graph->GetXaxis()->GetXmax(), 0);
73 
74  // The legend
75 
76  m_legend = new TLegend(0.75, 0.78, 0.98, 0.98);
77  m_legend->SetName("Confidence Levels");
78  m_legend->AddEntry(m_b_band_graph_1sigma, "-2lnQ #pm 1#sigma");
79  m_legend->AddEntry(m_b_band_graph_2sigma, "-2lnQ #pm 2#sigma");
80  m_legend->AddEntry(m_b_line_graph, "-2lnQ_{B}");
81  m_legend->AddEntry(m_sb_line_graph, "-2lnQ_{SB}");
82  if (m_data_line_graph != nullptr)
83  m_legend->AddEntry(m_data_line_graph, "-2lnQ_{Obs}");
84 
85  m_legend->SetFillColor(0);
86 
87  delete[] b_2rms;
88 }
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:91
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:97
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:103
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:85
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:100
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:88
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 = nullptr 
)

Constructor.

Definition at line 92 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 m_zero_line.

103  : StatisticalPlot(name, title, false) {
104  // bkg hypothesis line
105  m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
106  m_b_line_graph->SetLineWidth(2);
107  m_b_line_graph->SetLineStyle(2);
108  m_b_line_graph->SetFillColor(kWhite);
109 
110  // bkg hypothesis band 1 sigma
111  m_b_band_graph_1sigma = new TGraphAsymmErrors(n_points, x_vals, b_vals, nullptr, nullptr, b_down_bars1, b_up_bars1);
112  m_b_band_graph_1sigma->SetFillColor(kGreen);
113  m_b_band_graph_1sigma->SetLineColor(kGreen);
114  m_b_band_graph_1sigma->SetMarkerColor(kGreen);
115 
116  // bkg hypothesis band 2 sigma
117  m_b_band_graph_2sigma = new TGraphAsymmErrors(n_points, x_vals, b_vals, nullptr, nullptr, b_down_bars2, b_up_bars2);
118  m_b_band_graph_2sigma->SetFillColor(kYellow);
119  m_b_band_graph_2sigma->SetFillColor(kYellow);
120  m_b_band_graph_2sigma->SetLineColor(kYellow);
121  m_b_band_graph_2sigma->SetMarkerColor(kYellow);
122  m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");
123 
124  // sig+bkg hypothesis line
125  m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
126  m_sb_line_graph->SetLineWidth(2);
127  m_sb_line_graph->SetLineStyle(4);
128  m_sb_line_graph->SetLineColor(kRed);
129  m_sb_line_graph->SetFillColor(kWhite);
130 
131  // The points of the data
132  if (exp_vals != nullptr) {
133  m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
134  m_data_line_graph->SetLineWidth(2);
135  m_data_line_graph->SetFillColor(kWhite);
136  } else
137  m_data_line_graph = nullptr;
138 
139  // Line for 0
140  m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(), 0, m_b_line_graph->GetXaxis()->GetXmax(), 0);
141 
142  // The legend
143 
144  m_legend = new TLegend(0.75, 0.78, 0.98, 0.98);
145  m_legend->SetName("Confidence Levels");
146  m_legend->AddEntry(m_b_band_graph_1sigma, "-2lnQ #pm 1#sigma");
147  m_legend->AddEntry(m_b_band_graph_2sigma, "-2lnQ #pm 2#sigma");
148  m_legend->AddEntry(m_b_line_graph, "-2lnQ_{B}");
149  m_legend->AddEntry(m_sb_line_graph, "-2lnQ_{SB}");
150  if (m_data_line_graph != nullptr)
151  m_legend->AddEntry(m_data_line_graph, "-2lnQ_{Obs}");
152 
153  m_legend->SetFillColor(0);
154 }
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:91
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:97
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:103
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:85
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:100
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:88
StatisticalPlot(const char *name, const char *title, bool verbosity=true)
Constructor.
LEPBandPlot::~LEPBandPlot ( )
override

Destructor.

Definition at line 158 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 m_zero_line.

158  {
159  delete m_b_line_graph;
160  delete m_sb_line_graph;
161 
162  delete m_b_band_graph_1sigma;
163  delete m_b_band_graph_2sigma;
164 
165  if (m_data_line_graph != nullptr)
166  delete m_data_line_graph;
167 
168  delete m_zero_line;
169 
170  delete m_legend;
171 }
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:91
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:97
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:103
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:85
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:100
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:88

Member Function Documentation

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

Draw on canvas.

Implements StatisticalPlot.

Definition at line 191 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, runTheMatrix::opt, and StatisticalPlot::setCanvas().

191  {
192  setCanvas(new TCanvas(GetName(), GetTitle()));
193  getCanvas()->cd();
194 
195  getCanvas()->SetGridx();
196  getCanvas()->SetGridy();
197 
198  TString opt(options);
199  // Bands
200  if (opt.Contains("4") == 0) {
201  m_b_band_graph_2sigma->Draw("A3");
202  m_b_band_graph_1sigma->Draw("3");
203  } else {
204  m_b_band_graph_2sigma->Draw("A4");
205  m_b_band_graph_1sigma->Draw("4");
206  }
207 
208  // Lines
209  if (opt.Contains("4") == 0) {
210  m_b_line_graph->Draw("L");
211  m_sb_line_graph->Draw("L");
212  } else {
213  m_b_line_graph->Draw("C");
214  m_sb_line_graph->Draw("C");
215  }
216 
217  if (m_data_line_graph != nullptr)
218  m_data_line_graph->Draw("L");
219 
220  m_zero_line->Draw("Same");
221 
222  // Legend
223  m_legend->Draw("Same");
224 }
void setCanvas(TCanvas *new_canvas)
Set the canvas.
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:91
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:97
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:103
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:85
TLine * m_zero_line
The line at 0.
Definition: LEPBandPlot.h:100
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:88
TCanvas * getCanvas()
Get the canvas.
void LEPBandPlot::dumpToFile ( const char *  RootFileName,
const char *  options 
)
overridevirtual

All the objects are written to rootfile.

Implements StatisticalPlot.

Definition at line 228 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 indexGen::ofile.

228  {
229  TFile ofile(RootFileName, options);
230  ofile.cd();
231 
232  // Bands
233  m_b_band_graph_2sigma->Write("bkg_band_2sigma");
234  m_b_band_graph_1sigma->Draw("bkg_band_1sigma");
235 
236  // Lines
237  m_b_line_graph->Draw("bkg_line");
238  m_sb_line_graph->Draw("sigbkg_line");
239  if (m_data_line_graph)
240  m_data_line_graph->Draw("observed_line");
241 
242  // Legend
243  m_legend->Draw("IamTheLegend");
244 
245  ofile.Close();
246 }
TGraph * m_b_band_graph_1sigma
The b band 1 sigma.
Definition: LEPBandPlot.h:91
TGraph * m_sb_line_graph
The sb line.
Definition: LEPBandPlot.h:97
TLegend * m_legend
The legend.
Definition: LEPBandPlot.h:103
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94
TGraph * m_data_line_graph
The data line.
Definition: LEPBandPlot.h:85
TGraph * m_b_line_graph
The b line.
Definition: LEPBandPlot.h:88
void LEPBandPlot::print ( const char *  options = "")
overridevirtual

Print the relevant information.

Implements StatisticalPlot.

Definition at line 250 of file LEPBandPlot.cc.

References gather_cfg::cout.

250 { std::cout << "\nLEPBandPlot object " << GetName() << ":\n"; }
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 187 of file LEPBandPlot.cc.

References m_b_band_graph_2sigma.

187 { m_b_band_graph_2sigma->SetTitle(title); }
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94
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 179 of file LEPBandPlot.cc.

References m_b_band_graph_2sigma.

179 { m_b_band_graph_2sigma->GetXaxis()->SetTitle(title); }
TGraph * m_b_band_graph_2sigma
The b band 2 sigma.
Definition: LEPBandPlot.h:94

Member Data Documentation

TGraph* LEPBandPlot::m_b_band_graph_1sigma
private

The b band 1 sigma.

Definition at line 91 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 94 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 88 of file LEPBandPlot.h.

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

TGraph* LEPBandPlot::m_data_line_graph
private

The data line.

Definition at line 85 of file LEPBandPlot.h.

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

TLegend* LEPBandPlot::m_legend
private

The legend.

Definition at line 103 of file LEPBandPlot.h.

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

TGraph* LEPBandPlot::m_sb_line_graph
private

The sb line.

Definition at line 97 of file LEPBandPlot.h.

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

TLine* LEPBandPlot::m_zero_line
private

The line at 0.

Definition at line 100 of file LEPBandPlot.h.

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