Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <vector>
0002 #include <algorithm>
0003 #include "TMath.h"
0004 #include "macro/PlotHelpers.C"
0005 
0006 /////
0007 // Uncomment the following line to get more debuggin output
0008 // #define DEBUG
0009 
0010 void RecoMuonValHistoPublisher(const char* newFile = "NEW_FILE", const char* refFile = "REF_FILE") {
0011   cout << ">> Starting RecoMuonValHistoPublisher(" << newFile << "," << refFile << ")..." << endl;
0012 
0013   //====  To be replaced from python3 ====================
0014 
0015   const char* dataType = "DATATYPE";
0016   const char* refLabel("REF_LABEL, REF_RELEASE REFSELECTION");
0017   const char* newLabel("NEW_LABEL, NEW_RELEASE NEWSELECTION");
0018   const char* fastSim = "IS_FSIM";
0019 
0020   // ==== Initial settings and loads
0021   //gROOT->ProcessLine(".x HistoCompare_Tracks.C");
0022   //gROOT ->Reset();
0023   gROOT->SetBatch();
0024   gErrorIgnoreLevel = kWarning;  // Get rid of the info messages
0025 
0026   SetGlobalStyle();
0027 
0028   // ==== Some cleaning... is this needed?
0029   delete gROOT->GetListOfFiles()->FindObject(refFile);
0030   delete gROOT->GetListOfFiles()->FindObject(newFile);
0031 
0032   // ==== Opening files, moving to the right branch and getting the list of sub-branches
0033   cout << ">> Openning file, moving to the right branch and getting sub-branches..." << endl;
0034 
0035   cout << ">> Finding sources..." << endl;
0036   TFile* sfile = new TFile(newFile);
0037   TList* sl = getListOfBranches(dataType, sfile, "RecoMuonV");
0038   if (!sl) {
0039     cout << "ERROR: Could not find keys!!!" << endl;
0040     cerr << "ERROR: Could not find keys!!!" << endl;
0041     return;
0042   }
0043   TDirectory* sdir = gDirectory;
0044   for (unsigned int i = 0; i < sl->GetEntries(); i++)
0045     cout << "   + " << sl->At(i)->GetName() << endl;
0046 
0047   cout << ">> Finding references..." << endl;
0048   TFile* rfile = new TFile(refFile);
0049   TList* rl = getListOfBranches(dataType, rfile, "RecoMuonV");
0050   if (!rl) {
0051     cout << "ERROR: Could not find keys!!!" << endl;
0052     cerr << "ERROR: Could not find keys!!!" << endl;
0053     return;
0054   }
0055   TDirectory* rdir = gDirectory;
0056   for (unsigned int i = 0; i < sl->GetEntries(); i++)
0057     cout << "   + " << sl->At(i)->GetName() << endl;
0058 
0059   Float_t maxPT;
0060   TString File = newFile;
0061   if (File.Contains("SingleMuPt1000") || File.Contains("WpM") || File.Contains("ZpMM"))
0062     maxPT = 1400.;
0063   else if (File.Contains("SingleMuPt10")) {
0064     maxPT = 70.;
0065   } else if (File.Contains("SingleMuPt100")) {
0066     maxPT = 400.;
0067   } else
0068     maxPT = 400.;
0069 
0070   TIter iter_r(rl);
0071   TIter iter_s(sl);
0072   TKey* rKey = 0;
0073   TKey* sKey = 0;
0074   TString rcollname;
0075   TString scollname;
0076 
0077   while ((rKey = (TKey*)iter_r())) {
0078     TString myName = rKey->GetName();
0079 #ifdef DEBUG
0080     cout << "DEBUG: Checking key " << myName << endl;
0081 #endif
0082     rcollname = myName;
0083     sKey = (TKey*)iter_s();
0084     if (!sKey)
0085       continue;
0086     scollname = sKey->GetName();
0087     if ((rcollname != scollname) && (rcollname + "FS" != scollname) && (rcollname != scollname + "FS")) {
0088       cerr << "ERROR: Different collection names, please check: " << rcollname << " : " << scollname << endl;
0089       cout << "ERROR: Different collection names, please check: " << rcollname << " : " << scollname << endl;
0090       continue;
0091     }
0092 
0093     // ==== Now let's go for the plotting...
0094     cout << ">> Comparing plots in " << myName << "..." << endl;
0095     cerr << ">> Comparing plots in " << myName << "..." << endl;
0096     TString newDir("NEW_RELEASE/NEWSELECTION/NEW_LABEL/");
0097     newDir += myName;
0098     gSystem->mkdir(newDir, kTRUE);
0099     bool resolx = false;
0100     bool* resol = &resolx;
0101     bool logy[] = {false, false, false, false};
0102     bool doKolmo[] = {true, true, true, true};
0103     Double_t minx[] = {-1E100, -1E100, -1E100, 5., -1E100, -1E100};
0104     Double_t maxx[] = {-1E100, -1E100, -1E100, maxPT, -1E100, -1E100};
0105 
0106     Double_t norm[] = {0., 0., -999., -999., 0., 0.};  //Normalize to first histogram
0107 
0108     //===== reco muon distributions: GLB
0109     //TString baseh     = Form("RecoMuon_MuonAssoc_Glb%s/",fastSim);
0110 
0111     const char* plots1[] = {"RecoMuon_MuonAssoc_Glb/ErrPt",
0112                             "RecoMuon_MuonAssoc_Glb/ErrP",
0113                             "RecoMuon_MuonAssoc_Glb/ErrPt_vs_Eta_Sigma",
0114                             "RecoMuon_MuonAssoc_Glb/ErrPt_vs_Pt_Sigma"};
0115     const char* plotst1[] = {"GlobalMuon(GLB) #Delta p_{T}/p_{T}",
0116                              "GlobalMuon(GLB) #Delta p/p",
0117                              "GlobalMuon(GLB) #Delta p_{T}/p_{T} vs #sigma(#eta)",
0118                              "GlobalMuon(GLB) #Delta p_{T}/p_{T} vs #sigma(p_{T})"};
0119     Plot4Histograms(newDir + "/muonRecoGlb",
0120                     rdir,
0121                     sdir,
0122                     rcollname,
0123                     scollname,
0124                     "RecHistosGlb",
0125                     "Distributions for GlobalMuons (GLB)",
0126                     refLabel,
0127                     newLabel,
0128                     plots1,
0129                     plotst1,
0130                     logy,
0131                     doKolmo,
0132                     norm,
0133                     resol,
0134                     minx,
0135                     maxx);
0136 
0137     //==== efficiencies and fractions GLB
0138     const char* plots2[] = {"RecoMuon_MuonAssoc_Glb/EffP",
0139                             "RecoMuon_MuonAssoc_Glb/EffEta",
0140                             "RecoMuon_MuonAssoc_Glb/FractP",
0141                             "RecoMuon_MuonAssoc_Glb/FractEta"};
0142     const char* plotst2[] = {"GlobalMuon(GLB) #epsilon vs. p",
0143                              "GlobalMuon(GLB) #epsilon vs. #eta",
0144                              "GlobalMuon(GLB) fraction vs. p",
0145                              "GlobalMuon(GLB) fraction vs. #eta"};
0146     Double_t minx1[] = {5., -1E100, 5., -1E100, -1E100, -1E100};
0147     Double_t maxx1[] = {maxPT, -1E100, maxPT, -1E100, -1E100, -1E100};
0148     Double_t norm2[] = {-999., -999., -999., -999., -999., -999.};  //Normalize to first histogram
0149     Plot4Histograms(newDir + "/muonRecoGlbEff",
0150                     rdir,
0151                     sdir,
0152                     rcollname,
0153                     scollname,
0154                     "RecEffHistosGlb",
0155                     "Distributions for GlobalMuons (GLB), efficiencies and fractions",
0156                     refLabel,
0157                     newLabel,
0158                     plots2,
0159                     plotst2,
0160                     logy,
0161                     doKolmo,
0162                     norm2,
0163                     resol,
0164                     minx1,
0165                     maxx1);
0166 
0167     /*
0168     //===== reco muon distributions: GLBPF
0169     baseh             = Form("RecoMuon_MuonAssoc_GlbPF%s/",fastSim);
0170     const char* plots3[]  = {(baseh + "ErrPt").Data(), (baseh + "ErrP").Data(), 
0171                  (baseh + "ErrPt_vs_Eta_Sigma").Data(), (baseh + "ErrPt_vs_Pt_Sigma").Data()};   
0172     const char* plotst3[] = {"PFGlobalMuon(GLBPF) #Delta p_{T}/p_{T}", "PFGlobalMuon(GLBPF) #Delta p/p", 
0173                  "PFGlobalMuon(GLBPF) #Delta p_{T}/p_{T} vs #sigma(#eta)", "PFGlobalMuon(GLBPF) #Delta p_{T}/p_{T} vs #sigma(p_{T})"};
0174     Plot4Histograms(newDir + "/muonRecoGlbPF",
0175             rdir, sdir, 
0176             rcollname, scollname,
0177             "RecHistosGlbPF", "Distributions for PFGlobalMuons (GLBPF)",
0178             refLabel, newLabel,
0179             plots3, plotst3,
0180             logy, doKolmo, norm);
0181     
0182     
0183     //==== efficiencies and fractions GLBPF
0184     const char* plots4 [] = {(baseh + "EffP").Data(), (baseh + "EffEta").Data(), 
0185                  (baseh + "FractP").Data(), (baseh + "FractEta").Data()};   
0186     const char* plotst4[] = {"PFGlobalMuon(GLBPF) #epsilon vs. p", "PFGlobalMuon(GLBPF) #epsilon vs. #eta", 
0187                  "PFGlobalMuon(GLBPF) fraction vs. p", "PFGlobalMuon(GLBPF) fraction vs. #eta"};
0188     Plot4Histograms(newDir + "/muonRecoGlbPFEff",
0189             rdir, sdir, 
0190             rcollname, scollname,
0191             "RecEffHistosGlbPF", "Distributions for PFGlobalMuons (GLBPF), efficiencies and fractions",
0192             refLabel, newLabel,
0193             plots4, plotst4,
0194             logy, doKolmo, norm);
0195     */
0196 
0197     //===== reco muon distributions: STA
0198     //baseh             = Form("RecoMuon_MuonAssoc_Sta%s/",fastSim);
0199     const char* plots5[] = {"RecoMuon_MuonAssoc_Sta/ErrPt",
0200                             "RecoMuon_MuonAssoc_Sta/ErrP",
0201                             "RecoMuon_MuonAssoc_Sta/ErrPt_vs_Eta_Sigma",
0202                             "RecoMuon_MuonAssoc_Sta/ErrPt_vs_Pt_Sigma"};
0203     const char* plotst5[] = {"StandAloneMuon(STA) #Delta p_{T}/p_{T}",
0204                              "StandAloneMuon(STA) #Delta p/p",
0205                              "StandAloneMuon(STA) #Delta p_{T}/p_{T} vs #sigma(#eta)",
0206                              "StandAloneMuon(STA) #Delta p_{T}/p_{T} vs #sigma(p_{T})"};
0207     Plot4Histograms(newDir + "/muonRecoSta",
0208                     rdir,
0209                     sdir,
0210                     rcollname,
0211                     scollname,
0212                     "RecHistosSta",
0213                     "Distributions for StandAloneMuons (STA)",
0214                     refLabel,
0215                     newLabel,
0216                     plots5,
0217                     plotst5,
0218                     logy,
0219                     doKolmo,
0220                     norm,
0221                     resol,
0222                     minx,
0223                     maxx);
0224 
0225     //==== efficiencies and fractions STA
0226     const char* plots6[] = {"RecoMuon_MuonAssoc_Sta/EffP",
0227                             "RecoMuon_MuonAssoc_Sta/EffEta",
0228                             "RecoMuon_MuonAssoc_Sta/FractP",
0229                             "RecoMuon_MuonAssoc_Sta/FractEta"};
0230     const char* plotst6[] = {"StandAloneMuon(STA) #epsilon vs. p",
0231                              "StandAloneMuon(STA) #epsilon vs. #eta",
0232                              "StandAloneMuon(STA) fraction vs. p",
0233                              "StandAloneMuon(STA) fraction vs. #eta"};
0234     Plot4Histograms(newDir + "/muonRecoStaEff",
0235                     rdir,
0236                     sdir,
0237                     rcollname,
0238                     scollname,
0239                     "RecEffHistosSta",
0240                     "Distributions for StandAloneMuons (STA), efficiencies and fractions",
0241                     refLabel,
0242                     newLabel,
0243                     plots6,
0244                     plotst6,
0245                     logy,
0246                     doKolmo,
0247                     norm2,
0248                     resol,
0249                     minx1,
0250                     maxx1);
0251 
0252     //===== reco muon distributions: TRK
0253     //baseh             = Form("RecoMuon_MuonAssoc_Trk%s/",fastSim);
0254     const char* plots7[] = {"RecoMuon_MuonAssoc_Trk/ErrPt",
0255                             "RecoMuon_MuonAssoc_Trk/ErrP",
0256                             "RecoMuon_MuonAssoc_Trk/ErrPt_vs_Eta_Sigma",
0257                             "RecoMuon_MuonAssoc_Trk/ErrPt_vs_Pt_Sigma"};
0258     const char* plotst7[] = {"TrackerMuon(TRK) #Delta p_{T}/p_{T}",
0259                              "TrackerMuon(TRK) #Delta p/p",
0260                              "TrackerMuon(TRK) #Delta p_{T}/p_{T} vs #sigma(#eta)",
0261                              "TrackerMuon(TRK) #Delta p_{T}/p_{T} vs #sigma(p_{T})"};
0262     Plot4Histograms(newDir + "/muonRecoTrk",
0263                     rdir,
0264                     sdir,
0265                     rcollname,
0266                     scollname,
0267                     "RecHistosTrk",
0268                     "Distributions for TrackerMuons (TRK)",
0269                     refLabel,
0270                     newLabel,
0271                     plots7,
0272                     plotst7,
0273                     logy,
0274                     doKolmo,
0275                     norm,
0276                     resol,
0277                     minx,
0278                     maxx);
0279 
0280     //==== efficiencies and fractions TRK
0281     const char* plots8[] = {"RecoMuon_MuonAssoc_Trk/EffP",
0282                             "RecoMuon_MuonAssoc_Trk/EffEta",
0283                             "RecoMuon_MuonAssoc_Trk/FractP",
0284                             "RecoMuon_MuonAssoc_Trk/FractEta"};
0285     const char* plotst8[] = {"TrackerMuon(TRK) #epsilon vs. p",
0286                              "TrackerMuon(TRK) #epsilon vs. #eta",
0287                              "TrackerMuon(TRK) fraction vs. p",
0288                              "TrackerMuon(TRK) fraction vs. #eta"};
0289     Plot4Histograms(newDir + "/muonRecoTrkEff",
0290                     rdir,
0291                     sdir,
0292                     rcollname,
0293                     scollname,
0294                     "RecEffHistosTrk",
0295                     "Distributions for TrackerMuons (TRK), efficiencies and fractions",
0296                     refLabel,
0297                     newLabel,
0298                     plots8,
0299                     plotst8,
0300                     logy,
0301                     doKolmo,
0302                     norm2,
0303                     resol,
0304                     minx1,
0305                     maxx1);
0306 
0307     //
0308     //===== reco muon distributions: Tight Muons
0309     //
0310     //baseh             = Form("RecoMuon_MuonAssoc_Tgt%s/",fastSim);
0311     const char* plots9[] = {"RecoMuon_MuonAssoc_Tgt/ErrPt",
0312                             "RecoMuon_MuonAssoc_Tgt/ErrP",
0313                             "RecoMuon_MuonAssoc_Tgt/ErrPt_vs_Eta_Sigma",
0314                             "RecoMuon_MuonAssoc_Tgt/ErrPt_vs_Pt_Sigma"};
0315     const char* plotst9[] = {"Tight Muon #Delta p_{T}/p_{T}",
0316                              "Tight Muon #Delta p/p",
0317                              "Tight Muon #Delta p_{T}/p_{T} vs #sigma(#eta)",
0318                              "Tight Muon #Delta p_{T}/p_{T} vs #sigma(p_{T})"};
0319     Plot4Histograms(newDir + "/muonRecoTgt",
0320                     rdir,
0321                     sdir,
0322                     rcollname,
0323                     scollname,
0324                     "RecHistosTgt",
0325                     "Distributions for Tight Muons",
0326                     refLabel,
0327                     newLabel,
0328                     plots9,
0329                     plotst9,
0330                     logy,
0331                     doKolmo,
0332                     norm,
0333                     resol,
0334                     minx,
0335                     maxx);
0336 
0337     //==== efficiencies and fractions Tight Muons
0338     const char* plots10[] = {"RecoMuon_MuonAssoc_Tgt/EffP",
0339                              "RecoMuon_MuonAssoc_Tgt/EffEta",
0340                              "RecoMuon_MuonAssoc_Tgt/FractP",
0341                              "RecoMuon_MuonAssoc_Tgt/FractEta"};
0342     const char* plotst10[] = {"Tight Muon #epsilon vs. p",
0343                               "Tight Muon #epsilon vs. #eta",
0344                               "Tight Muon fraction vs. p",
0345                               "Tight Muon fraction vs. #eta"};
0346     Plot4Histograms(newDir + "/muonRecoTgtEff",
0347                     rdir,
0348                     sdir,
0349                     rcollname,
0350                     scollname,
0351                     "RecEffHistosTgt",
0352                     "Distributions for Tight Muons, efficiencies and fractions",
0353                     refLabel,
0354                     newLabel,
0355                     plots10,
0356                     plotst10,
0357                     logy,
0358                     doKolmo,
0359                     norm2,
0360                     resol,
0361                     minx1,
0362                     maxx1);
0363 
0364     //
0365     // Merge pdf histograms together into larger files, and name them based on the collection names
0366     //
0367     TString mergefile = "merged_recomuonval.pdf";          // File name where partial pdfs will be merged
0368     TString destfile = newDir + "/../" + myName + ".pdf";  // Destination file name
0369     TString gscommand = "gs -dBATCH -dNOPAUSE -q -sDEVICE=pdfwrite -sOutputFile=" + mergefile + " " + newDir +
0370                         "/muonRecoGlb.pdf " + newDir +
0371                         "/muonRecoGlbEff.pdf "
0372                         //      +newDir+"/muonRecoGlbPF.pdf "
0373                         //      +newDir+"/muonRecoGlbPFEff.pdf "
0374                         + newDir + "/muonRecoSta.pdf " + newDir + "/muonRecoStaEff.pdf " + newDir +
0375                         "/muonRecoTrk.pdf " + newDir + "/muonRecoTrkEff.pdf " + newDir + "/muonRecoTgt.pdf " + newDir +
0376                         "/muonRecoTgtEff.pdf ";
0377 
0378     cout << ">> Merging partial pdfs to " << mergefile << "..." << endl;
0379 #ifdef DEBUG
0380     cout << "DEBUG: ...with command \"" << gscommand << "\"" << endl;
0381 #endif
0382     gSystem->Exec(gscommand);
0383     cout << ">> Moving " << mergefile << " to " << destfile << "..." << endl;
0384     gSystem->Rename(mergefile, destfile);
0385     cout << "   ... Done" << endl;
0386 
0387     cout << ">> Deleting partial pdf files" << endl;
0388     gSystem->Exec("rm -r " + newDir + "/*.pdf");
0389 
0390   }  // end of "while loop"
0391 
0392   cout << ">> Removing the relval files from ROOT before closing..." << endl;
0393   gROOT->GetListOfFiles()->Remove(sfile);
0394   gROOT->GetListOfFiles()->Remove(rfile);
0395 
0396 #ifdef DEBUG
0397   cout << "DEBUG: Exiting!" << endl;
0398   cerr << "DEBUG: Exiting!" << endl;
0399 #endif
0400 }