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