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 void
00067 JetCollectionForEleHT::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00068 edm::ParameterSetDescription desc;
00069 desc.add<edm::InputTag>("HltElectronTag",edm::InputTag("triggerFilterObjectWithRefs"));
00070 desc.add<edm::InputTag>("SourceJetTag",edm::InputTag("caloJetCollection"));
00071 desc.add<double>("minDeltaR",0.5);
00072 descriptions.add("hltJetCollectionForEleHT",desc);
00073 }
00074
00075
00076
00077
00078
00079
00080
00081
00082 void
00083 JetCollectionForEleHT::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00084 {
00085 using namespace edm;
00086
00087 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
00088 iEvent.getByLabel(hltElectronTag,PrevFilterOutput);
00089
00090
00091 std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
00092 PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
00093 std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
00094 PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
00095
00096
00097 std::vector<TVector3> ElePs;
00098
00099 if(!clusCands.empty()){
00100 for(size_t candNr=0;candNr<clusCands.size();candNr++){
00101 TVector3 positionVector(
00102 clusCands[candNr]->superCluster()->position().x(),
00103 clusCands[candNr]->superCluster()->position().y(),
00104 clusCands[candNr]->superCluster()->position().z());
00105 ElePs.push_back(positionVector);
00106 }
00107 }else if(!eleCands.empty()){
00108 for(size_t candNr=0;candNr<eleCands.size();candNr++){
00109 TVector3 positionVector(
00110 eleCands[candNr]->superCluster()->position().x(),
00111 eleCands[candNr]->superCluster()->position().y(),
00112 eleCands[candNr]->superCluster()->position().z());
00113 ElePs.push_back(positionVector);
00114 }
00115 }
00116
00117 edm::Handle<reco::CaloJetCollection> theCaloJetCollectionHandle;
00118 iEvent.getByLabel(sourceJetTag, theCaloJetCollectionHandle);
00119 const reco::CaloJetCollection* theCaloJetCollection = theCaloJetCollectionHandle.product();
00120
00121 std::auto_ptr< reco::CaloJetCollection > theFilteredCaloJetCollection(new reco::CaloJetCollection);
00122
00123 bool isOverlapping;
00124
00125 for(unsigned int j=0; j<theCaloJetCollection->size(); j++) {
00126
00127 isOverlapping = false;
00128 for(unsigned int i=0; i<ElePs.size(); i++) {
00129
00130 TVector3 JetP((*theCaloJetCollection)[j].px(), (*theCaloJetCollection)[j].py(), (*theCaloJetCollection)[j].pz());
00131 double DR = ElePs[i].DeltaR(JetP);
00132
00133 if(DR<minDeltaR_) {
00134 isOverlapping = true;
00135 break;
00136 }
00137 }
00138
00139 if(!isOverlapping) theFilteredCaloJetCollection->push_back((*theCaloJetCollection)[j]);
00140 }
00141
00142
00143
00144 iEvent.put(theFilteredCaloJetCollection);
00145
00146 return;
00147
00148 }
00149
00150
00151 void
00152 JetCollectionForEleHT::beginJob()
00153 {
00154 }
00155
00156
00157 void
00158 JetCollectionForEleHT::endJob() {
00159 }
00160
00161
00162 DEFINE_FWK_MODULE(JetCollectionForEleHT);
00163