Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023
00024
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "HLTrigger/HLTfilters/interface/JetCollectionForEleHT.h"
00033
00034 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00035 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00036 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00037
00038 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidateIsolation.h"
00039 #include "DataFormats/Common/interface/AssociationMap.h"
00040 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00041
00042 #include "DataFormats/Common/interface/Handle.h"
00043
00044 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00045
00046 #include "DataFormats/HLTReco/interface/TriggerObject.h"
00047
00048 JetCollectionForEleHT::JetCollectionForEleHT(const edm::ParameterSet& iConfig):
00049 hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
00050 sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
00051 minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
00052 {
00053 produces<reco::CaloJetCollection>();
00054 }
00055
00056
00057
00058 JetCollectionForEleHT::~JetCollectionForEleHT()
00059 {
00060
00061
00062
00063
00064 }
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074 void
00075 JetCollectionForEleHT::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00076 {
00077 using namespace edm;
00078
00079 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
00080 iEvent.getByLabel(hltElectronTag,PrevFilterOutput);
00081
00082
00083 std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
00084 PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
00085 std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
00086 PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
00087
00088
00089 std::vector<TVector3> ElePs;
00090
00091 if(!clusCands.empty()){
00092 for(size_t candNr=0;candNr<clusCands.size();candNr++){
00093 TVector3 positionVector(
00094 clusCands[candNr]->superCluster()->position().x(),
00095 clusCands[candNr]->superCluster()->position().y(),
00096 clusCands[candNr]->superCluster()->position().z());
00097 ElePs.push_back(positionVector);
00098 }
00099 }else if(!eleCands.empty()){
00100 for(size_t candNr=0;candNr<eleCands.size();candNr++){
00101 TVector3 positionVector(
00102 eleCands[candNr]->superCluster()->position().x(),
00103 eleCands[candNr]->superCluster()->position().y(),
00104 eleCands[candNr]->superCluster()->position().z());
00105 ElePs.push_back(positionVector);
00106 }
00107 }
00108
00109 edm::Handle<reco::CaloJetCollection> theCaloJetCollectionHandle;
00110 iEvent.getByLabel(sourceJetTag, theCaloJetCollectionHandle);
00111 const reco::CaloJetCollection* theCaloJetCollection = theCaloJetCollectionHandle.product();
00112
00113 std::auto_ptr< reco::CaloJetCollection > theFilteredCaloJetCollection(new reco::CaloJetCollection);
00114
00115 bool isOverlapping;
00116
00117 for(unsigned int j=0; j<theCaloJetCollection->size(); j++) {
00118
00119 isOverlapping = false;
00120 for(unsigned int i=0; i<ElePs.size(); i++) {
00121
00122 TVector3 JetP((*theCaloJetCollection)[j].px(), (*theCaloJetCollection)[j].py(), (*theCaloJetCollection)[j].pz());
00123 double DR = ElePs[i].DeltaR(JetP);
00124
00125 if(DR<minDeltaR_) {
00126 isOverlapping = true;
00127 break;
00128 }
00129 }
00130
00131 if(!isOverlapping) theFilteredCaloJetCollection->push_back((*theCaloJetCollection)[j]);
00132 }
00133
00134
00135
00136 iEvent.put(theFilteredCaloJetCollection);
00137
00138 return;
00139
00140 }
00141
00142
00143 void
00144 JetCollectionForEleHT::beginJob()
00145 {
00146 }
00147
00148
00149 void
00150 JetCollectionForEleHT::endJob() {
00151 }
00152
00153
00154 DEFINE_FWK_MODULE(JetCollectionForEleHT);
00155