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