Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:56

0001 // Migrated to use DQMEDHarvester by: Jyothsna Rani Komaragiri, Oct 2014
0002 
0003 #include "DQMOffline/Trigger/interface/JetMETHLTOfflineClient.h"
0004 
0005 #include "FWCore/Framework/interface/Run.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009 
0010 #include "DQMServices/Core/interface/DQMStore.h"
0011 
0012 JetMETHLTOfflineClient::JetMETHLTOfflineClient(const edm::ParameterSet& iConfig) : conf_(iConfig) {
0013   debug_ = false;
0014   verbose_ = false;
0015 
0016   processname_ = iConfig.getParameter<std::string>("processname");
0017 
0018   hltTag_ = iConfig.getParameter<std::string>("hltTag");
0019   if (debug_)
0020     std::cout << hltTag_ << std::endl;
0021 
0022   dirName_ = iConfig.getParameter<std::string>("DQMDirName");
0023 }
0024 
0025 JetMETHLTOfflineClient::~JetMETHLTOfflineClient() = default;
0026 
0027 void JetMETHLTOfflineClient::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0028   ibooker.setCurrentFolder(dirName_);
0029 
0030   LogDebug("JetMETHLTOfflineClient") << "dqmEndJob" << std::endl;
0031   if (debug_)
0032     std::cout << "dqmEndJob" << std::endl;
0033 
0034   std::vector<MonitorElement*> hltMEs;
0035 
0036   // Look at all folders, go to the subfolder which includes the string "Eff"
0037   std::vector<std::string> fullPathHLTFolders = igetter.getSubdirs();
0038   for (auto& fullPathHLTFolder : fullPathHLTFolders) {
0039     // Move on only if the folder name contains "Eff" Or "Trigger Summary"
0040     if (debug_)
0041       std::cout << fullPathHLTFolder << std::endl;
0042     if ((fullPathHLTFolder.find("Eff") != std::string::npos)) {
0043       ibooker.setCurrentFolder(fullPathHLTFolder);
0044     } else {
0045       continue;
0046     }
0047 
0048     // Look at all subfolders, go to the subfolder which includes the string "Eff"
0049     std::vector<std::string> fullSubPathHLTFolders = igetter.getSubdirs();
0050     for (auto& fullSubPathHLTFolder : fullSubPathHLTFolders) {
0051       if (debug_)
0052         std::cout << fullSubPathHLTFolder << std::endl;
0053       ibooker.setCurrentFolder(fullSubPathHLTFolder);
0054 
0055       // Look at all MonitorElements in this folder
0056       hltMEs = igetter.getContents(fullSubPathHLTFolder);
0057       LogDebug("JetMETHLTOfflineClient") << "Number of MEs for this HLT path = " << hltMEs.size() << std::endl;
0058 
0059       for (unsigned int k = 0; k < hltMEs.size(); k++) {
0060         if (debug_)
0061           std::cout << hltMEs[k]->getName() << std::endl;
0062 
0063         //-----
0064         if ((hltMEs[k]->getName().find("ME_Numerator") != std::string::npos) &&
0065             hltMEs[k]->getName().find("ME_Numerator") == 0) {
0066           std::string name = hltMEs[k]->getName();
0067           name.erase(0, 12);  // Removed "ME_Numerator"
0068           if (debug_)
0069             std::cout << "==name==" << name << std::endl;
0070           //      if( name.find("EtaPhi") !=std::string::npos ) continue; // do not consider EtaPhi 2D plots
0071 
0072           for (unsigned int l = 0; l < hltMEs.size(); l++) {
0073             if (hltMEs[l]->getName() == "ME_Denominator" + name) {
0074               // found denominator too
0075               if (name.find("EtaPhi") != std::string::npos) {
0076                 TH2F* tNumerator = hltMEs[k]->getTH2F();
0077                 TH2F* tDenominator = hltMEs[l]->getTH2F();
0078 
0079                 std::string title = "Eff_" + hltMEs[k]->getTitle();
0080 
0081                 auto* teff = (TH2F*)tNumerator->Clone(title.c_str());
0082                 teff->Divide(tNumerator, tDenominator, 1, 1);
0083                 ibooker.book2D("ME_Eff_" + name, teff);
0084                 delete teff;
0085               } else {
0086                 TH1F* tNumerator = hltMEs[k]->getTH1F();
0087                 TH1F* tDenominator = hltMEs[l]->getTH1F();
0088 
0089                 std::string title = "Eff_" + hltMEs[k]->getTitle();
0090 
0091                 auto* teff = (TH1F*)tNumerator->Clone(title.c_str());
0092                 teff->Divide(tNumerator, tDenominator, 1, 1);
0093                 ibooker.book1D("ME_Eff_" + name, teff);
0094                 delete teff;
0095               }
0096             }  // Denominator
0097           }    // Loop-l
0098         }      // Numerator
0099 
0100       }  // Loop-k
0101     }    // fullSubPath
0102   }      // fullPath
0103 }
0104 
0105 #include "FWCore/Framework/interface/MakerMacros.h"
0106 DEFINE_FWK_MODULE(JetMETHLTOfflineClient);