CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/PhysicsTools/RooStatsCms/src/LEPBandPlot.cc

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: LEPBandPlot.cc,v 1.4 2009/05/15 09:55:59 dpiparo Exp $
00002 // Author: Danilo.Piparo@cern.ch   01/06/2008
00003 
00004 #include "assert.h"
00005 #include "math.h"
00006 
00007 #if (defined (STANDALONE) or defined (__CINT__) )
00008    #include "LEPBandPlot.h"
00009 #else
00010    #include "PhysicsTools/RooStatsCms/interface/LEPBandPlot.h"
00011 #endif
00012 #include "TGraphAsymmErrors.h"
00013 #include "TStyle.h"
00014 #include "TAxis.h"
00015 
00016 
00017 //For Cint
00018 #if (defined (STANDALONE) or defined (__CINT__) )
00019 ClassImp(LEPBandPlot)
00020 #endif
00021 /*----------------------------------------------------------------------------*/
00022 
00024 LEPBandPlot::LEPBandPlot(const char* name,
00025                          const char* title,
00026                          const int n_points,
00027                          double* x_vals,
00028                          double* sb_vals,
00029                          double* b_vals,
00030                          double* b_rms,
00031                          double* exp_vals):
00032     StatisticalPlot(name,title,false){
00033 
00034     // bkg hypothesis line
00035     m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
00036     m_b_line_graph->SetLineWidth(2);
00037     m_b_line_graph->SetLineStyle(2);
00038     m_b_line_graph->SetFillColor(kWhite);
00039 
00040 
00041     // bkg hypothesis band 1 sigma
00042     m_b_band_graph_1sigma = new TGraphErrors(n_points, x_vals, b_vals ,0, b_rms);
00043     m_b_band_graph_1sigma->SetFillColor(kGreen);
00044     m_b_band_graph_1sigma->SetLineColor(kGreen);
00045     m_b_band_graph_1sigma->SetMarkerColor(kGreen);
00046 
00047     // Make the band 2 times wider:
00048     double* b_2rms = new double[n_points];
00049     for (int i=0;i<n_points;++i)
00050         b_2rms[i]=2*b_rms[i];
00051 
00052     // bkg hypothesis band 2 sigma
00053     m_b_band_graph_2sigma = new TGraphErrors(n_points, x_vals, b_vals ,0, b_2rms);
00054     m_b_band_graph_2sigma->SetFillColor(kYellow);
00055     m_b_band_graph_2sigma->SetFillColor(kYellow);
00056     m_b_band_graph_2sigma->SetLineColor(kYellow);
00057     m_b_band_graph_2sigma->SetMarkerColor(kYellow);
00058     m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");
00059 
00060     // sig+bkg hypothesis line
00061     m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
00062     m_sb_line_graph->SetLineWidth(2);
00063     m_sb_line_graph->SetLineStyle(4);
00064     m_sb_line_graph->SetLineColor(kRed);
00065     m_sb_line_graph->SetFillColor(kWhite);
00066 
00067     // The points of the data
00068     if (exp_vals!=0){
00069         m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
00070         m_data_line_graph->SetLineWidth(2);
00071         m_data_line_graph->SetFillColor(kWhite);
00072         }
00073     else
00074         m_data_line_graph =NULL;
00075 
00076 
00077     // Line for 0
00078     m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(),0,
00079                             m_b_line_graph->GetXaxis()->GetXmax(),0);
00080 
00081     // The legend 
00082 
00083     m_legend = new TLegend(0.75,0.78,0.98,0.98);
00084     m_legend->SetName("Confidence Levels");
00085     m_legend->AddEntry(m_b_band_graph_1sigma,"-2lnQ #pm 1#sigma");
00086     m_legend->AddEntry(m_b_band_graph_2sigma,"-2lnQ #pm 2#sigma");
00087     m_legend->AddEntry(m_b_line_graph,"-2lnQ_{B}");
00088     m_legend->AddEntry(m_sb_line_graph,"-2lnQ_{SB}");
00089     if (m_data_line_graph!=NULL)
00090         m_legend->AddEntry(m_data_line_graph,"-2lnQ_{Obs}");
00091 
00092     m_legend->SetFillColor(0);
00093 
00094     delete[] b_2rms;
00095     }
00096 
00097 /*----------------------------------------------------------------------------*/
00098 
00099 LEPBandPlot::LEPBandPlot(const char* name,
00100                          const char* title,
00101                          const int n_points,
00102                          double* x_vals,
00103                          double* sb_vals,
00104                          double* b_vals,
00105                          double* b_up_bars1,
00106                          double* b_down_bars1,
00107                          double* b_up_bars2,
00108                          double* b_down_bars2,
00109                          double* exp_vals):
00110     StatisticalPlot(name,title,false){
00111 
00112     // bkg hypothesis line
00113     m_b_line_graph = new TGraph(n_points, x_vals, b_vals);
00114     m_b_line_graph->SetLineWidth(2);
00115     m_b_line_graph->SetLineStyle(2);
00116     m_b_line_graph->SetFillColor(kWhite);
00117 
00118 
00119 
00120     // bkg hypothesis band 1 sigma
00121     m_b_band_graph_1sigma = new TGraphAsymmErrors(n_points,
00122                                                   x_vals,
00123                                                   b_vals,
00124                                                   0,
00125                                                   0,
00126                                                   b_down_bars1,
00127                                                   b_up_bars1);
00128     m_b_band_graph_1sigma->SetFillColor(kGreen);
00129     m_b_band_graph_1sigma->SetLineColor(kGreen);
00130     m_b_band_graph_1sigma->SetMarkerColor(kGreen);
00131 
00132     // bkg hypothesis band 2 sigma
00133     m_b_band_graph_2sigma = new TGraphAsymmErrors(n_points,
00134                                                   x_vals,
00135                                                   b_vals,
00136                                                   0,
00137                                                   0,
00138                                                   b_down_bars2,
00139                                                   b_up_bars2);
00140     m_b_band_graph_2sigma->SetFillColor(kYellow);
00141     m_b_band_graph_2sigma->SetFillColor(kYellow);
00142     m_b_band_graph_2sigma->SetLineColor(kYellow);
00143     m_b_band_graph_2sigma->SetMarkerColor(kYellow);
00144     m_b_band_graph_2sigma->GetYaxis()->SetTitle("-2lnQ");
00145 
00146     // sig+bkg hypothesis line
00147     m_sb_line_graph = new TGraph(n_points, x_vals, sb_vals);
00148     m_sb_line_graph->SetLineWidth(2);
00149     m_sb_line_graph->SetLineStyle(4);
00150     m_sb_line_graph->SetLineColor(kRed);
00151     m_sb_line_graph->SetFillColor(kWhite);
00152 
00153     // The points of the data
00154     if (exp_vals!=0){
00155         m_data_line_graph = new TGraph(n_points, x_vals, exp_vals);
00156         m_data_line_graph->SetLineWidth(2);
00157         m_data_line_graph->SetFillColor(kWhite);
00158         }
00159     else
00160         m_data_line_graph =0;
00161 
00162     // Line for 0
00163     m_zero_line = new TLine(m_b_line_graph->GetXaxis()->GetXmin(),0,
00164                             m_b_line_graph->GetXaxis()->GetXmax(),0);
00165 
00166     // The legend 
00167 
00168     m_legend = new TLegend(0.75,0.78,0.98,0.98);
00169     m_legend->SetName("Confidence Levels");
00170     m_legend->AddEntry(m_b_band_graph_1sigma,"-2lnQ #pm 1#sigma");
00171     m_legend->AddEntry(m_b_band_graph_2sigma,"-2lnQ #pm 2#sigma");
00172     m_legend->AddEntry(m_b_line_graph,"-2lnQ_{B}");
00173     m_legend->AddEntry(m_sb_line_graph,"-2lnQ_{SB}");
00174     if (m_data_line_graph!=NULL)
00175         m_legend->AddEntry(m_data_line_graph,"-2lnQ_{Obs}");
00176 
00177     m_legend->SetFillColor(0);
00178     }
00179 
00180 /*----------------------------------------------------------------------------*/
00181 
00182 LEPBandPlot::~LEPBandPlot(){
00183 
00184     delete m_b_line_graph;
00185     delete m_sb_line_graph;
00186 
00187     delete m_b_band_graph_1sigma;
00188     delete m_b_band_graph_2sigma;
00189 
00190     if (m_data_line_graph!=NULL)
00191         delete m_data_line_graph;
00192 
00193     delete m_zero_line;
00194 
00195     delete m_legend;
00196 
00197     }
00198 
00199 /*----------------------------------------------------------------------------*/
00200 
00205 void LEPBandPlot::setXaxisTitle(const char* title){
00206         m_b_band_graph_2sigma->GetXaxis()->SetTitle(title);
00207         }
00208 
00209 /*----------------------------------------------------------------------------*/
00210 
00215 void LEPBandPlot::setTitle(const char* title){
00216     m_b_band_graph_2sigma->SetTitle(title);
00217     }
00218 
00219 /*----------------------------------------------------------------------------*/
00220 
00221 void LEPBandPlot::draw (const char* options){
00222 
00223     setCanvas(new TCanvas(GetName(),GetTitle()));
00224     getCanvas()->cd();
00225 
00226     getCanvas()->SetGridx();
00227     getCanvas()->SetGridy();
00228 
00229     TString opt(options);
00230     // Bands
00231     if (opt.Contains("4")==0){
00232         m_b_band_graph_2sigma->Draw("A3");
00233         m_b_band_graph_1sigma->Draw("3");
00234         }
00235     else{
00236         m_b_band_graph_2sigma->Draw("A4");
00237         m_b_band_graph_1sigma->Draw("4");
00238         }
00239 
00240     // Lines
00241     if (opt.Contains("4")==0){
00242         m_b_line_graph->Draw("L");
00243         m_sb_line_graph->Draw("L");
00244     }
00245     else{
00246         m_b_line_graph->Draw("C");
00247         m_sb_line_graph->Draw("C");
00248     }
00249 
00250     if (m_data_line_graph!=NULL)
00251         m_data_line_graph->Draw("L");
00252 
00253     m_zero_line->Draw("Same");
00254 
00255     // Legend
00256     m_legend->Draw("Same");
00257 
00258     }
00259 
00260 /*----------------------------------------------------------------------------*/
00261 
00262 void LEPBandPlot::dumpToFile (const char* RootFileName, const char* options){
00263 
00264     TFile ofile(RootFileName,options);
00265     ofile.cd();
00266 
00267     // Bands
00268     m_b_band_graph_2sigma->Write("bkg_band_2sigma");
00269     m_b_band_graph_1sigma->Draw("bkg_band_1sigma");
00270 
00271     // Lines
00272     m_b_line_graph->Draw("bkg_line");
00273     m_sb_line_graph->Draw("sigbkg_line");
00274     if (m_data_line_graph)
00275         m_data_line_graph->Draw("observed_line");
00276 
00277     // Legend
00278     m_legend->Draw("IamTheLegend");
00279 
00280     ofile.Close();
00281 
00282     }
00283 
00284 /*----------------------------------------------------------------------------*/
00285 
00286 void LEPBandPlot::print (const char* options){
00287     std::cout << "\nLEPBandPlot object " << GetName() << ":\n";
00288     }
00289 
00290 /*----------------------------------------------------------------------------*/
00291 // Automatically converted from the standalone version Wed Apr 15 11:36:34 2009