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 // Uncomment the following line for some extra debug information
0007 // #define DEBUG
0008 
0009 void RecoValHistoPublisher(const char* newFile = "NEW_FILE", const char* refFile = "REF_FILE") {
0010   cout << ">> Starting RecoValHistoPublisher(" << newFile << "," << refFile << ")..." << endl;
0011 
0012   //====  To be replaced from python3 ====================
0013 
0014   const char* dataType = "DATATYPE";
0015   const char* refLabel("REF_LABEL, REF_RELEASE REFSELECTION");
0016   const char* newLabel("NEW_LABEL, NEW_RELEASE NEWSELECTION");
0017 
0018   // ==== Initial settings and loads
0019   //gROOT->ProcessLine(".x HistoCompare_Tracks.C");
0020   //gROOT ->Reset();
0021   gROOT->SetBatch();
0022   gErrorIgnoreLevel = kWarning;  // Get rid of the info messages
0023 
0024   SetGlobalStyle();
0025 
0026   // ==== Some cleaning... is this needed?
0027   delete gROOT->GetListOfFiles()->FindObject(refFile);
0028   delete gROOT->GetListOfFiles()->FindObject(newFile);
0029 
0030   // ==== Opening files, moving to the right branch and getting the list of sub-branches
0031   cout << ">> Openning file, moving to the right branch and getting sub-branches..." << endl;
0032 
0033   cout << ">> Finding sources..." << endl;
0034   TFile* sfile = new TFile(newFile);
0035   TList* sl = getListOfBranches(dataType, sfile, "MuonRecoAnalyzer");
0036   if (!sl) {
0037     cout << "ERROR: Could not find keys!!!" << endl;
0038     cerr << "ERROR: Could not find keys!!!" << endl;
0039     return;
0040   }
0041   TDirectory* sdir = gDirectory;
0042   for (unsigned int i = 0; i < sl->GetEntries(); i++)
0043     cout << "   + " << sl->At(i)->GetName() << endl;
0044 
0045   cout << ">> Finding references..." << endl;
0046   TFile* rfile = new TFile(refFile);
0047   TList* rl = getListOfBranches(dataType, rfile, "MuonRecoAnalyzer");
0048   if (!rl) {
0049     cout << "ERROR: Could not find keys!!!" << endl;
0050     cerr << "ERROR: Could not find keys!!!" << endl;
0051     return;
0052   }
0053   TDirectory* rdir = gDirectory;
0054   for (unsigned int i = 0; i < sl->GetEntries(); i++)
0055     cout << "   + " << sl->At(i)->GetName() << endl;
0056 
0057   // Get the number of events for the normalization:
0058   TH1F *sevt, *revt;
0059   sdir->GetObject("RecoMuonV/RecoMuon_MuonAssoc_Glb/NMuon", sevt);
0060   rdir->GetObject("RecoMuonV/RecoMuon_MuonAssoc_Glb/NMuon", revt);
0061   /*
0062   if (sevt && revt) {
0063     if (revt->GetEntries()>0) 
0064       snorm = sevt->GetEntries()/revt->GetEntries();
0065     cout << "   + SEntries = " << sevt->GetEntries()
0066      << ", REntries = " << revt->GetEntries() 
0067      << ", Normalization = " << snorm << endl;
0068     cout << "   + SIntegral = " << sevt->Integral()
0069      << ", RIntegral = " << revt->Integral() 
0070      << ", Normalization = " << snorm << endl;
0071   }
0072   else {  
0073     cout << "WARNING: Missing normalization histos!" << endl; 
0074   }
0075   */
0076   Float_t maxPT;
0077   TString File = newFile;
0078 
0079   if (File.Contains("SingleMuPt10")) {
0080     maxPT = 70.;
0081   } else if (File.Contains("SingleMuPt100")) {
0082     maxPT = 400.;
0083   } else if (File.Contains("SingleMuPt1000") || File.Contains("WpM") || File.Contains("ZpMM"))
0084     maxPT = 1400.;
0085   else
0086     maxPT = 300.;
0087 
0088   //==== Iterate now over histograms and collections
0089   cout << ">> Iterating over histograms and collections..." << endl;
0090 
0091   TIter iter_r(rl);
0092   TIter iter_s(sl);
0093   TKey* rKey = 0;
0094   TKey* sKey = 0;
0095   TString rcollname;
0096   TString scollname;
0097   while ((rKey = (TKey*)iter_r())) {
0098     TString myName = rKey->GetName();
0099 #ifdef DEBUG
0100     cout << "DEBUG: Checking key " << myName << endl;
0101 #endif
0102     rcollname = myName;
0103     sKey = (TKey*)iter_s();
0104     if (!sKey)
0105       continue;
0106     scollname = sKey->GetName();
0107     if ((rcollname != scollname) && (rcollname + "FS" != scollname) && (rcollname != scollname + "FS")) {
0108       cerr << "ERROR: Different collection names, please check: " << rcollname << " : " << scollname << endl;
0109       cout << "ERROR: Different collection names, please check: " << rcollname << " : " << scollname << endl;
0110       continue;
0111     }
0112 
0113     // ==== Now let's go for the plotting...
0114     cout << ">> Comparing plots in " << myName << "..." << endl;
0115     cerr << ">> Comparing plots in " << myName << "..." << endl;
0116     TString newDir("NEW_RELEASE/NEWSELECTION/NEW_LABEL/");
0117     newDir += myName;
0118     gSystem->mkdir(newDir, kTRUE);
0119 
0120     bool resolx = false;
0121     bool* resol = &resolx;
0122     bool logy[] = {false, false, false, false};
0123     bool doKolmo[] = {true, true, true, true};
0124     Double_t minx[] = {-1E100, -1E100, 5., -1E100, -1E100, -1E100};
0125     Double_t maxx[] = {-1E100, -1E100, maxPT, -1E100, -1E100, -1E100};
0126     Double_t snorm[] = {0., 0., 0., 0., 0., 0.};
0127     //===== reco muon distributions: GLB_GLB
0128     const char* plots1[] = {"GlbMuon_Glb_eta", "GlbMuon_Glb_phi", "GlbMuon_Glb_pt", "GlbMuon_Glb_chi2OverDf"};
0129     const char* plotst1[] = {
0130         "GlobalMuon(GLB) #eta", "GlobalMuon(GLB) #phi", "GlobalMuon(GLB) pT", "GlobalMuon(GLB) #chi^{2}/ndf"};
0131     Plot4Histograms(newDir + "/muonReco1",
0132                     rdir,
0133                     sdir,
0134                     rcollname,
0135                     scollname,
0136                     "RecoHistos1",
0137                     "Distributions for GlobalMuons (GLB)",
0138                     refLabel,
0139                     newLabel,
0140                     plots1,
0141                     plotst1,
0142                     logy,
0143                     doKolmo,
0144                     snorm,
0145                     resol,
0146                     minx,
0147                     maxx);
0148 
0149     //===== reco muon distributions: GLB_STA
0150     const char* plots2[] = {"GlbMuon_Sta_eta", "GlbMuon_Sta_phi", "GlbMuon_Sta_pt", "GlbMuon_Sta_chi2OverDf"};
0151     const char* plotst2[] = {
0152         "GlobalMuon(STA) #eta", "GlobalMuon(STA) #phi", "GlobalMuon(STA) p_T", "GlobalMuon(STA) #chi^{2}/ndf"};
0153     Double_t minx1[] = {-1E100, -1E100, 5., -1E100, -1E100};
0154     Double_t maxx1[] = {-1E100, -1E100, maxPT, -1E100, -1E100};
0155 
0156     Plot4Histograms(newDir + "/muonReco2",
0157                     rdir,
0158                     sdir,
0159                     rcollname,
0160                     scollname,
0161                     "RecoHistos2",
0162                     "Distributions for GlobalMuons (STA)",
0163                     refLabel,
0164                     newLabel,
0165                     plots2,
0166                     plotst2,
0167                     logy,
0168                     doKolmo,
0169                     snorm,
0170                     resol,
0171                     minx1,
0172                     maxx1);
0173 
0174     //===== reco muon distributions: GLB_TK
0175     const char* plots3[] = {"GlbMuon_Tk_eta", "GlbMuon_Tk_phi", "GlbMuon_Tk_pt", "GlbMuon_Tk_chi2OverDf"};
0176     const char* plotst3[] = {
0177         "GlobalMuon(TK) #eta", "GlobalMuon(TK) #phi", "GlobalMuon(TK) pT", "GlobalMuon(TK) #chi^{2}/ndf"};
0178     Plot4Histograms(newDir + "/muonReco3",
0179                     rdir,
0180                     sdir,
0181                     rcollname,
0182                     scollname,
0183                     "RecoHistos3",
0184                     "Distributions for GlobalMuons (TK)",
0185                     refLabel,
0186                     newLabel,
0187                     plots3,
0188                     plotst3,
0189                     logy,
0190                     doKolmo,
0191                     snorm,
0192                     resol,
0193                     minx1,
0194                     maxx1);
0195 
0196     //===== reco muon distributions: STA
0197     const char* plots4[] = {"StaMuon_eta", "StaMuon_phi", "StaMuon_pt", "StaMuon_chi2OverDf"};
0198     const char* plotst4[] = {"StaMuon #eta", "StaMuon #phi", "StaMuon p_T", "StaMuon #chi^{2}/ndf"};
0199     Plot4Histograms(newDir + "/muonReco4",
0200                     rdir,
0201                     sdir,
0202                     rcollname,
0203                     scollname,
0204                     "RecoHistos4",
0205                     "Distributions for StandAlone Muons",
0206                     refLabel,
0207                     newLabel,
0208                     plots4,
0209                     plotst4,
0210                     logy,
0211                     doKolmo,
0212                     snorm,
0213                     resol,
0214                     minx1,
0215                     maxx1);
0216 
0217     //===== reco muon distributions: Tracker Muons
0218     const char* plots5[] = {"TkMuon_eta", "TkMuon_phi", "TkMuon_pt", "TkMuon_chi2OverDf"};
0219     const char* plotst5[] = {"TkMuon #eta", "TkMuon #phi", "TkMuon p_T", "TkMuon #chi^{2}/ndf"};
0220     Plot4Histograms(newDir + "/muonReco5",
0221                     rdir,
0222                     sdir,
0223                     rcollname,
0224                     scollname,
0225                     "RecoHistos5",
0226                     "Distributions for Tracker Muons",
0227                     refLabel,
0228                     newLabel,
0229                     plots5,
0230                     plotst5,
0231                     logy,
0232                     doKolmo,
0233                     snorm,
0234                     resol,
0235                     minx1,
0236                     maxx1);
0237 
0238     //// Merge pdf histograms together into larger files, and name them based on the collection names
0239     TString mergefile = "merged_reco.pdf";                 // File name where partial pdfs will be merged
0240     TString destfile = newDir + "/../" + myName + ".pdf";  // Destination file name
0241     TString gscommand = "gs -dBATCH -dNOPAUSE -q -sDEVICE=pdfwrite -sOutputFile=" + mergefile + " " + newDir +
0242                         "/muonReco1.pdf " + newDir + "/muonReco2.pdf " + newDir + "/muonReco3.pdf " + newDir +
0243                         "/muonReco4.pdf " + newDir + "/muonReco5.pdf ";
0244 
0245     cout << ">> Merging partial pdfs to " << mergefile << "..." << endl;
0246 #ifdef DEBUG
0247     cout << "DEBUG: ...with command \"" << gscommand << "\"" << endl;
0248 #endif
0249     gSystem->Exec(gscommand);
0250     cout << ">> Moving " << mergefile << " to " << destfile << "..." << endl;
0251     gSystem->Rename(mergefile, destfile);
0252 
0253     cout << ">> Deleting partial pdf files" << endl;
0254     gSystem->Exec("rm -r " + newDir + "/*.pdf");
0255     cout << "   ... Done" << endl;
0256 
0257   }  // end of "while loop"
0258 
0259   cout << ">> Removing the relval files from ROOT before closing..." << endl;
0260   gROOT->GetListOfFiles()->Remove(sfile);
0261   gROOT->GetListOfFiles()->Remove(rfile);
0262 
0263 #ifdef DEBUG
0264   cout << "DEBUG: Exiting!" << endl;
0265   cerr << "DEBUG: Exiting!" << endl;
0266 #endif
0267 }