CMS 3D CMS Logo

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