CMS 3D CMS Logo

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