CMS 3D CMS Logo

DQMHistNormalizer.cc
Go to the documentation of this file.
1 
8 // framework & common header files
15 
16 //DQM services
20 
21 //Regexp handling
22 #include "classlib/utils/RegexpMatch.h"
23 #include "classlib/utils/Regexp.h"
24 
25 #include <string>
26 #include <vector>
27 #include <map>
28 
29 using namespace std;
30 
32 public:
35 
36  explicit DQMHistNormalizer(const edm::ParameterSet&);
37  ~DQMHistNormalizer() override;
38  void analyze(const edm::Event&, const edm::EventSetup&) override;
39  void endJob() override {}
40  void endRun(const edm::Run& r, const edm::EventSetup& c) override;
41 
42 private:
43  lat::Regexp* buildRegex(const string& expr);
44  vector<string> plotNamesToNormalize_; //root name used by all the plots that must be normalized
45  string reference_;
46 };
47 
49  : plotNamesToNormalize_(cfg.getParameter<std::vector<string> >("plotNamesToNormalize")),
50  reference_(cfg.getParameter<string>("reference")) {
51  //std::cout << "<DQMHistNormalizer::DQMHistNormalizer>:" << std::endl;
52 }
53 
55  //--- nothing to be done yet
56 }
57 
59  //--- nothing to be done yet
60 }
61 
62 lat::Regexp* DQMHistNormalizer::buildRegex(const string& expr) {
63  lat::Regexp* rx = nullptr;
64  try {
65  rx = new lat::Regexp(expr, 0, lat::Regexp::Wildcard);
66  rx->study();
67  } catch (lat::Error& e) {
68  raiseDQMError("DQMStore", "Invalid regular expression '%s': %s", expr.c_str(), e.explain().c_str());
69  }
70  return rx;
71 }
72 
74  //std::cout << "<DQMHistNormalizer::endJob>:" << std::endl;
75 
76  //--- check that DQMStore service is available
77  if (!edm::Service<DQMStore>().isAvailable()) {
78  edm::LogError("endJob") << " Failed to access dqmStore --> histograms will NOT be plotted !!";
79  return;
80  }
81 
83 
84  vector<MonitorElement*> allOurMEs = dqmStore.getAllContents("RecoTauV/");
85  lat::Regexp* refregex = buildRegex("*RecoTauV/*/" + reference_);
86  vector<lat::Regexp*> toNormRegex;
87  for (std::vector<string>::const_iterator toNorm = plotNamesToNormalize_.begin();
88  toNorm != plotNamesToNormalize_.end();
89  ++toNorm)
90  toNormRegex.push_back(buildRegex("*RecoTauV/*/" + *toNorm));
91 
92  map<string, MonitorElement*> refsMap;
93  vector<MonitorElement*> toNormElements;
94 
95  for (vector<MonitorElement*>::const_iterator element = allOurMEs.begin(); element != allOurMEs.end(); ++element) {
96  string pathname = (*element)->getFullname();
97  //cout << pathname << endl;
98  //Matches reference
99  if (refregex->match(pathname)) {
100  //cout << "Matched to ref" << endl;
101  string dir = pathname.substr(0, pathname.rfind("/"));
102  if (refsMap.find(dir) != refsMap.end()) {
103  edm::LogInfo("DQMHistNormalizer")
104  << "DQMHistNormalizer::endRun: Warning! found multiple normalizing references for dir: " << dir << "!";
105  edm::LogInfo("DQMHistNormalizer") << " " << (refsMap[dir])->getFullname();
106  edm::LogInfo("DQMHistNormalizer") << " " << pathname;
107  } else {
108  refsMap[dir] = *element;
109  }
110  }
111 
112  //Matches targets
113  for (vector<lat::Regexp*>::const_iterator reg = toNormRegex.begin(); reg != toNormRegex.end(); ++reg) {
114  if ((*reg)->match(pathname)) {
115  //cout << "Matched to target" << endl;
116  toNormElements.push_back(*element);
117  //cout << "Filled the collection" << endl;
118  }
119  }
120  }
121 
122  delete refregex;
123  for (vector<lat::Regexp*>::const_iterator reg = toNormRegex.begin(); reg != toNormRegex.end(); ++reg)
124  delete *reg;
125 
126  for (vector<MonitorElement*>::const_iterator matchingElement = toNormElements.begin();
127  matchingElement != toNormElements.end();
128  ++matchingElement) {
129  string meName = (*matchingElement)->getFullname();
130  string dir = meName.substr(0, meName.rfind("/"));
131 
132  if (refsMap.find(dir) == refsMap.end()) {
133  edm::LogInfo("DQMHistNormalizer") << "DQMHistNormalizer::endRun: Error! normalizing references for " << meName
134  << " not found! Skipping...";
135  continue;
136  }
137 
138  float norm = refsMap[dir]->getTH1()->GetEntries();
139  TH1* hist = (*matchingElement)->getTH1();
140  if (norm != 0.) {
141  if (!hist->GetSumw2N())
142  hist->Sumw2();
143  hist->Scale(1 / norm); //use option "width" to divide the bin contents and errors by the bin width?
144  } else {
145  edm::LogInfo("DQMHistNormalizer") << "DQMHistNormalizer::endRun: Error! Normalization failed in "
146  << hist->GetTitle() << "!";
147  }
148 
149  } // for(vector<MonitorElement *>::const_iterator matchingElement = matchingElemts.begin(); matchingElement = matchingElemts.end(); ++matchingElement)
150 }
151 
edm::ErrorSummaryEntry Error
~DQMHistNormalizer() override
void endRun(const edm::Run &r, const edm::EventSetup &c) override
vector< string > plotNamesToNormalize_
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< MonitorElement * > getAllContents(std::string const &path, uint32_t runNumber=0, uint32_t lumi=0) const
Definition: DQMStore.cc:1616
lat::Regexp * buildRegex(const string &expr)
void endJob() override
void analyze(const edm::Event &, const edm::EventSetup &) override
DQMHistNormalizer(const edm::ParameterSet &)
dqm::legacy::MonitorElement MonitorElement
dqm::legacy::DQMStore DQMStore
Definition: Run.h:45
void raiseDQMError(const char *context, const char *fmt,...)
Definition: DQMError.cc:10