Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class DQMHistNormalizer
0002  *
0003  *  Class to produce efficiency histograms by dividing nominator by denominator histograms
0004  *
0005  *  \author Christian Veelken, UC Davis
0006  */
0007 
0008 // framework & common header files
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 
0016 //DQM services
0017 #include "DQMServices/Core/interface/DQMStore.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 
0020 //Regexp handling
0021 #include <regex>
0022 #include <string>
0023 #include <vector>
0024 #include <map>
0025 
0026 using namespace std;
0027 
0028 namespace {
0029   // Three implementations were tested: char-by-char (this version),
0030   // using std::string::find + std::string::replace and std::regex_replace.
0031   // First one takes ~60 ns per iteration, second one ~85 ns,
0032   // and the regex implementation takes nearly 1 us
0033   std::string globToRegex(const std::string& s) {
0034     std::string out;
0035     out.reserve(s.size());
0036     for (auto ch : s) {
0037       if (ch == '*') {
0038         out.push_back('.');
0039       }
0040       out.push_back(ch);
0041     }
0042     return out;
0043   }
0044 }  // namespace
0045 
0046 class DQMHistNormalizer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
0047 public:
0048   typedef dqm::legacy::DQMStore DQMStore;
0049   typedef dqm::legacy::MonitorElement MonitorElement;
0050 
0051   explicit DQMHistNormalizer(const edm::ParameterSet&);
0052   ~DQMHistNormalizer() override;
0053   void analyze(const edm::Event&, const edm::EventSetup&) override;
0054   void beginRun(const edm::Run& r, const edm::EventSetup& c) override {}
0055   void endRun(const edm::Run& r, const edm::EventSetup& c) override;
0056 
0057 private:
0058   vector<string> plotNamesToNormalize_;  //root name used by all the plots that must be normalized
0059   string reference_;
0060 };
0061 
0062 DQMHistNormalizer::DQMHistNormalizer(const edm::ParameterSet& cfg)
0063     : plotNamesToNormalize_(cfg.getParameter<std::vector<string> >("plotNamesToNormalize")),
0064       reference_(cfg.getParameter<string>("reference")) {
0065   usesResource("DQMStore");
0066   //std::cout << "<DQMHistNormalizer::DQMHistNormalizer>:" << std::endl;
0067 }
0068 
0069 DQMHistNormalizer::~DQMHistNormalizer() {
0070   //--- nothing to be done yet
0071 }
0072 
0073 void DQMHistNormalizer::analyze(const edm::Event&, const edm::EventSetup&) {
0074   //--- nothing to be done yet
0075 }
0076 
0077 void DQMHistNormalizer::endRun(const edm::Run& r, const edm::EventSetup& c) {
0078   //std::cout << "<DQMHistNormalizer::endJob>:" << std::endl;
0079 
0080   //--- check that DQMStore service is available
0081   if (!edm::Service<DQMStore>().isAvailable()) {
0082     edm::LogError("endJob") << " Failed to access dqmStore --> histograms will NOT be plotted !!";
0083     return;
0084   }
0085 
0086   DQMStore& dqmStore = (*edm::Service<DQMStore>());
0087 
0088   vector<MonitorElement*> allOurMEs = dqmStore.getAllContents("RecoTauV/");
0089   std::regex refregex = std::regex(".*RecoTauV/.*/" + globToRegex(reference_), std::regex::nosubs);
0090   vector<std::regex> toNormRegex;
0091   for (std::vector<string>::const_iterator toNorm = plotNamesToNormalize_.begin();
0092        toNorm != plotNamesToNormalize_.end();
0093        ++toNorm)
0094     toNormRegex.emplace_back(".*RecoTauV/.*/" + globToRegex(*toNorm), std::regex::nosubs);
0095 
0096   map<string, MonitorElement*> refsMap;
0097   vector<MonitorElement*> toNormElements;
0098   std::smatch path_match;
0099 
0100   for (vector<MonitorElement*>::const_iterator element = allOurMEs.begin(); element != allOurMEs.end(); ++element) {
0101     string pathname = (*element)->getFullname();
0102     //cout << pathname << endl;
0103     //Matches reference
0104     if (std::regex_match(pathname, path_match, refregex)) {
0105       //cout << "Matched to ref" << endl;
0106       string dir = pathname.substr(0, pathname.rfind('/'));
0107       if (refsMap.find(dir) != refsMap.end()) {
0108         edm::LogInfo("DQMHistNormalizer")
0109             << "DQMHistNormalizer::endRun: Warning! found multiple normalizing references for dir: " << dir << "!";
0110         edm::LogInfo("DQMHistNormalizer") << "     " << (refsMap[dir])->getFullname();
0111         edm::LogInfo("DQMHistNormalizer") << "     " << pathname;
0112       } else {
0113         refsMap[dir] = *element;
0114       }
0115     }
0116 
0117     //Matches targets
0118     for (const auto& reg : toNormRegex) {
0119       if (std::regex_match(pathname, path_match, reg)) {
0120         //cout << "Matched to target" << endl;
0121         toNormElements.push_back(*element);
0122         //cout << "Filled the collection" << endl;
0123       }
0124     }
0125   }
0126 
0127   for (vector<MonitorElement*>::const_iterator matchingElement = toNormElements.begin();
0128        matchingElement != toNormElements.end();
0129        ++matchingElement) {
0130     string meName = (*matchingElement)->getFullname();
0131     string dir = meName.substr(0, meName.rfind('/'));
0132 
0133     if (refsMap.find(dir) == refsMap.end()) {
0134       edm::LogInfo("DQMHistNormalizer") << "DQMHistNormalizer::endRun: Error! normalizing references for " << meName
0135                                         << " not found! Skipping...";
0136       continue;
0137     }
0138 
0139     float norm = refsMap[dir]->getTH1()->GetEntries();
0140     TH1* hist = (*matchingElement)->getTH1();
0141     if (norm != 0.) {
0142       if (!hist->GetSumw2N())
0143         hist->Sumw2();
0144       hist->Scale(1 / norm);  //use option "width" to divide the bin contents and errors by the bin width?
0145     } else {
0146       edm::LogInfo("DQMHistNormalizer") << "DQMHistNormalizer::endRun: Error! Normalization failed in "
0147                                         << hist->GetTitle() << "!";
0148     }
0149 
0150   }  //    for(vector<MonitorElement *>::const_iterator matchingElement = matchingElemts.begin(); matchingElement = matchingElemts.end(); ++matchingElement)
0151 }
0152 
0153 DEFINE_FWK_MODULE(DQMHistNormalizer);