Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:13

0001 ///////////////////////////////////////////////////////////////////////////////
0002 //
0003 //   plotEffiAll(approve, ratio, savePlot);
0004 //      Plots reconstructionefficiency as a function of p, pt, eta, phi
0005 //
0006 //   approve  (bool)   To include "CMS Preliminary" or not (true)
0007 //   ratio    (bool)   If the ratio with the first input to be plotted (true)
0008 //   savePlot (int)    Saving the plot or not: (no:gif:jpg:pdf -1:0:1:2) (2)
0009 //
0010 //   plotCompareAll(cdir1, cdir2, cvers, cfil1, cfil2, ctype1, ctype2,
0011 //                  postfix, ratio, logy, save, norm);
0012 //       Plots comparison plots of calorimetric quatities from native/vecgeom
0013 //       versions
0014 //
0015 //   cdir1   (std::string) Name of the native directory ("10.7.r06.g4")
0016 //   cdir2   (std::string) Name of the vecgeom directory ("10.7.r06.vg")
0017 //   cvers   (std::string) Name of the input version ("10.7.r06 MinBias")
0018 //   cfil1   (std::string) Name of the first input file ("minbias.root")
0019 //   cfil2   (std::string) Name of the second input file ("minbias.root")
0020 //   ctype1  (std::string) Name of the first type ("Native")
0021 //   ctype2  (std::string) Name of the second type ("VecGeom v1.1.16")
0022 //   postfix (std::string) Tag to be given to the canvas name ("MBR")
0023 //   ratio   (bool)        If the ratio to be plotted (false)
0024 //   logy    (bool)        If the scale is log or linear (true)
0025 //   save    (bool)        If the canvas is to be saved (true)
0026 //   norm    (bool)        If the plots to be normalized (true)
0027 //
0028 ///////////////////////////////////////////////////////////////////////////////
0029 #include "TCanvas.h"
0030 #include "TDirectory.h"
0031 #include "TF1.h"
0032 #include "TFile.h"
0033 #include "TFitResult.h"
0034 #include "TGraph.h"
0035 #include "TGraphAsymmErrors.h"
0036 #include "TH1D.h"
0037 #include "TH2D.h"
0038 #include "THStack.h"
0039 #include "TLegend.h"
0040 #include "TMath.h"
0041 #include "TProfile.h"
0042 #include "TPaveStats.h"
0043 #include "TPaveText.h"
0044 #include "TROOT.h"
0045 #include "TString.h"
0046 #include "TStyle.h"
0047 
0048 #include <iostream>
0049 #include <iomanip>
0050 #include <vector>
0051 #include <string>
0052 
0053 TH1D* getEffi(TFile* file, std::string varname, unsigned int ifl);
0054 TCanvas* plotEffi(int type, bool approve);
0055 TCanvas* plotEffiRatio(int type, bool approve);
0056 void plotEffiAll(bool approve = true, bool ratio = true, int savePlot = 2);
0057 void plotCompare(const char* infile1,
0058                  const char* text1,
0059                  const char* infile2,
0060                  const char* text2,
0061                  int type1 = -1,
0062                  int type2 = -1,
0063                  int type3 = -1,
0064                  bool ratio = false,
0065                  bool logy = true,
0066                  std::string postfix = "",
0067                  bool save = false,
0068                  bool normalize = false);
0069 void plotCompareAll(std::string cdir1 = "11.0.r06.g4",
0070                     std::string cdir2 = "11.0.r06.vg",
0071                     std::string cvers = "11.0.r06 MinBias",
0072                     std::string cfil1 = "minbias.root",
0073                     std::string cfil2 = "minbias.root",
0074                     std::string ctype1 = "Native",
0075                     std::string ctype2 = "VecGeom v1.1.20",
0076                     std::string postfix = "MBR",
0077                     bool ratio = false,
0078                     bool logy = true,
0079                     bool save = true,
0080                     bool norm = true);
0081 
0082 int markStyle[7] = {23, 22, 24, 25, 21, 33, 20};
0083 int colors[7] = {2, 4, 6, 7, 46, 3, 1};
0084 int lineStyle[7] = {1, 2, 3, 4, 1, 2, 3};
0085 
0086 /*
0087 const unsigned int nmodels = 4;
0088 std::string filem[nmodels] = {"pikp/FBE7r00vMixStudyHLT.root",
0089                   "pikp/FBE7p01vMixStudyHLT.root",
0090                   "pikp/FBE7p02vMixStudyHLT.root",
0091                   "pikp/FBE7p03vMixStudyHLT.root"};
0092 std::string typem[nmodels] = {"10.7  FTFP_BERT_EMM",
0093                   "10.7.p01 FTFP_BERT_EMM",
0094                   "10.7.p02 FTFP_BERT_EMM",
0095                   "10.7.p03 FTFP_BERT_EMM"};
0096 */
0097 const unsigned int nmodels = 5;
0098 std::string filem[nmodels] = {"pikp/FBE4p3vMixStudyHLT.root",
0099                               "pikp/FBE7p02vMixStudyHLT.root",
0100                               "pikp/FBE110p01vMixStudyHLT.root",
0101                               "pikp/FBE110r04MixStudyHLT.root",
0102                               "pikp/FBE110r06vMixStudyHLT.root"};
0103 std::string typem[nmodels] = {"10.4.p03 FTFP_BERT_EMM",
0104                               "10.7.p01 FTFP_BERT_EMM",
0105                               "11.0.p01 FTFP_BERT_EMM",
0106                               "11.0.ref04 FTFP_BERT_EMM",
0107                               "11.0.ref06 FTFP_BERT_EMM"};
0108 
0109 TH1D* getEffi(TFile* file, std::string varname, unsigned int ifl) {
0110   char name[100];
0111   sprintf(name, "h_%s_All_0", varname.c_str());
0112   TH1D* hist1 = (TH1D*)file->FindObjectAny(name);
0113   sprintf(name, "h_%s_All_1", varname.c_str());
0114   TH1D* hist2 = (TH1D*)file->FindObjectAny(name);
0115   if (hist1 && hist2) {
0116     sprintf(name, "h_%s_Effy_%d", varname.c_str(), ifl);
0117     int nbins = hist1->GetNbinsX();
0118     double xl = hist1->GetBinLowEdge(1);
0119     double xh = hist1->GetBinLowEdge(nbins) + hist1->GetBinWidth(nbins);
0120     TH1D* hist = new TH1D(name, hist1->GetTitle(), nbins, xl, xh);
0121     for (int i = 1; i < nbins; ++i) {
0122       double den = hist1->GetBinContent(i);
0123       double val = (den > 0) ? (hist2->GetBinContent(i)) / den : 0;
0124       double err = (den > 0) ? (hist1->GetBinError(i)) * (val / den) : 0;
0125       hist->SetBinContent(i, val);
0126       hist->SetBinError(i, err);
0127     }
0128     return hist;
0129   } else {
0130     return 0;
0131   }
0132 }
0133 
0134 TCanvas* plotEffi(int type, bool approve) {
0135   std::string varnam[4] = {"pt", "p", "eta", "phi"};
0136   std::string xvtitl[4] = {"p_{T} (GeV)", "p (GeV)", "#eta", "#phi"};
0137   bool irng[4] = {true, true, true, false};
0138   double xlowr[4] = {0.0, 0.0, -2.2, -3.1415926};
0139   double xtopr[4] = {20.0, 20.0, 2.2, 3.1415926};
0140 
0141   TCanvas* c(0);
0142   if (type < 0 || type > 3)
0143     type = 0;
0144   TObjArray histArr;
0145   for (unsigned k = 0; k < nmodels; ++k) {
0146     TFile* file = TFile::Open(filem[k].c_str());
0147     TH1D* hist = getEffi(file, varnam[type], k);
0148     if (hist) {
0149       hist->GetXaxis()->SetTitle(xvtitl[type].c_str());
0150       hist->GetYaxis()->SetTitle("Efficiency");
0151       if (irng[type])
0152         hist->GetXaxis()->SetRangeUser(xlowr[type], xtopr[type]);
0153       histArr.AddLast(hist);
0154     }
0155   }
0156   if (histArr.GetEntries() > 0) {
0157     gStyle->SetCanvasBorderMode(0);
0158     gStyle->SetCanvasColor(kWhite);
0159     gStyle->SetPadColor(kWhite);
0160     gStyle->SetFillColor(kWhite);
0161     gStyle->SetOptTitle(kFALSE);
0162     gStyle->SetPadBorderMode(0);
0163     gStyle->SetCanvasBorderMode(0);
0164     gStyle->SetOptStat(0);
0165 
0166     char cname[50];
0167     sprintf(cname, "c_%sEff", varnam[type].c_str());
0168     c = new TCanvas(cname, cname, 500, 500);
0169     gPad->SetTopMargin(0.10);
0170     gPad->SetBottomMargin(0.10);
0171     gPad->SetLeftMargin(0.15);
0172     gPad->SetRightMargin(0.025);
0173 
0174     TLegend* legend = new TLegend(0.30, 0.15, 0.975, 0.30);
0175     TPaveText* text1 = new TPaveText(0.05, 0.94, 0.35, 0.99, "brNDC");
0176     legend->SetBorderSize(1);
0177     legend->SetFillColor(kWhite);
0178     char texts[200];
0179     sprintf(texts, "CMS Preliminary");
0180     text1->AddText(texts);
0181     THStack* Hs = new THStack("hs2", " ");
0182     for (int i = 0; i < histArr.GetEntries(); i++) {
0183       TH1D* h = (TH1D*)histArr[i];
0184       h->SetLineColor(colors[i]);
0185       h->SetLineWidth(2);
0186       h->SetMarkerSize(markStyle[i]);
0187       Hs->Add(h, "hist sames");
0188       legend->AddEntry(h, typem[i].c_str(), "l");
0189     }
0190     Hs->Draw("nostack");
0191     c->Update();
0192     Hs->GetHistogram()->GetXaxis()->SetTitle(xvtitl[type].c_str());
0193     Hs->GetHistogram()->GetXaxis()->SetLabelSize(0.035);
0194     Hs->GetHistogram()->GetYaxis()->SetTitleOffset(1.6);
0195     Hs->GetHistogram()->GetYaxis()->SetTitle("Track Reconstruction Efficiency");
0196     Hs->GetHistogram()->GetYaxis()->SetRangeUser(0.0, 1.2);
0197     if (irng[type])
0198       Hs->GetHistogram()->GetXaxis()->SetRangeUser(xlowr[type], xtopr[type]);
0199     c->Modified();
0200     c->Update();
0201     legend->Draw("");
0202     if (approve)
0203       text1->Draw("same");
0204     c->Modified();
0205     c->Update();
0206   }
0207   return c;
0208 }
0209 
0210 TCanvas* plotEffiRatio(int type, bool approve) {
0211   std::string varnam[4] = {"pt", "p", "eta", "phi"};
0212   std::string xvtitl[4] = {"p_{T} (GeV)", "p (GeV)", "#eta", "#phi"};
0213   bool irng[4] = {true, true, true, false};
0214   double xlowr[4] = {0.0, 0.0, -2.2, -3.1415926};
0215   double xtopr[4] = {20.0, 20.0, 2.2, 3.1415926};
0216 
0217   TCanvas* c(0);
0218   if (type < 0 || type > 3)
0219     type = 0;
0220   TObjArray histArr;
0221   TH1D* hist0(0);
0222   for (unsigned k = 0; k < nmodels; ++k) {
0223     TFile* file = TFile::Open(filem[k].c_str());
0224     TH1D* hist = getEffi(file, varnam[type], k);
0225     if (hist) {
0226       if (hist0 == nullptr) {
0227         hist0 = hist;
0228       } else {
0229         int nbin = hist->GetNbinsX();
0230         int npt1(0), npt2(0);
0231         double sumNum(0), sumDen(0);
0232         for (int i = 1; i < nbin; ++i) {
0233           double val1 = hist->GetBinContent(i);
0234           double err1 = (val1 > 0) ? ((hist->GetBinError(i)) / val1) : 0;
0235           double val2 = hist0->GetBinContent(i);
0236           double err2 = (val2 > 0) ? ((hist0->GetBinError(i)) / val2) : 0;
0237           double ratio = (val2 > 0) ? (val1 / val2) : 0;
0238           double drat = ratio * sqrt(err1 * err1 + err2 * err2);
0239           if ((((hist->GetBinLowEdge(i)) >= xlowr[type]) &&
0240                ((hist->GetBinLowEdge(i) + hist->GetBinWidth(i)) <= xtopr[type])) ||
0241               (!irng[type])) {
0242             ++npt1;
0243             if (val2 > 0) {
0244               double temp1 = (ratio > 1.0) ? 1.0 / ratio : ratio;
0245               double temp2 = (ratio > 1.0) ? drat / (ratio * ratio) : drat;
0246               sumNum += (fabs(1 - temp1) / (temp2 * temp2));
0247               sumDen += (1.0 / (temp2 * temp2));
0248               ++npt2;
0249             }
0250           }
0251           hist->SetBinContent(i, ratio);
0252           hist->SetBinError(i, drat);
0253         }
0254         sumNum = (sumDen > 0) ? (sumNum / sumDen) : 0;
0255         sumDen = (sumDen > 0) ? 1.0 / sqrt(sumDen) : 0;
0256         std::cout << "Get Ratio of mean for " << npt2 << ":" << npt1 << ":" << nbin << " points: Mean " << sumNum
0257                   << " +- " << sumDen << " from " << filem[k] << std::endl;
0258         hist->GetXaxis()->SetTitle(xvtitl[type].c_str());
0259         hist->GetYaxis()->SetTitle("Efficiency Ratio");
0260         hist->GetYaxis()->SetRangeUser(0.8, 1.2);
0261         if (irng[type])
0262           hist->GetXaxis()->SetRangeUser(xlowr[type], xtopr[type]);
0263         histArr.AddLast(hist);
0264       }
0265     }
0266   }
0267 
0268   if (histArr.GetEntries() > 0) {
0269     gStyle->SetCanvasBorderMode(0);
0270     gStyle->SetCanvasColor(kWhite);
0271     gStyle->SetPadColor(kWhite);
0272     gStyle->SetFillColor(kWhite);
0273     gStyle->SetOptTitle(kFALSE);
0274     gStyle->SetPadBorderMode(0);
0275     gStyle->SetCanvasBorderMode(0);
0276     gStyle->SetOptStat(0);
0277 
0278     char cname[50];
0279     sprintf(cname, "c_%sEffRat", varnam[type].c_str());
0280     c = new TCanvas(cname, cname, 500, 500);
0281     gPad->SetTopMargin(0.10);
0282     gPad->SetBottomMargin(0.10);
0283     gPad->SetLeftMargin(0.15);
0284     gPad->SetRightMargin(0.025);
0285 
0286     TLegend* legend = new TLegend(0.30, 0.15, 0.975, 0.30);
0287     TPaveText* text1 = new TPaveText(0.05, 0.94, 0.35, 0.99, "brNDC");
0288     legend->SetBorderSize(1);
0289     legend->SetFillColor(kWhite);
0290     char texts[200];
0291     sprintf(texts, "CMS Preliminary");
0292     text1->AddText(texts);
0293     THStack* Hs = new THStack("hs2", " ");
0294     for (int i = 0; i < histArr.GetEntries(); i++) {
0295       TH1D* h = (TH1D*)histArr[i];
0296       h->SetLineColor(colors[i + 1]);
0297       h->SetLineWidth(2);
0298       h->SetMarkerSize(markStyle[i + 1]);
0299       Hs->Add(h, "hist sames");
0300       legend->AddEntry(h, typem[i + 1].c_str(), "l");
0301     }
0302     Hs->Draw("nostack");
0303     c->Update();
0304     Hs->GetHistogram()->GetXaxis()->SetTitle(xvtitl[type].c_str());
0305     Hs->GetHistogram()->GetXaxis()->SetLabelSize(0.035);
0306     Hs->GetHistogram()->GetYaxis()->SetTitleOffset(1.6);
0307     Hs->GetHistogram()->GetYaxis()->SetTitle("Track Reconstruction Efficiency Ratio");
0308     Hs->GetHistogram()->GetYaxis()->SetRangeUser(0.8, 1.2);
0309     if (irng[type])
0310       Hs->GetHistogram()->GetXaxis()->SetRangeUser(xlowr[type], xtopr[type]);
0311     c->Modified();
0312     c->Update();
0313     legend->Draw("");
0314     if (approve)
0315       text1->Draw("same");
0316     c->Modified();
0317     c->Update();
0318   }
0319   return c;
0320 }
0321 
0322 void plotEffiAll(bool approve, bool ratio, int savePlot) {
0323   for (int var = 0; var <= 4; ++var) {
0324     TCanvas* c = (ratio) ? plotEffiRatio(var, approve) : plotEffi(var, approve);
0325     if (c != 0 && savePlot >= 0 && savePlot < 3) {
0326       std::string ext[3] = {"eps", "gif", "pdf"};
0327       char name[200];
0328       sprintf(name, "%s.%s", c->GetName(), ext[savePlot].c_str());
0329       c->Print(name);
0330     }
0331   }
0332 }
0333 
0334 void plotCompare(const char* infile1,
0335                  const char* text1,
0336                  const char* infile2,
0337                  const char* text2,
0338                  int type1,
0339                  int type2,
0340                  int type3,
0341                  bool ratio,
0342                  bool logy,
0343                  std::string postfix,
0344                  bool save,
0345                  bool normalize) {
0346   int ndets[4] = {1, 9, 9, 15};
0347   int types[4] = {7, 8, 2, 3};
0348   std::string htype0[7] = {"PtInc", "EneInc", "EtaInc", "PhiInc", "HitLow", "HitHigh", "HitMu"};
0349   int rebin0[7] = {1, 1, 1, 1, 10, 10, 10};
0350   std::string htype1[8] = {"Hit", "Time", "Edep", "Etot", "TimeAll", "EdepEM", "EdepHad", "EtotG"};
0351   int rebin1[8] = {10, 10, 100, 100, 10, 50, 50, 50};
0352   std::string htype2[2] = {"EdepCal", "EdepCalT"};
0353   int rebin2[2] = {50, 50};
0354   std::string htype3[3] = {"HitTk", "TimeTk", "EdepTk"};
0355   int rebin3[3] = {10, 10, 50};
0356   const double ymin(10.0);
0357 
0358   gStyle->SetCanvasBorderMode(0);
0359   gStyle->SetCanvasColor(kWhite);
0360   gStyle->SetPadColor(kWhite);
0361   gStyle->SetFillColor(kWhite);
0362   gStyle->SetOptTitle(0);
0363   gStyle->SetOptFit(0);
0364   if (ratio || normalize)
0365     gStyle->SetOptStat(0);
0366   else
0367     gStyle->SetOptStat(1110);
0368 
0369   int itmin1 = (type1 >= 0) ? type1 : 0;
0370   int itmax1 = (type1 >= 0) ? type1 : 3;
0371   TFile* file1 = new TFile(infile1);
0372   TFile* file2 = new TFile(infile2);
0373   std::cout << "File1: " << infile1 << ":" << file1 << " File2: " << infile2 << ":" << file2 << std::endl;
0374   if (file1 != 0 && file2 != 0) {
0375     for (int it1 = itmin1; it1 <= itmax1; ++it1) {
0376       int itmin2 = (type2 >= 0) ? type2 : 0;
0377       int itmax2 = (type2 >= 0) ? type2 : ((type1 == 1) ? 3 : types[it1] - 1);
0378       int itmin3 = (type3 >= 0) ? type3 : 0;
0379       int itmax3 = (type3 >= 0) ? type3 : ndets[it1] - 1;
0380       for (int it2 = itmin2; it2 <= itmax2; ++it2) {
0381         int rebin(1);
0382         if (it1 == 0)
0383           rebin = rebin0[it2];
0384         else if (it1 == 1)
0385           rebin = rebin1[it2];
0386         else if (it1 == 2)
0387           rebin = rebin2[it2];
0388         else if (it1 == 3)
0389           rebin = rebin3[it2];
0390         for (int it3 = itmin3; it3 <= itmax3; ++it3) {
0391           if (type1 == 1 && (it3 == 1 || it3 == 2))
0392             continue;
0393           char name[20], namec[22];
0394           if (it1 == 0)
0395             sprintf(name, "%s", htype0[it2].c_str());
0396           else if (it1 == 1)
0397             sprintf(name, "%s%d", htype1[it2].c_str(), it3);
0398           else if (it1 == 2)
0399             sprintf(name, "%s%d", htype2[it2].c_str(), it3);
0400           else
0401             sprintf(name, "%s%d", htype3[it2].c_str(), it3);
0402           TH1D* hist[2];
0403           hist[0] = (TH1D*)file1->FindObjectAny(name);
0404           hist[1] = (TH1D*)file2->FindObjectAny(name);
0405           std::cout << name << " Hist " << hist[0] << ":" << hist[1] << std::endl;
0406           if (hist[0] != 0 && hist[1] != 0) {
0407             sprintf(namec, "c_%s%s", name, postfix.c_str());
0408             TCanvas* pad = new TCanvas(namec, namec, 500, 500);
0409             pad->SetRightMargin(0.10);
0410             pad->SetTopMargin(0.10);
0411             if (logy)
0412               pad->SetLogy();
0413             double xedge = (ratio) ? 0.88 : 0.64;
0414             TLegend* legend = new TLegend(0.12, 0.79, xedge, 0.89);
0415             legend->SetFillColor(kWhite);
0416             if (ratio) {
0417               int nbin = hist[0]->GetNbinsX();
0418               int nbinX = nbin / rebin;
0419               sprintf(namec, "%sR", name);
0420               TH1D* hist0 = (TH1D*)gROOT->FindObjectAny(namec);
0421               if (hist0)
0422                 hist0->Delete();
0423               double xlow = hist[0]->GetBinLowEdge(1);
0424               double xhigh = hist[0]->GetBinLowEdge(nbin) + hist[0]->GetBinWidth(nbin);
0425               hist0 = new TH1D(namec, hist[0]->GetTitle(), nbinX, xlow, xhigh);
0426               int npt1(0), npt2(0);
0427               double sumNum(0), sumDen(0);
0428               xhigh = hist0->GetBinLowEdge(1) + hist0->GetBinWidth(1);
0429               double fact = (hist[1]->GetEntries()) / (hist[0]->GetEntries());
0430               bool empty(false);
0431               for (int i = 1; i < nbinX; ++i) {
0432                 double val1(0), err1(0), val2(0), err2(0);
0433                 for (int j = 0; j < rebin; ++j) {
0434                   int i1 = (i - 1) * rebin + j + 1;
0435                   val1 += hist[0]->GetBinContent(i1);
0436                   val2 += hist[1]->GetBinContent(i1);
0437                   err1 += ((hist[0]->GetBinError(i1)) * (hist[0]->GetBinError(i1)));
0438                   err2 += ((hist[1]->GetBinError(i1)) * (hist[1]->GetBinError(i1)));
0439                 }
0440                 err1 = (val1 > 0) ? (sqrt(err1) / val1) : 0;
0441                 err2 = (val2 > 0) ? (sqrt(err2) / val2) : 0;
0442                 double ratio = (((val1 > 0) && (val2 > 0)) ? fact * (val1 / val2) : 1);
0443                 double drat = (((val1 > 0) && (val2 > 0)) ? (ratio * sqrt(err1 * err1 + err2 * err2)) : 0);
0444                 ++npt1;
0445                 if ((val1 > ymin) && (val2 > ymin)) {
0446                   double temp1 = (ratio > 1.0) ? 1.0 / ratio : ratio;
0447                   double temp2 = (ratio > 1.0) ? drat / (ratio * ratio) : drat;
0448                   sumNum += (fabs(1 - temp1) / (temp2 * temp2));
0449                   sumDen += (1.0 / (temp2 * temp2));
0450                   ++npt2;
0451                 }
0452                 hist0->SetBinContent(i, ratio);
0453                 hist0->SetBinError(i, drat);
0454                 if (i <= 10)
0455                   xhigh = (hist0->GetBinLowEdge(i) + hist0->GetBinWidth(i));
0456                 else if (val1 > 0 && val2 > 0 && (!empty))
0457                   xhigh = (hist0->GetBinLowEdge(i) + hist0->GetBinWidth(i));
0458                 else if ((val1 <= 0 || val2 <= 0) && i > 10)
0459                   empty = true;
0460               }
0461               sumNum = (sumDen > 0) ? (sumNum / sumDen) : 0;
0462               sumDen = (sumDen > 0) ? 1.0 / sqrt(sumDen) : 0;
0463               std::cout << "Get Ratio of mean for " << hist[0]->GetTitle() << " with " << npt2 << ":" << npt1 << ":"
0464                         << nbinX << " points: Mean " << sumNum << " +- " << sumDen << std::endl;
0465               if (it1 == 0)
0466                 sprintf(namec, "Ratio for %s", htype0[it2].c_str());
0467               else if (it1 == 1)
0468                 sprintf(namec, "Ratio for %s", htype1[it2].c_str());
0469               else if (it1 == 2)
0470                 sprintf(namec, "Ratio for %s", htype2[it2].c_str());
0471               else
0472                 sprintf(namec, "Ratio for %s", htype3[it2].c_str());
0473               hist0->SetMarkerStyle(markStyle[0]);
0474               hist0->SetMarkerColor(colors[0]);
0475               hist0->SetLineStyle(lineStyle[0]);
0476               hist0->SetLineColor(colors[0]);
0477               hist0->GetXaxis()->SetTitle(hist[0]->GetTitle());
0478               hist0->GetXaxis()->SetRangeUser(xlow, xhigh);
0479               hist0->GetYaxis()->SetTitle(namec);
0480               hist0->GetYaxis()->SetRangeUser(0.0, 2.0);
0481               hist0->GetYaxis()->SetTitleOffset(1.3);
0482               sprintf(namec, "%s vs %s", text1, text2);
0483               legend->AddEntry(hist0, namec, "lp");
0484               hist0->Draw();
0485               pad->Update();
0486               TLine* line = new TLine(xlow, 1.0, xhigh, 1.0);
0487               line->SetLineWidth(2);
0488               line->SetLineStyle(2);
0489               line->Draw("same");
0490               TPaveText* text0 = new TPaveText(0.12, 0.12, 0.65, 0.17, "brNDC");
0491               char texts[200];
0492               sprintf(texts, "Mean Deviation = %5.3f #pm %5.3f", sumNum, sumDen);
0493               text0->SetFillColor(kWhite);
0494               text0->AddText(texts);
0495               text0->Draw("same");
0496             } else {
0497               double ymax(0.90), dy(0.08);
0498               for (int ih = 0; ih < 2; ++ih) {
0499                 hist[ih]->GetXaxis()->SetTitle(hist[ih]->GetTitle());
0500                 hist[ih]->SetMarkerStyle(markStyle[ih]);
0501                 hist[ih]->SetMarkerColor(colors[ih]);
0502                 hist[ih]->SetLineStyle(lineStyle[ih]);
0503                 hist[ih]->SetLineColor(colors[ih]);
0504                 hist[ih]->SetLineWidth(2);
0505                 hist[ih]->GetYaxis()->SetTitleOffset(1.20);
0506                 if (rebin > 1)
0507                   hist[ih]->Rebin(rebin);
0508                 if (ih == 0) {
0509                   legend->AddEntry(hist[ih], text1, "lp");
0510                   if (normalize)
0511                     hist[ih]->DrawNormalized("hist");
0512                   else
0513                     hist[ih]->Draw();
0514                 } else {
0515                   legend->AddEntry(hist[ih], text2, "lp");
0516                   if (normalize)
0517                     hist[ih]->DrawNormalized("sames hist");
0518                   else
0519                     hist[ih]->Draw("sames");
0520                 }
0521                 pad->Update();
0522                 TPaveStats* st1 = (TPaveStats*)hist[ih]->GetListOfFunctions()->FindObject("stats");
0523                 if (st1 != NULL) {
0524                   st1->SetLineColor(colors[ih]);
0525                   st1->SetTextColor(colors[ih]);
0526                   st1->SetY1NDC(ymax - dy);
0527                   st1->SetY2NDC(ymax);
0528                   st1->SetX1NDC(0.65);
0529                   st1->SetX2NDC(0.90);
0530                   ymax -= dy;
0531                 }
0532                 pad->Update();
0533               }
0534             }
0535             legend->Draw("same");
0536             pad->Update();
0537             if (save) {
0538               sprintf(name, "%s.pdf", pad->GetName());
0539               pad->Print(name);
0540             }
0541           }
0542         }
0543       }
0544     }
0545   }
0546 }
0547 
0548 void plotCompareAll(std::string cdir1,
0549                     std::string cdir2,
0550                     std::string cvers,
0551                     std::string cfil1,
0552                     std::string cfil2,
0553                     std::string ctype1,
0554                     std::string ctype2,
0555                     std::string postfix,
0556                     bool ratio,
0557                     bool logy,
0558                     bool save,
0559                     bool norm) {
0560   char infile1[200], infile2[200], text1[200], text2[200];
0561   sprintf(infile1, "%s/%s", cdir1.c_str(), cfil1.c_str());
0562   sprintf(infile2, "%s/%s", cdir2.c_str(), cfil2.c_str());
0563   sprintf(text1, "%s (%s)", cvers.c_str(), ctype1.c_str());
0564   sprintf(text2, "%s (%s)", cvers.c_str(), ctype2.c_str());
0565   plotCompare(infile1, text1, infile2, text2, 1, -1, 0, ratio, logy, postfix, save, norm);
0566   plotCompare(infile1, text1, infile2, text2, 1, -1, 3, ratio, logy, postfix, save, norm);
0567   plotCompare(infile1, text1, infile2, text2, 1, -1, 4, ratio, logy, postfix, save, norm);
0568   plotCompare(infile1, text1, infile2, text2, 1, -1, 5, ratio, logy, postfix, save, norm);
0569   plotCompare(infile1, text1, infile2, text2, 1, -1, 6, ratio, logy, postfix, save, norm);
0570   plotCompare(infile1, text1, infile2, text2, 1, -1, 7, ratio, logy, postfix, save, norm);
0571   plotCompare(infile1, text1, infile2, text2, 1, -1, 8, ratio, logy, postfix, save, norm);
0572 }