CMS 3D CMS Logo

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 TNamed

List of all members.

Public Member Functions

void draw (const char *options="")
 Draw on canvas.
void dumpToFile (const char *RootFileName, const char *options)
 All the objects are written to rootfile.
 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.
 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 print (const char *options="")
 Print the relevant information.
void setTitle (const char *title)
 Set the title of the plot.
void setXaxisTitle (const char *title)
 Set the title of the x axis.
 ~LEPBandPlot ()
 Destructor.

Private Attributes

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

Detailed Description

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

Revision:
1.5
Date:
2009/09/23 19:41:22
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.

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 99 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.

                                          :
    StatisticalPlot(name,title,false){

    // bkg hypothesis line
    m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
    m_b_line_graph->SetLineWidth(2);
    m_b_line_graph->SetLineStyle(2);
    m_b_line_graph->SetFillColor(kWhite);



    // bkg hypothesis band 1 sigma
    m_b_band_graph_1sigma = new TGraphAsymmErrors(n_points,
                                                  x_vals,
                                                  b_vals,
                                                  0,
                                                  0,
                                                  b_down_bars1,
                                                  b_up_bars1);
    m_b_band_graph_1sigma->SetFillColor(kGreen);
    m_b_band_graph_1sigma->SetLineColor(kGreen);
    m_b_band_graph_1sigma->SetMarkerColor(kGreen);

    // bkg hypothesis band 2 sigma
    m_b_band_graph_2sigma = new TGraphAsymmErrors(n_points,
                                                  x_vals,
                                                  b_vals,
                                                  0,
                                                  0,
                                                  b_down_bars2,
                                                  b_up_bars2);
    m_b_band_graph_2sigma->SetFillColor(kYellow);
    m_b_band_graph_2sigma->SetFillColor(kYellow);
    m_b_band_graph_2sigma->SetLineColor(kYellow);
    m_b_band_graph_2sigma->SetMarkerColor(kYellow);
    m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");

    // sig+bkg hypothesis line
    m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
    m_sb_line_graph->SetLineWidth(2);
    m_sb_line_graph->SetLineStyle(4);
    m_sb_line_graph->SetLineColor(kRed);
    m_sb_line_graph->SetFillColor(kWhite);

    // The points of the data
    if (exp_vals!=0){
        m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
        m_data_line_graph->SetLineWidth(2);
        m_data_line_graph->SetFillColor(kWhite);
        }
    else
        m_data_line_graph =0;

    // Line for 0
    m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(),0,
                            m_b_line_graph->GetXaxis()->GetXmax(),0);

    // The legend 

    m_legend = new TLegend(0.75,0.78,0.98,0.98);
    m_legend->SetName("Confidence Levels");
    m_legend->AddEntry(m_b_band_graph_1sigma,"-2lnQ #pm 1#sigma");
    m_legend->AddEntry(m_b_band_graph_2sigma,"-2lnQ #pm 2#sigma");
    m_legend->AddEntry(m_b_line_graph,"-2lnQ_{B}");
    m_legend->AddEntry(m_sb_line_graph,"-2lnQ_{SB}");
    if (m_data_line_graph!=NULL)
        m_legend->AddEntry(m_data_line_graph,"-2lnQ_{Obs}");

    m_legend->SetFillColor(0);
    }
LEPBandPlot::~LEPBandPlot ( )

Member Function Documentation

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

Draw on canvas.

Implements StatisticalPlot.

Definition at line 221 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().

                                          {

    setCanvas(new TCanvas(GetName(),GetTitle()));
    getCanvas()->cd();

    getCanvas()->SetGridx();
    getCanvas()->SetGridy();

    TString opt(options);
    // Bands
    if (opt.Contains("4")==0){
        m_b_band_graph_2sigma->Draw("A3");
        m_b_band_graph_1sigma->Draw("3");
        }
    else{
        m_b_band_graph_2sigma->Draw("A4");
        m_b_band_graph_1sigma->Draw("4");
        }

    // Lines
    if (opt.Contains("4")==0){
        m_b_line_graph->Draw("L");
        m_sb_line_graph->Draw("L");
    }
    else{
        m_b_line_graph->Draw("C");
        m_sb_line_graph->Draw("C");
    }

    if (m_data_line_graph!=NULL)
        m_data_line_graph->Draw("L");

    m_zero_line->Draw("Same");

    // Legend
    m_legend->Draw("Same");

    }
void LEPBandPlot::dumpToFile ( const char *  RootFileName,
const char *  options 
) [virtual]

All the objects are written to rootfile.

Implements StatisticalPlot.

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

                                                                          {

    TFile ofile(RootFileName,options);
    ofile.cd();

    // Bands
    m_b_band_graph_2sigma->Write("bkg_band_2sigma");
    m_b_band_graph_1sigma->Draw("bkg_band_1sigma");

    // Lines
    m_b_line_graph->Draw("bkg_line");
    m_sb_line_graph->Draw("sigbkg_line");
    if (m_data_line_graph)
        m_data_line_graph->Draw("observed_line");

    // Legend
    m_legend->Draw("IamTheLegend");

    ofile.Close();

    }
void LEPBandPlot::print ( const char *  options = "") [virtual]

Print the relevant information.

Implements StatisticalPlot.

Definition at line 286 of file LEPBandPlot.cc.

References gather_cfg::cout.

                                           {
    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 215 of file LEPBandPlot.cc.

References m_b_band_graph_2sigma.

                                           {
    m_b_band_graph_2sigma->SetTitle(title);
    }
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 205 of file LEPBandPlot.cc.

References m_b_band_graph_2sigma.

                                                {
        m_b_band_graph_2sigma->GetXaxis()->SetTitle(title);
        }

Member Data Documentation

The b band 1 sigma.

Definition at line 95 of file LEPBandPlot.h.

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

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().