Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:51

0001 //gcc `root-config --cflags --libs` CalculateElectronHLTEff.cpp -o CalculateElectronHLTEff.exe
0002 
0003 //gcc -pthread -m32 -I/afs/cern.ch/cms/sw/slc4_ia32_gcc345/lcg/root/5.14.00g-CMS18b//include -L/afs/cern.ch/cms/sw/slc4_ia32_gcc345/lcg/root/5.14.00g-CMS18b//lib -lCore -lCint -lHist -lGraf -lGraf3d -lGpad -lTree -lRint -lPostscript -lMatrix -lPhysics -pthread -lm -ldl -rdynamic CalculateElectronHLTEff.cpp -o CalculateElectronHLTEff.exe
0004 
0005 #include <iostream>
0006 #include <string.h>
0007 #include <iomanip>
0008 #include<fstream>
0009 #include <math.h>
0010 
0011 #include "TFile.h"
0012 #include "TH1F.h"
0013 
0014 using namespace std;
0015 int main(int argc, char* argv[])
0016 {
0017   if(argc < 3 ){
0018     cout<<"Usage: ./CalculateEgammaHLTEff.exe histofile1.root histofile2.root"<<endl;
0019     return -1;
0020   }
0021   //cout<<"Ciao"<<endl;
0022   TFile file1(argv[1]);
0023   TFile file2(argv[2]);
0024   ofstream outfile ("test.tex");
0025   outfile<<"\\documentclass{article}"<<endl;
0026   outfile<<"\\usepackage{times}"<<endl;
0027   outfile<<"\\usepackage{color}"<<endl;
0028   outfile<<"\\usepackage{amssymb}"<<endl;
0029   outfile<<"\\begin{document}"<<endl;
0030   const int npaths = 3;
0031   string paths[npaths]={"singleElectron","singleElectronRelaxed","singleElectronLargeWindow"};
0032   string Histo_paths[npaths];
0033   for(int i=0;i<npaths;i++){Histo_paths[i]=paths[i]+"DQM/total eff";}
0034   //Single electron path
0035   TH1F* h1=0, *h2=0;
0036   TH1F* hf1=0, *hf2=0;
0037 
0038   for(int path=0; path<npaths; path++){
0039    
0040     cout<<"@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@"<<endl;
0041     cout<<"@@@@@@@@@@@@              "<< paths[path]<<"             \t    @@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@"<<endl;
0042     cout<<"@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@"<<endl;
0043     hf1 = (TH1F*)file1.Get(Histo_paths[path].c_str());
0044     hf2 = (TH1F*)file2.Get(Histo_paths[path].c_str());
0045     if(hf1 != 0 && hf2 !=0 ){
0046       int Maxbin = hf1->GetNbinsX();
0047       if (hf2->GetNbinsX() != Maxbin){
0048     cout<<"AAAA histos for path "<< paths[path] << " have different number of bins"<<endl;
0049     continue;
0050       }
0051       string nm1=paths[path]+"_1";
0052       string nm2=paths[path]+"_2";
0053       h1= new TH1F(nm1.c_str(),nm1.c_str(),Maxbin,0,Maxbin);
0054       h2= new TH1F(nm2.c_str(),nm2.c_str(),Maxbin,0,Maxbin);
0055       
0056       for(int bin=0;bin < Maxbin-1; bin++){
0057     h1->SetBinContent(bin+2,hf1->GetBinContent(bin));
0058     h1->GetXaxis()->SetBinLabel(bin+2,hf1->GetXaxis()->GetBinLabel(bin));
0059     h2->SetBinContent(bin+2,hf2->GetBinContent(bin));
0060     h2->GetXaxis()->SetBinLabel(bin+2,hf2->GetXaxis()->GetBinLabel(bin));
0061     
0062       }
0063       h1->SetBinContent(1,hf1->GetBinContent(Maxbin-1));
0064       h1->SetBinContent(2,hf1->GetBinContent(Maxbin));
0065       h1->GetXaxis()->SetBinLabel(1,"Total");
0066       h1->GetXaxis()->SetBinLabel(2,"Gen");
0067 
0068       h2->SetBinContent(1,hf2->GetBinContent(Maxbin-1));
0069       h2->SetBinContent(2,hf2->GetBinContent(Maxbin));
0070       h2->GetXaxis()->SetBinLabel(1,"Total");
0071       h2->GetXaxis()->SetBinLabel(2,"Gen");
0072 
0073       outfile<<"\\begin{tabular}{ | c || c | c | }"<<endl;
0074       outfile<<"\\hline"<<endl;
0075       outfile<<"\\multicolumn{3}{|c|}{"<<paths[path]<<"}\\\\"<<endl;
0076       outfile<<"\\hline"<<endl;
0077       outfile<<" module & 2\\_1\\_0\\_pre6 & 2\\_1\\_4 \\\\"<<endl;
0078       outfile<<"\\hline"<<endl;
0079       outfile<<" Total Events & "<< h1->GetBinContent(1) <<" & " <<h2->GetBinContent(1)<<" \\\\"<<endl;
0080       outfile<<"\\hline"<<endl;
0081       
0082       cout<<"Total events: "<<h1->GetBinContent(1)<<"  and  "<<h2->GetBinContent(1)<<endl;
0083       for(int bin=2; bin < Maxbin+1; bin++){
0084     //name=h->GetXaxis()->GetBinLabel(bin);
0085     string lab1(h1->GetXaxis()->GetBinLabel(bin));
0086     string lab2(h2->GetXaxis()->GetBinLabel(bin));
0087     if(lab1 != lab2 ){
0088       cout<<"AAAA histos for path "<< paths[path] << " have different labes for bin "<<bin<<" : "<<lab1<<" || "<<lab2<<endl;
0089       //  break;
0090     }
0091     float den1 = h1->GetBinContent(bin-1);
0092     float den2 = h2->GetBinContent(bin-1);
0093     if(den1 >0 && den2 >0 ){
0094       float eff1=h1->GetBinContent(bin)/den1;
0095       float sigma1 = sqrt(eff1*(1-eff1)/den1);
0096 
0097       float eff2=h2->GetBinContent(bin)/den2;
0098       float sigma2 = sqrt(eff2*(1-eff2)/den2);
0099       
0100       float thr = 0.015;//could be tuned on sigma1,2
0101       if(fabs(eff1-eff2)<thr){
0102         cout<<lab1<<" Eff1: "<<fixed<<setprecision(2)<<eff1*100<<" +/- "<<sigma1*100<<" %"<<" ; Eff2: "<<eff2*100<<" +/- "<<sigma2*100<<" %"<<endl; 
0103         outfile<<lab1<< " & "<<fixed<<setprecision(2)<<  eff1*100<<" $\\pm$ "<< sigma1*100 <<" & "<<eff2*100 <<" $\\pm$ "<< sigma2*100<<" \\\\"<<endl;
0104       }
0105       else{
0106         cout<<"AAAABBBB "<<lab1<<" Eff1: "<<fixed<<setprecision(2)<<eff1*100<<" +/- "<<sigma1*100<<" %"<<" ; Eff2: "<<eff2*100<<" +/- "<<sigma2*100<<" %"<<endl;    
0107         outfile<<"\\textcolor{red}{\\textbf{"<<lab1<< "}} & \\textcolor{red}{\\textbf{"<<fixed<<setprecision(2)<<  eff1*100<<" $\\pm$ "<< sigma1*100 <<"}} & \\textcolor{red}{\\textbf{"<<eff2*100 <<" $\\pm$ "<< sigma2*100<<"}} \\\\"<<endl;
0108       }
0109     }
0110     else{cout<<"AAAA no event found in: "<< h1->GetXaxis()->GetBinLabel(bin-1)<<endl;}
0111       
0112 
0113       }
0114     
0115       outfile<<"\\hline"<<endl;
0116       float L1acc1=h1->GetBinContent(3);
0117       float L1acc2=h2->GetBinContent(3);
0118       if(L1acc1 >0 && L1acc2>0 ){
0119 
0120     float eff1=h1->GetBinContent(Maxbin)/L1acc1;
0121     float sigma1 = sqrt(eff1*(1-eff1)/L1acc1);
0122 
0123     float eff2=h2->GetBinContent(Maxbin)/L1acc2;
0124     float sigma2 = sqrt(eff2*(1-eff2)/L1acc2);
0125      float thr = 0.015;//could be tuned on sigma1,2
0126       if(fabs(eff1-eff2)<thr){
0127         cout<<"HLT/L1: "<<fixed<<setprecision(2)<< eff1*100 <<" +/- "<<sigma1*100<<" %"<< " |||| " << eff2*100 <<" +/- "<<sigma2*100<<" %"<<endl;
0128         outfile<<"HLT/L1  & "<<fixed<<setprecision(2)<<  eff1*100<<" $\\pm$ "<< sigma1*100 <<" & "<<eff2*100 <<" $\\pm$ "<< sigma2*100<<" \\\\"<<endl;
0129       }
0130       else{
0131         cout<<"AAAABBBB HLT/L1: "<<fixed<<setprecision(2)<< eff1*100 <<" +/- "<<sigma1*100<<" %"<< " |||| " << eff2*100 <<" +/- "<<sigma2*100<<" %"<<endl;
0132          outfile<<"\\textcolor{red}{\\textbf{HLT/L1}} & \\textcolor{red}{\\textbf{"<<fixed<<setprecision(2)<<  eff1*100<<" $\\pm$ "<< sigma1*100 <<"}} & \\textcolor{red}{\\textbf{"<<eff2*100 <<" $\\pm$ "<< sigma2*100<<"}} \\\\"<<endl;
0133       }
0134       }
0135       else{
0136     cout<<"AAAA No L1 accepted"<<endl;
0137       }
0138       outfile<<"\\hline"<<endl;
0139       float Genacc1=h1->GetBinContent(2);
0140       float Genacc2=h2->GetBinContent(2);
0141       if(Genacc1 >0 && Genacc2>0 ){
0142 
0143     float eff1=h1->GetBinContent(Maxbin)/Genacc1;
0144     float sigma1 = sqrt(eff1*(1-eff1)/Genacc1);
0145 
0146     float eff2=h2->GetBinContent(Maxbin)/Genacc2;
0147     float sigma2 = sqrt(eff2*(1-eff2)/Genacc2);
0148      float thr = 0.015;//could be tuned on sigma1,2
0149       if(fabs(eff1-eff2)<thr){
0150         cout<<"HLT/Gen: "<<fixed<<setprecision(2)<< eff1*100 <<" +/- "<<sigma1*100<<" %"<< " |||| " << eff2*100 <<" +/- "<<sigma2*100<<" %"<<endl;
0151         outfile<<"HLT/Gen  & "<<fixed<<setprecision(2)<<  eff1*100<<" $\\pm$ "<< sigma1*100 <<" & "<<eff2*100 <<" $\\pm$ "<< sigma2*100<<" \\\\"<<endl;
0152       }
0153       else{
0154         cout<<"AAAABBBB HLT/Gen: "<<fixed<<setprecision(2)<< eff1*100 <<" +/- "<<sigma1*100<<" %"<< " |||| " << eff2*100 <<" +/- "<<sigma2*100<<" %"<<endl;
0155          outfile<<"\\textcolor{red}{\\textbf{HLT/Gen}} & \\textcolor{red}{\\textbf{"<<fixed<<setprecision(2)<<  eff1*100<<" $\\pm$ "<< sigma1*100 <<"}} & \\textcolor{red}{\\textbf{"<<eff2*100 <<" $\\pm$ "<< sigma2*100<<"}} \\\\"<<endl;
0156       }
0157       }
0158       else{
0159     cout<<"AAAA No Gen accepted"<<endl;
0160       }
0161       outfile<<"\\hline"<<endl;
0162       outfile<<"\\end{tabular}"<<endl;
0163       outfile<<"\\newline"<<endl;
0164       outfile<<"\\vspace{ 0.8 cm}"<<endl;
0165        outfile<<"\\newline"<<endl;
0166       cout<<"########################################################################################"<<endl<<endl;
0167     }  
0168   }
0169   outfile<<"\\end{document}"<<endl;
0170   system("latex test.tex");
0171   system("dvips test.dvi -o a.ps");
0172   return 1;
0173 }