Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:04

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