CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DQMOffline/Trigger/src/TopElectronHLTOfflineSource.cc

Go to the documentation of this file.
00001 #include "DQMOffline/Trigger/interface/TopElectronHLTOfflineSource.h"
00002 
00003 #include "DQMServices/Core/interface/MonitorElement.h"
00004 
00005 #include "FWCore/Framework/interface/Run.h"
00006 
00007 #include <boost/algorithm/string.hpp>
00008 
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
00011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00012 #include "DataFormats/Common/interface/TriggerResults.h"
00013 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00014 #include "DQMOffline/Trigger/interface/TopElectronHLTOfflineSource.h"
00015 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00017 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00018 #include "DataFormats/Common/interface/ValueMap.h"
00019 #include "DataFormats/Math/interface/deltaR.h"
00020 
00021 #include "FWCore/ServiceRegistry/interface/Service.h"
00022 
00023 //using namespace egHLT;
00024 
00025 TopElectronHLTOfflineSource::TopElectronHLTOfflineSource(const edm::ParameterSet& conf) :
00026   beamSpot_(conf.getParameter<edm::InputTag>("beamSpot")) {
00027 
00028         dbe_ = edm::Service<DQMStore>().operator->();
00029         
00030         if (!dbe_) 
00031         {
00032                 edm::LogInfo("TopElectronHLTOfflineSource") << "unable to get DQMStore service?";
00033         }
00034         
00035         if(conf.getUntrackedParameter<bool>("DQMStore", false)) 
00036         {
00037           if(!dbe_) dbe_->setVerbose(0);
00038         }
00039 
00040         dirName_ = conf.getParameter<std::string>("DQMDirName");
00041         
00042         electronIdNames_ = conf.getParameter<std::vector<std::string> >("electronIdNames");
00043         hltTag_ = conf.getParameter<std::string>("hltTag");
00044         superTriggerNames_ = conf.getParameter<std::vector<std::string> >("superTriggerNames");
00045         electronTriggerNames_ = conf.getParameter<std::vector<std::string> >("electronTriggerNames");
00046         
00047         triggerResultsLabel_ = conf.getParameter<edm::InputTag>("triggerResultsLabel");
00048         triggerSummaryLabel_ = conf.getParameter<edm::InputTag>("triggerSummaryLabel");
00049         electronLabel_ = conf.getParameter<edm::InputTag>("electronCollection");
00050         primaryVertexLabel_ = conf.getParameter<edm::InputTag>("primaryVertexCollection");
00051         triggerJetFilterLabel_  = conf.getParameter<edm::InputTag>("triggerJetFilterLabel");
00052         triggerElectronFilterLabel_ = conf.getParameter<edm::InputTag>("triggerElectronFilterLabel");
00053 
00054         excludeCloseJets_ = conf.getParameter<bool>("excludeCloseJets");
00055         requireTriggerMatch_ = conf.getParameter<bool>("requireTriggerMatch");
00056         electronMinEt_ = conf.getParameter<double>("electronMinEt");
00057         electronMaxEta_ = conf.getParameter<double>("electronMaxEta");
00058                 
00059         addExtraId_ = conf.getParameter<bool>("addExtraId");
00060         extraIdCutsSigmaEta_ = conf.getParameter<double>("extraIdCutsSigmaEta");
00061         extraIdCutsSigmaPhi_ = conf.getParameter<double>("extraIdCutsSigmaPhi");
00062         extraIdCutsDzPV_ = conf.getParameter<double>("extraIdCutsDzPV");
00063 
00064 }
00065 TopElectronHLTOfflineSource::~TopElectronHLTOfflineSource()
00066 {
00067 }
00068 
00069 void TopElectronHLTOfflineSource::beginJob()
00070 {
00071   if(!dbe_) return;
00072         dbe_->setCurrentFolder(dirName_);
00073         for (size_t i = 0; i < superTriggerNames_.size(); ++i)
00074         {
00075                 eleMEs_.push_back(EleMEs(dbe_, electronIdNames_, addExtraId_, superTriggerNames_[i]));
00076                 for (size_t j = 0; j < electronTriggerNames_.size(); ++j)
00077                 {
00078                         eleMEs_.push_back(EleMEs(dbe_, electronIdNames_, addExtraId_, superTriggerNames_[i]+"_"+electronTriggerNames_[j]));
00079                         //std::cout <<superTriggerNames_[i]+"_"+electronTriggerNames_[j]<<std::endl;
00080                         
00081                 }
00082         }
00083         //std::cout <<"done"<<std::endl;
00084 }
00085 void TopElectronHLTOfflineSource::setupHistos(std::vector<EleMEs> topEleHists)
00086 {
00087         for (size_t i = 0; i < eleMEs_.size(); ++i)
00088         {
00089                 topEleHists.push_back(eleMEs_[i]);
00090         }
00091 }
00092 void TopElectronHLTOfflineSource::endJob()
00093 {
00094 }
00095 void TopElectronHLTOfflineSource::beginRun(const edm::Run& run, const edm::EventSetup& c)
00096 {
00097   hltConfigValid_=hltConfig_.init(run,c,hltTag_,hltConfigChanged_);
00098 }
00099 void TopElectronHLTOfflineSource::endRun(const edm::Run& run, const edm::EventSetup& c)
00100 {
00101 }
00102 
00103 void TopElectronHLTOfflineSource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00104 {
00105         if(!dbe_) return;
00106         // ---- Get Trigger Decisions for all triggers under investigation ----
00107         edm::Handle<edm::TriggerResults> hltResults;
00108         if(!iEvent.getByLabel(triggerResultsLabel_, hltResults) || !hltResults.product()) return; //bail if we didnt get trigger results
00109       
00110 
00111         
00112         if (!hltConfigValid_) return;
00113         
00114         std::vector<bool> superTriggerAccepts;
00115         std::vector<bool> electronTriggerAccepts;
00116         
00117         for (size_t i = 0; i < superTriggerNames_.size(); ++i)
00118         {
00119                 unsigned int triggerIndex( hltConfig_.triggerIndex(superTriggerNames_[i]) );
00120                 bool accept = false;
00121                 
00122                 if (triggerIndex < hltResults->size())
00123                 {
00124                         accept = hltResults->accept(triggerIndex);
00125                 }
00126                 
00127                 superTriggerAccepts.push_back(accept);
00128         }
00129         
00130         for (size_t i = 0; i < electronTriggerNames_.size(); ++i)
00131         {
00132                 unsigned int triggerIndex( hltConfig_.triggerIndex(electronTriggerNames_[i]) );
00133                 bool accept = false;
00134                 
00135                 if (triggerIndex < hltResults->size())
00136                 {
00137                         accept = hltResults->accept(triggerIndex);
00138                 } 
00139                 
00140                 electronTriggerAccepts.push_back(accept);
00141         }
00142         
00143         // get reconstructed electron collection
00144         if(!iEvent.getByLabel(electronLabel_, eleHandle_) || !eleHandle_.product()) return;
00145         
00146         // Get Trigger Event, providing the information about trigger objects   
00147         if(!iEvent.getByLabel(triggerSummaryLabel_, triggerEvent_) || !triggerEvent_.product()) return; 
00148         
00149         edm::Handle<reco::VertexCollection> vertexHandle;
00150         if(!iEvent.getByLabel(primaryVertexLabel_, vertexHandle) || !vertexHandle.product()) return;
00151 
00152         reco::Vertex::Point vertexPoint(0., 0., 0.);
00153         if (vertexHandle.product()->size() != 0)
00154         {
00155                 const reco::Vertex& theVertex = vertexHandle.product()->front();
00156                 vertexPoint = theVertex.position();
00157         }
00158         else
00159         {
00160                 edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00161                 if(!iEvent.getByLabel(beamSpot_, recoBeamSpotHandle) ||  !recoBeamSpotHandle.product()) return;
00162 
00163                 vertexPoint = recoBeamSpotHandle->position();
00164         }
00165         
00166         trigger::size_type jetFilterPos = triggerEvent_->filterIndex(triggerJetFilterLabel_);
00167         std::vector<const trigger::TriggerObject*> triggerJets;
00168         
00169         if (jetFilterPos != triggerEvent_->sizeFilters())
00170         {
00171                 for (size_t i = 0; i < triggerEvent_->filterKeys(jetFilterPos).size(); ++i)
00172                 {
00173                         size_t objNr = triggerEvent_->filterKeys(jetFilterPos)[i];
00174                         if(objNr<triggerEvent_->sizeObjects()){
00175                           triggerJets.push_back(& triggerEvent_->getObjects()[objNr]);
00176                         }
00177                 }       
00178         }
00179         
00180         trigger::size_type eleFilterPos = triggerEvent_->filterIndex(triggerElectronFilterLabel_);
00181         std::vector<const trigger::TriggerObject*> triggerElectrons;
00182  
00183         if (eleFilterPos != triggerEvent_->sizeFilters())
00184         {
00185                 for (size_t i = 0; i < triggerEvent_->filterKeys(eleFilterPos).size(); ++i)
00186                 {
00187                         size_t objNr = triggerEvent_->filterKeys(eleFilterPos)[i];
00188                         if(objNr<triggerEvent_->sizeObjects()){
00189                           triggerElectrons.push_back(& triggerEvent_->getObjects()[objNr]);
00190                         }
00191                 }
00192         }
00193         
00194         const reco::GsfElectronCollection& eles = *eleHandle_;
00195         
00196         for(size_t eleNr=0; eleNr < eles.size(); ++eleNr)
00197         {
00198                 
00199                 const reco::GsfElectron& ele = eles[eleNr];
00200                 
00201                 // electron selection
00202                 
00203                 if(ele.et() > electronMinEt_ && std::abs(ele.eta()) < electronMaxEta_)
00204                 {
00205                         size_t index = 0;
00206                         for (size_t i = 0; i < superTriggerNames_.size(); ++i)
00207                         {
00208                                 if (superTriggerAccepts[i])
00209                                         fill(eleMEs_[index], iEvent, eleNr, triggerJets, triggerElectrons, vertexPoint);
00210                                 index++;
00211                                 
00212                                 for (size_t j = 0; j < electronTriggerNames_.size(); ++j)
00213                                 {
00214                                         if (superTriggerAccepts[i] && electronTriggerAccepts[j]) 
00215                                                 fill(eleMEs_[index], iEvent, eleNr, triggerJets, triggerElectrons, vertexPoint);
00216                                         index++;
00217                                 }
00218                         }
00219                 }
00220         }
00221 }
00222 
00223 void TopElectronHLTOfflineSource::EleMEs::setup(DQMStore* dbe, const std::vector<std::string>& eleIdNames, bool addExtraId, const std::string& name)
00224 { 
00225         for (size_t i = 0; i < eleIdNames.size(); ++i)
00226         {
00227                 eleIdNames_.push_back(eleIdNames[i]);
00228                 if (addExtraId)
00229                         eleIdNames_.push_back(eleIdNames[i]+"extraId");
00230         }
00231         
00232         addMESets(name);
00233         
00234         for (size_t i = 0; i < eleMESets_.size(); ++i)
00235         {
00236                 setupMESet(eleMESets_[i], dbe, name+"_"+fullName(i));
00237                 LogDebug("TopElectronHLTOfflineSource") << "Booked MonitorElement with name " << name;
00238         }
00239 }
00240 void TopElectronHLTOfflineSource::EleMEs::setupMESet(EleMESet& eleSet, DQMStore* dbe, const std::string& name)
00241 {
00242         eleSet.ele_et = dbe->book1D("ele_"+name+"_et", "ele_"+name+"_et", 50, 0., 500.);
00243         eleSet.ele_eta = dbe->book1D("ele_"+name+"_eta", "ele_"+name+"_eta", 50, -2.5, 2.5);
00244         eleSet.ele_phi = dbe->book1D("ele_"+name+"_phi","ele_"+name+"_phi", 50, -3.1416, 3.1416);
00245         eleSet.ele_isolEm = dbe->book1D("ele_"+name+"_isolEm", "ele_"+name+"_isolEm", 50, -0.05, 3.);
00246         eleSet.ele_isolHad = dbe->book1D("ele_"+name+"_isolHad", "ele_"+name+"_isolHad", 50, -0.05, 5.);
00247         eleSet.ele_minDeltaR = dbe->book1D("ele_"+name+"_minDeltaR", "ele_"+name+"_minDeltaR", 50, 0., 1.);
00248         eleSet.global_n30jets = dbe->book1D("ele_"+name+"_global_n30jets", "ele_"+name+"_global_n30jets", 10, -0.5, 9.5);
00249         eleSet.global_sumEt = dbe->book1D("ele_"+name+"_global_sumEt", "ele_"+name+"_global_sumEt", 50, 0., 1000.);
00250         eleSet.ele_gsftrack_etaError = dbe->book1D("ele_"+name+"_gsftrack_etaError", "ele_"+name+"_gsftrack_etaError", 50, 0., 0.005);
00251         eleSet.ele_gsftrack_phiError = dbe->book1D("ele_"+name+"_gsftrack_phiError", "ele_"+name+"_gsftrack_phiError", 50, 0., 0.005);
00252         eleSet.ele_gsftrack_numberOfValidHits = dbe->book1D("ele_"+name+"_gsftrack_numberOfValidHits", "ele_"+name+"_gsftrack_numberOfValidHits", 25, -0.5, 24.5);
00253         eleSet.ele_gsftrack_dzPV = dbe->book1D("ele_"+name+"_gsftrack_dzPV", "ele_"+name+"_gsftrack_dzPV", 50, 0., 0.2);
00254 }
00255 
00256 void TopElectronHLTOfflineSource::EleMEs::addMESets(const std::string& name)
00257 {
00258         eleMENames_.push_back("EB");
00259         eleMENames_.push_back("EE");
00260         name_ = name;
00261         for (size_t i=0; i < eleIdNames_.size() * eleMENames_.size(); ++i)
00262         {
00263                 eleMESets_.push_back(EleMESet());
00264         }
00265 }
00266 
00267 void TopElectronHLTOfflineSource::fill(EleMEs& eleMEs, const edm::Event& iEvent, size_t eleIndex, const std::vector<const trigger::TriggerObject*>& triggerJets, const std::vector<const trigger::TriggerObject*>& triggerElectrons, const reco::Vertex::Point& vertexPoint)
00268 {
00269         const reco::GsfElectron& ele = (*eleHandle_)[eleIndex];
00270         
00271         float dzPV = std::abs(ele.gsfTrack()->dz(vertexPoint));
00272         
00273         bool isTriggerMatched = false;
00274         for (size_t i = 0; i < triggerElectrons.size(); ++i)
00275         {
00276                 if (deltaR(*(triggerElectrons[i]), ele.p4()) < 0.3)
00277                          isTriggerMatched = true;
00278         }
00279         
00280         if (requireTriggerMatch_ && !isTriggerMatched)
00281                 return;
00282         
00283         // Calculate minimum deltaR to closest jet and sumEt (all jets)
00284         float minDeltaR = 999.;
00285         float sumEt = 0.;
00286         
00287         for (size_t jetNr = 0; jetNr < triggerJets.size(); ++jetNr)
00288         {
00289                 const trigger::TriggerObject& jet = *triggerJets[jetNr];
00290 
00291                 sumEt += jet.et();
00292                 
00293                 float dr = deltaR(jet, ele.p4());
00294                 
00295                 if (!excludeCloseJets_ && dr < minDeltaR)
00296                         minDeltaR = dr;
00297                 if (excludeCloseJets_ && dr > 0.1 && dr < minDeltaR)
00298                         minDeltaR = dr;
00299         }
00300         
00301         for (size_t j = 0; j < eleMEs.eleIdNames().size(); ++j)
00302         {
00303                 bool eId = true;
00304 
00305                 edm::Handle<edm::ValueMap<float> > eIdMapHandle;
00306                 iEvent.getByLabel(eleMEs.eleIdNames()[j], eIdMapHandle);
00307                 const edm::ValueMap<float>& eIdMap = *eIdMapHandle;
00308                 eId = eIdMap[edm::Ref<reco::GsfElectronCollection>(eleHandle_, eleIndex)];
00309                 
00310                 bool extraId = true;
00311                 if (addExtraId_)
00312                 {
00313                         if (ele.gsfTrack()->etaError() > extraIdCutsSigmaEta_)
00314                                 extraId = false;
00315                         if (ele.gsfTrack()->phiError() > extraIdCutsSigmaPhi_)
00316                                 extraId = false;
00317                         if (dzPV > extraIdCutsDzPV_)
00318                                 extraId = false;
00319                 }
00320                 
00321                 for (size_t i = 0; i < eleMEs.eleMENames().size(); ++i)
00322                 {
00323                         if (eId && eleMEs.eleMENames()[i] == "EB" && ele.isEB()&& !ele.isGap())
00324                                 eleMEs.fill(eleMEs.getMESet(i, j), ele, minDeltaR, sumEt, triggerJets.size(), dzPV);
00325                         if (eId && eleMEs.eleMENames()[i] == "EE" && ele.isEE()&& !ele.isGap())
00326                                 eleMEs.fill(eleMEs.getMESet(i, j), ele, minDeltaR, sumEt, triggerJets.size(), dzPV);
00327                         if (addExtraId_)
00328                         {
00329                                 if (eId && extraId && eleMEs.eleMENames()[i] == "EB" && ele.isEB()&& !ele.isGap())
00330                                         eleMEs.fill(eleMEs.getMESet(i, j+1), ele, minDeltaR, sumEt, triggerJets.size(), dzPV);
00331                                 if (eId && extraId && eleMEs.eleMENames()[i] == "EE" && ele.isEE()&& !ele.isGap())
00332                                         eleMEs.fill(eleMEs.getMESet(i, j+1), ele, minDeltaR, sumEt, triggerJets.size(), dzPV);
00333                         }
00334                 }
00335                 if (addExtraId_)
00336                         ++j;
00337 
00338         }
00339 }
00340 void TopElectronHLTOfflineSource::EleMEs::fill(EleMESet& eleMESet, const reco::GsfElectron& ele, float minDeltaR, float sumEt, int n30jets, float dzPV)
00341 {
00342         LogDebug("TopElectronHLTOfflineSource") << "filling the histos with " << ele.et();
00343 
00344         eleMESet.ele_et->Fill(ele.et());
00345         eleMESet.ele_eta->Fill(ele.eta());
00346         eleMESet.ele_phi->Fill(ele.phi());
00347         eleMESet.ele_isolEm->Fill(ele.dr03EcalRecHitSumEt());
00348         eleMESet.ele_isolHad->Fill(ele.dr03HcalTowerSumEt());
00349         eleMESet.ele_minDeltaR->Fill(minDeltaR);
00350         eleMESet.global_n30jets->Fill(n30jets);
00351         eleMESet.global_sumEt->Fill(sumEt);
00352         eleMESet.ele_gsftrack_etaError->Fill(ele.gsfTrack()->etaError());
00353         eleMESet.ele_gsftrack_phiError->Fill(ele.gsfTrack()->phiError());
00354         eleMESet.ele_gsftrack_numberOfValidHits->Fill(ele.gsfTrack()->numberOfValidHits());
00355         eleMESet.ele_gsftrack_dzPV->Fill(dzPV);         
00356 }