CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/DQMOffline/PFTau/plugins/PFCandidateDQMAnalyzer.cc

Go to the documentation of this file.
00001 #include "DQMOffline/PFTau/plugins/PFCandidateDQMAnalyzer.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include "FWCore/Framework/interface/Event.h"
00006 #include "DataFormats/Common/interface/Handle.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 
00010 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00011 #include "DataFormats/Candidate/interface/Candidate.h"
00012 
00013 #include "DQMServices/Core/interface/DQMStore.h"
00014 #include "DQMServices/Core/interface/MonitorElement.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 //
00017 // -- Constructor
00018 //
00019 PFCandidateDQMAnalyzer::PFCandidateDQMAnalyzer(const edm::ParameterSet& parameterSet)  
00020   
00021 {
00022   pSet_                = parameterSet;
00023   inputLabel_          = pSet_.getParameter<edm::InputTag>("InputCollection");
00024   matchLabel_          = pSet_.getParameter<edm::InputTag>("MatchCollection");
00025   benchmarkLabel_      = pSet_.getParameter<std::string>("BenchmarkLabel"); 
00026 
00027   pfCandidateMonitor_.setParameters(parameterSet);  
00028   
00029 }
00030 //
00031 // -- BeginJob
00032 //
00033 void PFCandidateDQMAnalyzer::beginJob() {
00034 
00035   Benchmark::DQM_ = edm::Service<DQMStore>().operator->();
00036   // part of the following could be put in the base class
00037   std::string path = "ParticleFlow/" + benchmarkLabel_;
00038   Benchmark::DQM_->setCurrentFolder(path.c_str());
00039   edm::LogInfo("PFCandidateDQMAnalyzer") << " PFCandidateDQMAnalyzer::beginJob " << "Histogram Folder path set to "<< path;
00040   pfCandidateMonitor_.setup(pSet_);  
00041   nBadEvents_ = 0;
00042 }
00043 //
00044 // -- Analyze
00045 //
00046 void PFCandidateDQMAnalyzer::analyze(edm::Event const& iEvent, 
00047                                       edm::EventSetup const& iSetup) {
00048   edm::Handle< edm::View<reco::Candidate> > candCollection;
00049   iEvent.getByLabel( inputLabel_, candCollection);
00050 
00051   edm::Handle< edm::View<reco::Candidate> > matchedCandCollection;
00052   iEvent.getByLabel( matchLabel_, matchedCandCollection);
00053 
00054   float maxRes = 0.0;
00055   float minRes = 99.99;
00056   if (candCollection.isValid() && matchedCandCollection.isValid()) {
00057     pfCandidateMonitor_.fill( *candCollection, *matchedCandCollection, minRes, maxRes);
00058     edm::ParameterSet skimPS = pSet_.getParameter<edm::ParameterSet>("SkimParameter");
00059     if ( (skimPS.getParameter<bool>("switchOn")) &&  
00060          (nBadEvents_ <= skimPS.getParameter<int32_t>("maximumNumberToBeStored")) ) {
00061       if ( minRes < skimPS.getParameter<double>("lowerCutOffOnResolution")) {
00062         nBadEvents_++; 
00063         storeBadEvents(iEvent,minRes);
00064       } else if (maxRes > skimPS.getParameter<double>("upperCutOffOnResolution")) {
00065         nBadEvents_++;
00066         storeBadEvents(iEvent,maxRes);
00067       }
00068     }
00069   }
00070 }
00071 void PFCandidateDQMAnalyzer::storeBadEvents(edm::Event const& iEvent, float& val) {
00072   unsigned int runNb  = iEvent.id().run();
00073   unsigned int evtNb  = iEvent.id().event();
00074   unsigned int lumiNb = iEvent.id().luminosityBlock();
00075   
00076   std::string path = "ParticleFlow/" + benchmarkLabel_ + "/BadEvents";
00077   Benchmark::DQM_->setCurrentFolder(path.c_str());
00078   std::ostringstream eventid_str;
00079   eventid_str << runNb << "_"<< evtNb << "_" << lumiNb;
00080   MonitorElement* me = Benchmark::DQM_->get(path + "/" + eventid_str.str());
00081   if (me) me->Reset();
00082   else me = Benchmark::DQM_->bookFloat(eventid_str.str());
00083   me->Fill(val);  
00084 }
00085 //
00086 // -- EndJob
00087 // 
00088 void PFCandidateDQMAnalyzer::endJob() {
00089 }
00090 #include "FWCore/Framework/interface/MakerMacros.h"
00091 DEFINE_FWK_MODULE (PFCandidateDQMAnalyzer) ;