CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/Validation/RecoTau/plugins/DQMHistNormalizer.cc

Go to the documentation of this file.
00001 
00010 // framework & common header files
00011 #include "FWCore/Framework/interface/EDAnalyzer.h"
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/MakerMacros.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 
00018 
00019 //DQM services
00020 #include "DQMServices/Core/interface/DQMStore.h"
00021 #include "FWCore/ServiceRegistry/interface/Service.h"
00022 #include "DQMServices/Core/interface/MonitorElement.h"
00023 #include "DQMServices/Core/src/DQMError.h"
00024 
00025 
00026 //Regexp handling
00027 #include "classlib/utils/RegexpMatch.h"
00028 #include "classlib/utils/Regexp.h"
00029 
00030 
00031 #include <string>
00032 #include <vector>
00033 #include <map>
00034 
00035 using namespace std;
00036 
00037 
00038 class DQMHistNormalizer : public edm::EDAnalyzer
00039 {
00040 
00041  public:
00042   explicit DQMHistNormalizer(const edm::ParameterSet&);
00043   virtual ~DQMHistNormalizer();
00044   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00045   virtual void endJob(){}
00046   virtual void endRun(const edm::Run& r, const edm::EventSetup& c);
00047 
00048 private:
00049   lat::Regexp* buildRegex(const string & expr);
00050   vector<string> plotNamesToNormalize_; //root name used by all the plots that must be normalized
00051   string reference_;
00052 };
00053 
00054 DQMHistNormalizer::DQMHistNormalizer(const edm::ParameterSet& cfg):
00055   plotNamesToNormalize_(cfg.getParameter< std::vector<string> >("plotNamesToNormalize")),
00056   reference_(cfg.getParameter< string >("reference"))
00057 {
00058   //std::cout << "<DQMHistNormalizer::DQMHistNormalizer>:" << std::endl;
00059 }
00060 
00061 DQMHistNormalizer::~DQMHistNormalizer() 
00062 {
00063 //--- nothing to be done yet
00064 }
00065 
00066 void DQMHistNormalizer::analyze(const edm::Event&, const edm::EventSetup&)
00067 {
00068 //--- nothing to be done yet
00069 }
00070 
00071 lat::Regexp* DQMHistNormalizer::buildRegex(const string & expr)
00072 {
00073   lat::Regexp* rx = 0;
00074   try
00075   {
00076     rx = new lat::Regexp(expr, 0, lat::Regexp::Wildcard);
00077     rx->study();
00078   }
00079   catch (lat::Error &e)
00080   {
00081     raiseDQMError("DQMStore", "Invalid regular expression '%s': %s",
00082                   expr.c_str(), e.explain().c_str());
00083   }
00084   return rx;
00085 }
00086 
00087 void DQMHistNormalizer::endRun(const edm::Run& r, const edm::EventSetup& c)
00088 {
00089   //std::cout << "<DQMHistNormalizer::endJob>:" << std::endl;
00090 
00091 //--- check that DQMStore service is available
00092   if ( !edm::Service<DQMStore>().isAvailable() ) {
00093     edm::LogError ("endJob") << " Failed to access dqmStore --> histograms will NOT be plotted !!";
00094     return;
00095   }
00096 
00097   DQMStore& dqmStore = (*edm::Service<DQMStore>());
00098 
00099   vector<MonitorElement *> allOurMEs = dqmStore.getAllContents("RecoTauV/");
00100   lat::Regexp* refregex = buildRegex("*RecoTauV/*/" + reference_);
00101   vector<lat::Regexp*> toNormRegex;
00102   for ( std::vector<string>::const_iterator toNorm = plotNamesToNormalize_.begin(); toNorm != plotNamesToNormalize_.end(); ++toNorm )
00103     toNormRegex.push_back( buildRegex("*RecoTauV/*/" + *toNorm) );
00104 
00105   map<string, MonitorElement *> refsMap;
00106   vector<MonitorElement *> toNormElements;
00107 
00108   for(vector<MonitorElement *>::const_iterator element = allOurMEs.begin(); element != allOurMEs.end(); ++element){
00109     string pathname = (*element)->getFullname();
00110     //cout << pathname << endl;
00111     //Matches reference
00112     if(refregex->match(pathname)){
00113       //cout << "Matched to ref" << endl;
00114       string dir = pathname.substr(0, pathname.rfind("/"));
00115       if(refsMap.find(dir) != refsMap.end()){
00116         edm::LogInfo("DQMHistNormalizer")<<"DQMHistNormalizer::endRun: Warning! found multiple normalizing references for dir: "<<dir<<"!";
00117         edm::LogInfo("DQMHistNormalizer")<<"     " << (refsMap[dir])->getFullname();
00118         edm::LogInfo("DQMHistNormalizer")<<"     " << pathname;
00119       }
00120       else{
00121         refsMap[dir] = *element;
00122       }
00123     }
00124 
00125     //Matches targets
00126     for(vector<lat::Regexp*>::const_iterator reg = toNormRegex.begin(); reg != toNormRegex.end(); ++reg){
00127       if((*reg)->match(pathname)){
00128         //cout << "Matched to target" << endl;
00129         toNormElements.push_back(*element);
00130         //cout << "Filled the collection" << endl;
00131       }
00132     }
00133   }
00134 
00135   delete refregex;
00136   for(vector<lat::Regexp*>::const_iterator reg = toNormRegex.begin(); reg != toNormRegex.end(); ++reg)
00137     delete *reg;
00138 
00139   for(vector<MonitorElement *>::const_iterator matchingElement = toNormElements.begin(); matchingElement != toNormElements.end(); ++matchingElement){
00140     string meName = (*matchingElement)->getFullname();
00141     string dir = meName.substr(0, meName.rfind("/"));
00142     
00143     if(refsMap.find(dir) == refsMap.end()){
00144       edm::LogInfo("DQMHistNormalizer")<<"DQMHistNormalizer::endRun: Error! normalizing references for "<<meName<<" not found! Skipping...";
00145       continue;
00146     }
00147       
00148     float norm = refsMap[dir]->getTH1()->GetEntries();
00149     TH1* hist = (*matchingElement)->getTH1();
00150     if ( norm != 0. ) {
00151       if( !hist->GetSumw2N() ) hist->Sumw2();
00152       hist->Scale(1/norm);//use option "width" to divide the bin contents and errors by the bin width?
00153     }else{
00154       edm::LogInfo("DQMHistNormalizer")<<"DQMHistNormalizer::endRun: Error! Normalization failed in "<<hist->GetTitle()<<"!";
00155     }
00156       
00157   }//    for(vector<MonitorElement *>::const_iterator matchingElement = matchingElemts.begin(); matchingElement = matchingElemts.end(); ++matchingElement)
00158 }
00159 
00160 DEFINE_FWK_MODULE(DQMHistNormalizer);
00161 
00162 
00163