Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:58

0001 #ifndef UTIL
0002 #define UTIL
0003 #include <utility>
0004 #endif
0005 
0006 #ifndef MAP
0007 #define MAP
0008 #include <map>
0009 #endif
0010 
0011 #ifndef VECTOR
0012 #define VECTOR
0013 #include <vector>
0014 #endif
0015 
0016 template <class T>
0017 inline std::string to_string (const T& t)
0018 {
0019   std::stringstream ss; 
0020   ss << t;
0021   return ss.str();
0022 }
0023 
0024 
0025 void getEfficiency(float efficientEvents, float allEvents, std::vector<float> &efficiencyResult);
0026 void getType(int iE, int iS, int iR, float & verticalScale);
0027 
0028 void computeEfficiencies(){
0029   //TFile* anaFile = TFile::Open("efficiencies.root", "RECREATE"); // my output file
0030   //cout << tmp[0] << endl;
0031   char *file_name = "cscHists_bigger.root";
0032   TFile *f1=
0033     (TFile*)gROOT->GetListOfFiles()->FindObject(file_name);
0034   if (!f1){
0035     TFile *f1 = new TFile(file_name);
0036   }
0037   // All files in a vector
0038   std::vector< TFile * > DataFiles;
0039   DataFiles.push_back(f1);
0040   
0041   string mytitle;
0042   char *myTitle;
0043   
0044   const Int_t nx = 36;
0045   const Int_t ny = 16;
0046   
0047   char *chambers[nx]  = {"01","02","03","04","05","06","07","08","09","10",
0048              "11","12","13","14","15","16","17","18","19","20",
0049              "21","22","23","24","25","26","27","28","29","30",
0050              "31","32","33","34","35","36"};
0051   char *types[ny] = {"ME-41","ME-32","ME-31","ME-22","ME-21","ME-13","ME-12","ME-11",
0052              "ME+11","ME+12","ME+13","ME+21","ME+22","ME+31","ME+32","ME+41"};
0053   
0054   
0055   TH2F *h_alctEfficiency = new TH2F("h_alctEfficiency","ALCT efficiency; chamber number  ",
0056                     nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0057   TH2F *h_clctEfficiency = new TH2F("h_clctEfficiency","CLCT efficiency; chamber number  ",
0058                     nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0059   TH2F *h_corrlctEfficiency = new TH2F("h_corrlctEfficiency","CorrLCT efficiency; chamber number  ",
0060                        nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0061   TH2F *h_rhEfficiency = new TH2F("h_rhEfficiency","RecHit efficiency; chamber number ",
0062                   nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0063   TH2F *h_simrhEfficiency = new TH2F("h_simrhEfficiency","RecHit efficiency (simhit based); chamber number ",
0064                   nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0065   TH2F *h_segEfficiency = new TH2F("h_segEfficiency","Segment efficiency; chamber number  ",
0066                    nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0067   
0068   TH2F *h_stripEfficiency = new TH2F("h_stripEfficiency","Strip efficiency; chamber number  ",
0069                      nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0070   TH2F *h_wireEfficiency = new TH2F("h_wireEfficiency","Wire group efficiency; chamber number  ",
0071                     nx,0. + 0.5,float(nx) + 0.5, ny,0. + 0.5,float(ny) + 0.5);
0072 
0073   for(int i=0;i<nx;++i){
0074     h_alctEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0075     h_clctEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0076     h_corrlctEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0077     h_rhEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0078     h_simrhEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0079     h_segEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0080     h_stripEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0081     h_wireEfficiency->GetXaxis()->SetBinLabel(i+1,chambers[i]);
0082   }
0083   for(int i=0;i<ny;++i){
0084     h_alctEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0085     h_clctEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0086     h_corrlctEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0087     h_rhEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0088     h_simrhEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0089     h_segEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0090     h_stripEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0091     h_wireEfficiency->GetYaxis()->SetBinLabel(i+1,types[i]);
0092   }
0093   
0094   //TH1F *h_alct_theta_Efficiency = new TH1F("alct_theta_Efficiency",
0095   //                   "ALCT efficiency vs local theta;local theta;efficiency",
0096   //                   );
0097   
0098   TFile* anaFile = TFile::Open("efficiencies.root", "RECREATE"); // my output file
0099   TCanvas *c1 = new TCanvas("c1", "canvas 1",16,30,700,500);
0100   gStyle->SetOptStat(0);
0101   gStyle->SetErrorX(0);
0102   gPad->SetFillColor(0);
0103   gStyle->SetPalette(1);
0104   
0105   TH1F *data_p1;
0106   TH1F *data_p2;
0107   TH1F *data_p3;
0108   TH1F *data_p4;
0109 
0110   TH1F * h_Chi2_ndf;
0111   TH1F * h_hitsInSegm;
0112   TH1F * h_Residual_Segm;
0113 
0114 
0115 
0116   TH1F * hV_All_alct_dydz[ny];
0117   TH1F * hV_Eff_alct_dydz[ny];
0118   TH1F * hV_All_clct_dxdz[ny];
0119   TH1F * hV_Eff_clct_dxdz[ny];
0120 
0121   TH1F * h_All_alct_dydz;
0122   TH1F * h_Eff_alct_dydz;
0123 
0124   TH1F * h_attachment_AllHits;
0125   TH1F * h_attachment_EffHits;
0126 
0127   bool nextLoop[ny];
0128   for(int i=0;i<ny;++i){
0129     nextLoop[i] = false;
0130   }
0131   //std::cout<<" ny = "<<ny<<" size = "<<  nextLoop.size()<<std::endl;
0132 
0133   string dirName;
0134   char * charName;
0135   char * charName_2;
0136   char * charName_3;
0137   char * charName_4;
0138   char * charName_5;
0139   string histo_1, histo_2, histo_3, histo_4, histo_5;
0140   int firstSt_it = 0;
0141   //TH1F *sum_histo;
0142   int nEfficient = 0;
0143   int nAll = 0;
0144   std::vector<float> efficiencyResult(2);
0145   float verticalScale;
0146   int iterations = 0;
0147   for(int iE=1;iE<3;++iE){
0148     for(int iS=1;iS<5;++iS){
0149       dirName = Form("Stations__E%d_S%d",iE,iS);
0150 
0151       histo_1 = dirName + "/" + Form("segmentChi2_ndf_St%d",iS);
0152       histo_2 = dirName + "/" + Form("hitsInSegment_St%d",iS);
0153       histo_3 = dirName + "/" + Form("ResidualSegments_St%d",iS);
0154       charName = histo_1.c_str();
0155       charName_2 = histo_2.c_str();
0156       charName_3 = histo_3.c_str();
0157       data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0158       data_p2 =(TH1F*)(DataFiles[0]->Get(charName_2));
0159       data_p3 =(TH1F*)(DataFiles[0]->Get(charName_3));
0160 
0161       if(!firstSt_it){
0162     h_Chi2_ndf = (TH1F *)data_p1->Clone();
0163     
0164     h_hitsInSegm = (TH1F *)data_p2->Clone();
0165     h_Residual_Segm = (TH1F *)data_p3->Clone();
0166     //sum_histo = (TH1F *)data_p1->Clone();
0167     
0168       }
0169       else{
0170     
0171     h_Chi2_ndf->Add(data_p1);
0172     h_hitsInSegm->Add(data_p2);
0173     h_Residual_Segm->Add(data_p3);
0174     //sum_histo->Add(data_p1);
0175     
0176       }
0177       
0178       ++firstSt_it;
0179       for(int iR=1;iR<4;++iR){
0180     if(1!=iS && iR>2){
0181       continue;
0182     }
0183     else if(2==iR && 4==iS){
0184       continue;
0185     }
0186     getType(iE,iS,iR, verticalScale);
0187     //std::cout<<" iE = "<<iE<< " iS = "<<iS<<" iR = "<<iR<<" verticalScale = "<<verticalScale<<std::endl;
0188     for(int iC=1;iC<37;++iC){
0189       if(1!=iS && 1==iR && iC >18){
0190         continue;
0191       }
0192       // Attachment 
0193       histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/AllSingleHits_Ch%d",iE,iS,iR,iC,iC);
0194       histo_2 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/InefficientSingleHits_Ch%d",iE,iS,iR,iC,iC);
0195       charName = histo_1.c_str();
0196       data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0197       charName_2 = histo_2.c_str();
0198       data_p2 =(TH1F*)(DataFiles[0]->Get(charName_2));
0199       if(!iterations){
0200         h_attachment_AllHits=(TH1F*)data_p1->Clone();
0201         h_attachment_EffHits=(TH1F*)data_p1->Clone();
0202         h_attachment_EffHits->Add(data_p2, -1.);
0203       }
0204       else{
0205         h_attachment_AllHits->Add(data_p1);
0206         h_attachment_EffHits->Add(data_p1);
0207         h_attachment_EffHits->Add(data_p2, -1.);
0208       }
0209       // RecHits
0210       histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/EfficientRechits_good_Ch%d",iE,iS,iR,iC,iC);
0211       charName = histo_1.c_str();
0212       data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0213       std::cout<<" chamber name : "<<charName<<std::endl;
0214       nEfficient = 0;
0215       for(int iL=1;iL<7;++iL){
0216         nEfficient +=data_p1->GetBinContent(iL+1);
0217       }
0218       nAll = 6 * data_p1->GetBinContent(8+1);
0219       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0220       std::cout<<" Rechit eff = "<<efficiencyResult[0]<<" +-"<<efficiencyResult[1]<<std::endl;
0221       h_rhEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0222       h_rhEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]); 
0223           //(simhit based)
0224           histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/Sim_Simhits_Ch%d",iE,iS,iR,iC,iC);
0225           histo_2 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/Sim_Rechits_Ch%d",iE,iS,iR,iC,iC);
0226           charName = histo_1.c_str();
0227           charName_2 = histo_2.c_str();
0228           data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0229           data_p2 =(TH1F*)(DataFiles[0]->Get(charName_2));
0230       std::cout<<" chamber name : "<<charName<<std::endl;
0231           nEfficient = 0;
0232           nAll =0;
0233           for(int iL=1;iL<7;++iL){
0234             nAll += data_p1->GetBinContent(iL+1);
0235             nEfficient += data_p2->GetBinContent(iL+1);
0236           }
0237           getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0238           h_simrhEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0239           h_simrhEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0240 
0241 
0242       // Segments
0243       histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/digiAppearanceCount_Ch%d",iE,iS,iR,iC,iC);
0244       charName = histo_1.c_str();
0245       data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0246       nEfficient = data_p1->GetBinContent(2);
0247       nAll = data_p1->GetBinContent(1) + nEfficient;
0248       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0249       h_segEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0250       h_segEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0251       // LCT
0252       nEfficient = data_p1->GetBinContent(4);
0253       nAll = data_p1->GetBinContent(3) + nEfficient;
0254       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0255       h_alctEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0256       h_alctEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0257       //
0258       nEfficient = data_p1->GetBinContent(6);
0259       nAll = data_p1->GetBinContent(5) + nEfficient;
0260       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0261       h_clctEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0262       h_clctEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0263       //
0264       nEfficient = data_p1->GetBinContent(8);
0265       nAll = data_p1->GetBinContent(7) + nEfficient;
0266       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0267       h_corrlctEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0268       h_corrlctEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0269       //
0270       
0271       // Strips
0272       histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/EfficientStrips_Ch%d",iE,iS,iR,iC,iC);
0273       charName = histo_1.c_str();
0274       data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0275       nEfficient = 0;
0276       
0277       for(int iBin = 2;iBin<8;++iBin){
0278         nEfficient +=  data_p1->GetBinContent(iBin);
0279       }
0280       nAll = 6*data_p1->GetBinContent(10);    
0281       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0282       h_stripEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0283       h_stripEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0284       // Wires
0285       histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/EfficientWireGroups_Ch%d",iE,iS,iR,iC,iC);
0286       charName = histo_1.c_str();
0287       data_p1 =(TH1F*)(DataFiles[0]->Get(charName));
0288       nEfficient = 0;
0289       for(int iBin = 2;iBin<8;++iBin){
0290         nEfficient +=  data_p1->GetBinContent(iBin);
0291       }
0292       nAll = 6*data_p1->GetBinContent(10);
0293       getEfficiency(float(nEfficient), float(nAll), efficiencyResult);
0294       h_wireEfficiency->SetBinContent(iC,int(verticalScale+0.5), efficiencyResult[0]);
0295       h_wireEfficiency->SetBinError(iC,int(verticalScale+0.5), efficiencyResult[1]);
0296       
0297       //ALCT
0298       histo_1 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/EfficientALCT_dydz_Ch%d",iE,iS,iR,iC,iC);
0299       charName_2 = histo_1.c_str();
0300       histo_2 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/InefficientALCT_dydz_Ch%d",iE,iS,iR,iC,iC);
0301       charName_3 = histo_2.c_str();
0302       data_p1 =(TH1F*)(DataFiles[0]->Get(charName_2));
0303       data_p2 =(TH1F*)(DataFiles[0]->Get(charName_3));
0304 
0305       histo_4 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/EfficientCLCT_dxdz_Ch%d",iE,iS,iR,iC,iC);
0306       charName_4 = histo_4.c_str();
0307       histo_5 = Form("Chambers__E%d_S%d_R%d_Chamber_%d/InefficientCLCT_dxdz_Ch%d",iE,iS,iR,iC,iC);
0308       charName_5 = histo_5.c_str();
0309       data_p3 =(TH1F*)(DataFiles[0]->Get(charName_4));
0310       data_p4 =(TH1F*)(DataFiles[0]->Get(charName_5));
0311 
0312       int type = int(verticalScale+0.1);
0313       if(!nextLoop[type]){
0314         nextLoop[type] = true;
0315         hV_All_alct_dydz[int(verticalScale+0.5)-1]=(TH1F*)data_p1->Clone();
0316         hV_All_alct_dydz[int(verticalScale+0.5)-1]->Add(data_p2);
0317         hV_Eff_alct_dydz[int(verticalScale+0.5)-1]=(TH1F*)data_p1->Clone();
0318 
0319         hV_All_clct_dxdz[int(verticalScale+0.5)-1]=(TH1F*)data_p3->Clone();
0320         hV_All_clct_dxdz[int(verticalScale+0.5)-1]->Add(data_p4);
0321         hV_Eff_clct_dxdz[int(verticalScale+0.5)-1]=(TH1F*)data_p3->Clone();
0322 
0323       }
0324       else{
0325         hV_All_alct_dydz[int(verticalScale+0.5)-1]->Add(data_p1);
0326         hV_All_alct_dydz[int(verticalScale+0.5)-1]->Add(data_p2);
0327         hV_Eff_alct_dydz[int(verticalScale+0.5)-1]->Add(data_p1);
0328 
0329         hV_All_clct_dxdz[int(verticalScale+0.5)-1]->Add(data_p3);
0330         hV_All_clct_dxdz[int(verticalScale+0.5)-1]->Add(data_p4);
0331         hV_Eff_clct_dxdz[int(verticalScale+0.5)-1]->Add(data_p3);
0332 
0333       }
0334 
0335       ++iterations;
0336     }
0337       }
0338     }
0339   }
0340 
0341   TGraphAsymmErrors* h_attachment_efficiency =  new TGraphAsymmErrors(h_attachment_AllHits->GetNbinsX());
0342   h_attachment_efficiency->BayesDivide(h_attachment_EffHits,h_attachment_AllHits);
0343   mytitle = "Attachment efficiency (all chambers); layer number; efficiency ";
0344   myTitle = mytitle.c_str();
0345   h_attachment_efficiency->SetTitle(myTitle);
0346   mytitle = "h_attachment_efficiency";
0347   myTitle = mytitle.c_str();
0348   h_attachment_efficiency->SetName(myTitle);
0349   h_attachment_efficiency->Write();
0350   h_attachment_EffHits->Delete();
0351   h_attachment_AllHits->Delete();
0352 
0353   TGraphAsymmErrors* gV_effALCT_dydz[ny];
0354   TGraphAsymmErrors* gV_effCLCT_dxdz[ny];
0355 
0356   for(int i = 0;i<ny;++i){
0357     if(!i){
0358       h_All_alct_dydz=(TH1F*)hV_All_alct_dydz[i]->Clone();
0359       h_Eff_alct_dydz=(TH1F*)hV_Eff_alct_dydz[i]->Clone();
0360     }
0361     else{
0362       h_All_alct_dydz->Add(hV_All_alct_dydz[i]);
0363       h_Eff_alct_dydz->Add(hV_Eff_alct_dydz[i]);
0364     }
0365     gV_effALCT_dydz[i] = new TGraphAsymmErrors(hV_Eff_alct_dydz[i]->GetNbinsX());
0366     gV_effALCT_dydz[i]->BayesDivide(hV_Eff_alct_dydz[i],hV_All_alct_dydz[i]);
0367     mytitle = "Efficiency - ALCT vs dydz : " + to_string(types[i]);
0368     mytitle += "; local dydz (ME 3 and 4 flipped); efficiency ";
0369     myTitle = mytitle.c_str();
0370     gV_effALCT_dydz[i]->SetTitle(myTitle);
0371     mytitle = "g_effALCT_dydz_" + to_string(types[i]);
0372     myTitle = mytitle.c_str();
0373     gV_effALCT_dydz[i]->SetName(myTitle);
0374     hV_Eff_alct_dydz[i]->Delete();
0375     hV_All_alct_dydz[i]->Delete();
0376     gV_effALCT_dydz[i]->Write();
0377 
0378 
0379     gV_effCLCT_dxdz[i] = new TGraphAsymmErrors(hV_Eff_clct_dxdz[i]->GetNbinsX());
0380     gV_effCLCT_dxdz[i]->BayesDivide(hV_Eff_clct_dxdz[i],hV_All_clct_dxdz[i]);
0381     mytitle = "Efficiency - CLCT vs dxdz : " + to_string(types[i]);
0382     mytitle += "; dxdz (local); efficiency ";
0383     myTitle = mytitle.c_str();
0384     gV_effCLCT_dxdz[i]->SetTitle(myTitle);
0385     mytitle = "g_effCLCT_dxdz_" + to_string(types[i]);
0386     myTitle = mytitle.c_str();
0387     gV_effCLCT_dxdz[i]->SetName(myTitle);
0388     hV_Eff_clct_dxdz[i]->Delete();
0389     hV_All_clct_dxdz[i]->Delete();
0390     gV_effCLCT_dxdz[i]->Write();
0391   }
0392 
0393   TGraphAsymmErrors* g_effALCT_dydz;
0394   g_effALCT_dydz = new TGraphAsymmErrors(h_Eff_alct_dydz->GetNbinsX());
0395   g_effALCT_dydz->BayesDivide(h_Eff_alct_dydz,h_All_alct_dydz);
0396   mytitle = "Efficiency - ALCT vs dydz : All chambers";
0397   mytitle += "; local dydz (ME 3 and 4 flipped); efficiency ";
0398   myTitle = mytitle.c_str();
0399   g_effALCT_dydz->SetTitle(myTitle);
0400   mytitle = "g_effALCT_dydz_AllCh";
0401   myTitle = mytitle.c_str();
0402   g_effALCT_dydz->SetName(myTitle);
0403   h_Eff_alct_dydz->Delete();
0404   h_All_alct_dydz->Delete();
0405   g_effALCT_dydz->Write();
0406   
0407   
0408   //sum_All_alct_theta->Delete();
0409   //sum_Eff_alct_theta->Delete();
0410   h_rhEfficiency->Draw("colz");
0411   
0412   //f1->Close();
0413   //TFile* anaFile = TFile::Open("efficiencies.root", "RECREATE"); // my output file
0414   //sum_histo->Write();
0415   //h_Chi2_ndf->Write();
0416   //h_hitsInSegm->Write();
0417   //h_Residual_Segm->Write();
0418   mytitle = "Chi2/ndf (segments) : All chambers; chi2/ndf; entries";
0419   myTitle = mytitle.c_str();
0420   h_Chi2_ndf->SetName(myTitle);
0421   mytitle = "h_Chi2_ndf_AllCh";
0422   myTitle = mytitle.c_str();
0423   h_Chi2_ndf->SetName(myTitle);
0424 
0425   mytitle = "Rechits in a segment : All chambers;number of rechits;entries";
0426   myTitle = mytitle.c_str();
0427   h_hitsInSegm->SetTitle(myTitle);
0428   mytitle = "h_hitsInSegm_AllCh";
0429   myTitle = mytitle.c_str();
0430   h_hitsInSegm->SetName(myTitle);
0431 
0432   mytitle = "Residuals (track to segment) : All chambers; resirual, cm;entries";
0433   myTitle = mytitle.c_str();
0434   h_Residual_Segm->SetTitle(myTitle);
0435   mytitle = "h_ResidInSegm_AllCh";
0436   myTitle = mytitle.c_str();
0437   h_Residual_Segm->SetName(myTitle);
0438 
0439   h_alctEfficiency->Write();
0440   h_clctEfficiency->Write();
0441   h_corrlctEfficiency->Write();
0442 
0443   h_rhEfficiency->Write();
0444   h_simrhEfficiency->Write();
0445   h_segEfficiency->Write();
0446 
0447   h_stripEfficiency->Write();
0448   h_wireEfficiency->Write();
0449 
0450   anaFile->Write();
0451   anaFile->Close();
0452   std::cout<<" Efficiency calculations done. Open the file efficiency.root to look at the efficiency histograms."<<std::endl;
0453   //f1->Close();
0454 }
0455 
0456 void getEfficiency(float efficientEvents, float allEvents, std::vector<float> &efficiencyResult){
0457   //---- Efficiency with binomial error
0458   float efficiency = 0.;
0459   float effError = 0.;
0460   if(fabs(allEvents)>0.000000001){
0461     efficiency = efficientEvents/allEvents;
0462     if(efficientEvents<allEvents){
0463       effError = sqrt( (1.-efficiency)*efficiency/allEvents );
0464     }
0465     else{
0466       double effTemp = (allEvents -1)/allEvents;
0467       if(allEvents<=1){
0468     effError = 1;
0469       }
0470       else{
0471     effError = sqrt( (1.-effTemp)*effTemp/allEvents );
0472       }
0473     }
0474   }
0475   efficiencyResult.clear();
0476   efficiencyResult.push_back(efficiency);
0477   efficiencyResult.push_back(effError);
0478 }
0479 
0480 void getType(int iE, int iS, int iR, float & verticalScale){
0481   if (1==iS){
0482     if(4==iR){
0483       //verticalScale = 0.5;
0484     }
0485     else if(1==iR){
0486       verticalScale = 0.5;
0487     }
0488     else if(2==iR){
0489       verticalScale = 1.5;
0490     }
0491     else if(3==iR){
0492       verticalScale = 2.5;
0493     }
0494   }
0495   else if (2==iS){
0496     if(1==iR){
0497       verticalScale = 3.5;
0498     }
0499     else if(2==iR){
0500       verticalScale = 4.5;
0501     }
0502   }
0503   else if(3==iS){
0504     if(1==iR){
0505       verticalScale = 5.5;
0506     }
0507     else if(2==iR){
0508       verticalScale = 6.5;
0509     }
0510   }
0511   else if( 4==iS){
0512     if(1==iR){
0513       verticalScale = 7.5;
0514     }
0515   }
0516   if(2==iE){
0517     verticalScale = - verticalScale;
0518   } 
0519   verticalScale +=8.5;
0520   //std:cout<<" IN: iE"<<iE<<" iS = "<<iS<<" iR = "<<iR<<" verticalScale = "<<verticalScale<<std::endl;
0521 }