CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/DQM/HLTEvF/plugins/FourVectorHLT.cc

Go to the documentation of this file.
00001 // $Id: FourVectorHLT.cc,v 1.13 2010/02/17 22:53:25 wdd Exp $
00002 // See header file for information. 
00003 #include "FWCore/Framework/interface/EDAnalyzer.h"
00004 #include "DataFormats/Common/interface/Handle.h"
00005 #include "FWCore/Framework/interface/Run.h"
00006 #include "FWCore/Framework/interface/MakerMacros.h"
00007 
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "DQM/HLTEvF/interface/FourVectorHLT.h"
00010 
00011 #include "DataFormats/HLTReco/interface/TriggerObject.h"
00012 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00013 
00014 #include "DQMServices/Core/interface/MonitorElement.h"
00015 
00016 using namespace edm;
00017 
00018 FourVectorHLT::FourVectorHLT(const edm::ParameterSet& iConfig):
00019   resetMe_(true),  currentRun_(-99)
00020 {
00021   LogDebug("FourVectorHLT") << "constructor...." ;
00022 
00023   dbe_ = Service < DQMStore > ().operator->();
00024   if ( ! dbe_ ) {
00025     LogWarning("Status") << "unable to get DQMStore service?";
00026   }
00027   if (iConfig.getUntrackedParameter < bool > ("DQMStore", false)) {
00028     dbe_->setVerbose(0);
00029   }
00030   
00031   
00032   dirname_="HLT/FourVectorHLT" ;
00033   
00034   if (dbe_ != 0 ) {
00035     LogDebug("Status") << "Setting current directory to " << dirname_;
00036     dbe_->setCurrentFolder(dirname_);
00037   }
00038   
00039   
00040   // plotting paramters
00041   ptMin_ = iConfig.getUntrackedParameter<double>("ptMin",0.);
00042   ptMax_ = iConfig.getUntrackedParameter<double>("ptMax",200.);
00043   nBins_ = iConfig.getUntrackedParameter<unsigned int>("Nbins",50);
00044   
00045   plotAll_ = iConfig.getUntrackedParameter<bool>("plotAll", false);
00046 
00047   // this is the list of paths to look at.
00048   std::vector<edm::ParameterSet> filters = 
00049     iConfig.getParameter<std::vector<edm::ParameterSet> >("filters");
00050   for(std::vector<edm::ParameterSet>::iterator 
00051         filterconf = filters.begin() ; filterconf != filters.end(); 
00052       filterconf++) {
00053     std::string me  = filterconf->getParameter<std::string>("name");
00054     int objectType = filterconf->getParameter<unsigned int>("type");
00055     float ptMin = filterconf->getUntrackedParameter<double>("ptMin");
00056     float ptMax = filterconf->getUntrackedParameter<double>("ptMax");
00057     hltPaths_.push_back(PathInfo(me, objectType, ptMin, ptMax));
00058   }
00059   if ( hltPaths_.size() && plotAll_) {
00060     // these two ought to be mutually exclusive....
00061     LogWarning("Configuration") << "Using both plotAll and a list. "
00062       "list will be ignored." ;
00063     hltPaths_.clear();
00064   }
00065   triggerSummaryLabel_ = 
00066     iConfig.getParameter<edm::InputTag>("triggerSummaryLabel");
00067 }
00068 
00069 
00070 FourVectorHLT::~FourVectorHLT()
00071 {
00072  
00073    // do anything here that needs to be done at desctruction time
00074    // (e.g. close files, deallocate resources etc.)
00075 
00076 }
00077 
00078 
00079 //
00080 // member functions
00081 //
00082 
00083 // ------------ method called to for each event  ------------
00084 void
00085 FourVectorHLT::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00086 {
00087   using namespace edm;
00088   using namespace trigger;
00089   ++nev_;
00090   LogDebug("Status")<< "analyze" ;
00091   
00092   edm::Handle<TriggerEvent> triggerObj;
00093   iEvent.getByLabel(triggerSummaryLabel_,triggerObj); 
00094   if(!triggerObj.isValid()) { 
00095     edm::LogInfo("Status") << "Summary HLT object (TriggerEvent) not found, "
00096       "skipping event"; 
00097     return;
00098   }
00099   
00100   const trigger::TriggerObjectCollection & toc(triggerObj->getObjects());
00101 
00102   if ( plotAll_ ) {
00103     for ( size_t ia = 0; ia < triggerObj->sizeFilters(); ++ ia) {
00104       std::string fullname = triggerObj->filterTag(ia).encode();
00105       // the name can have in it the module label as well as the process and
00106       // other labels - strip 'em
00107       std::string name;
00108       size_t p = fullname.find_first_of(':');
00109       if ( p != std::string::npos) {
00110         name = fullname.substr(0, p);
00111       }
00112       else {
00113         name = fullname;
00114       }
00115 
00116       LogDebug("Parameter") << "filter " << ia << ", full name = " << fullname
00117                             << ", p = " << p 
00118                             << ", abbreviated = " << name ;
00119       
00120       PathInfoCollection::iterator pic =  hltPaths_.find(name);
00121       if ( pic == hltPaths_.end() ) {
00122         // doesn't exist - add it
00123         MonitorElement *et(0), *eta(0), *phi(0), *etavsphi(0);
00124         
00125         std::string histoname(name+"_et");
00126         LogDebug("Status") << "new histo with name "<<  histoname ;
00127         dbe_->setCurrentFolder(dirname_);
00128         std::string title(name+" E_{T}");
00129         et =  dbe_->book1D(histoname.c_str(),
00130                           title.c_str(),nBins_, 0, 100);
00131       
00132         histoname = name+"_eta";
00133         title = name+" #eta";
00134         eta =  dbe_->book1D(histoname.c_str(),
00135                            title.c_str(),nBins_,-2.7,2.7);
00136       
00137         histoname = name+"_phi";
00138         title = name+" #phi";
00139         phi =  dbe_->book1D(histoname.c_str(),
00140                            title.c_str(),nBins_,-3.14,3.14);
00141       
00142       
00143         histoname = name+"_etaphi";
00144         title = name+" #eta vs #phi";
00145         etavsphi =  dbe_->book2D(histoname.c_str(),
00146                                 title.c_str(),
00147                                 nBins_,-2.7,2.7,
00148                                 nBins_,-3.14, 3.14);
00149       
00150         // no idea how to get the bin boundries in this mode
00151         PathInfo e(name,0, et, eta, phi, etavsphi, ptMin_,ptMax_);
00152         hltPaths_.push_back(e);  
00153         pic = hltPaths_.begin() + hltPaths_.size()-1;
00154       }
00155       const trigger::Keys & k = triggerObj->filterKeys(ia);
00156       for (trigger::Keys::const_iterator ki = k.begin(); ki !=k.end(); ++ki ) {
00157         LogDebug("Parameters") << "pt, eta, phi = " 
00158                                << toc[*ki].pt() << ", "
00159                                << toc[*ki].eta() << ", "
00160                                << toc[*ki].phi() ;
00161         pic->getEtHisto()->Fill(toc[*ki].pt());
00162         pic->getEtaHisto()->Fill(toc[*ki].eta());
00163         pic->getPhiHisto()->Fill(toc[*ki].phi());
00164         pic->getEtaVsPhiHisto()->Fill(toc[*ki].eta(), toc[*ki].phi());
00165       }  
00166 
00167     }
00168 
00169   }
00170   else { // not plotAll_
00171     for(PathInfoCollection::iterator v = hltPaths_.begin();
00172         v!= hltPaths_.end(); ++v ) {
00173       const int index = triggerObj->filterIndex(v->getName());
00174       if ( index >= triggerObj->sizeFilters() ) {
00175         continue; // not in this event
00176       }
00177       LogDebug("Status") << "filling ... " ;
00178       const trigger::Keys & k = triggerObj->filterKeys(index);
00179       for (trigger::Keys::const_iterator ki = k.begin(); ki !=k.end(); ++ki ) {
00180         v->getEtHisto()->Fill(toc[*ki].pt());
00181         v->getEtaHisto()->Fill(toc[*ki].eta());
00182         v->getPhiHisto()->Fill(toc[*ki].phi());
00183         v->getEtaVsPhiHisto()->Fill(toc[*ki].eta(), toc[*ki].phi());
00184       }  
00185     }
00186   }
00187 }
00188 
00189 
00190 // -- method called once each job just before starting event loop  --------
00191 void 
00192 FourVectorHLT::beginJob()
00193 {
00194   nev_ = 0;
00195   DQMStore *dbe = 0;
00196   dbe = Service<DQMStore>().operator->();
00197   
00198   if (dbe) {
00199     dbe->setCurrentFolder(dirname_);
00200     dbe->rmdir(dirname_);
00201   }
00202   
00203   
00204   if (dbe) {
00205     dbe->setCurrentFolder(dirname_);
00206 
00207     if ( ! plotAll_ ) {
00208       for(PathInfoCollection::iterator v = hltPaths_.begin();
00209           v!= hltPaths_.end(); ++v ) {
00210         MonitorElement *et, *eta, *phi, *etavsphi=0;
00211         std::string histoname(v->getName()+"_et");
00212         std::string title(v->getName()+" E_t");
00213         et =  dbe->book1D(histoname.c_str(),
00214                           title.c_str(),nBins_,
00215                           v->getPtMin(),
00216                           v->getPtMax());
00217       
00218         histoname = v->getName()+"_eta";
00219         title = v->getName()+" #eta";
00220         eta =  dbe->book1D(histoname.c_str(),
00221                            title.c_str(),nBins_,-2.7,2.7);
00222 
00223         histoname = v->getName()+"_phi";
00224         title = v->getName()+" #phi";
00225         phi =  dbe->book1D(histoname.c_str(),
00226                            histoname.c_str(),nBins_,-3.14,3.14);
00227  
00228 
00229         histoname = v->getName()+"_etaphi";
00230         title = v->getName()+" #eta vs #phi";
00231         etavsphi =  dbe->book2D(histoname.c_str(),
00232                                 title.c_str(),
00233                                 nBins_,-2.7,2.7,
00234                                 nBins_,-3.14, 3.14);
00235       
00236         v->setHistos( et, eta, phi, etavsphi);
00237       } 
00238     } // ! plotAll_ - for plotAll we discover it during the event
00239   }
00240 }
00241 
00242 // - method called once each job just after ending the event loop  ------------
00243 void 
00244 FourVectorHLT::endJob() 
00245 {
00246    LogInfo("Status") << "endJob: analyzed " << nev_ << " events";
00247    return;
00248 }
00249 
00250 
00251 // BeginRun
00252 void FourVectorHLT::beginRun(const edm::Run& run, const edm::EventSetup& c)
00253 {
00254   LogDebug("Status") << "beginRun, run " << run.id();
00255 }
00256 
00258 void FourVectorHLT::endRun(const edm::Run& run, const edm::EventSetup& c)
00259 {
00260   LogDebug("Status") << "endRun, run " << run.id();
00261 }