File indexing completed on 2024-04-06 12:18:33
0001 #include <string>
0002 #include <vector>
0003
0004 #include "TVector3.h"
0005
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0008 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0009 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0012 #include "HLTrigger/JetMET/interface/HLTJetCollectionsForElePlusJets.h"
0013 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0014
0015 template <typename T>
0016 HLTJetCollectionsForElePlusJets<T>::HLTJetCollectionsForElePlusJets(const edm::ParameterSet& iConfig)
0017 : hltElectronTag(iConfig.getParameter<edm::InputTag>("HltElectronTag")),
0018 sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
0019
0020
0021
0022 minDeltaR_(iConfig.getParameter<double>("minDeltaR"))
0023
0024
0025
0026 {
0027 typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
0028 TCollectionVector;
0029 m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
0030 m_theJetToken = consumes<std::vector<T>>(sourceJetTag);
0031 produces<TCollectionVector>();
0032 }
0033
0034 template <typename T>
0035 HLTJetCollectionsForElePlusJets<T>::~HLTJetCollectionsForElePlusJets() {
0036
0037
0038 }
0039
0040 template <typename T>
0041 void HLTJetCollectionsForElePlusJets<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0042 edm::ParameterSetDescription desc;
0043 desc.add<edm::InputTag>("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
0044 desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
0045
0046
0047
0048 desc.add<double>("minDeltaR", 0.5);
0049
0050
0051
0052 descriptions.add(defaultModuleLabel<HLTJetCollectionsForElePlusJets<T>>(), desc);
0053 }
0054
0055
0056
0057
0058
0059
0060 template <typename T>
0061 void HLTJetCollectionsForElePlusJets<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0062 using namespace edm;
0063 using namespace std;
0064
0065 typedef vector<T> TCollection;
0066 typedef Ref<TCollection> TRef;
0067
0068 typedef edm::RefVector<TCollection> TRefVector;
0069
0070 typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
0071 TCollectionVector;
0072
0073 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0074 iEvent.getByToken(m_theElectronToken, PrevFilterOutput);
0075
0076
0077
0078 std::vector<edm::Ref<reco::RecoEcalCandidateCollection>> clusCands;
0079 PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
0080
0081 std::vector<edm::Ref<reco::ElectronCollection>> eleCands;
0082 PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
0083
0084 trigger::VRphoton photonCands;
0085 PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
0086
0087
0088 std::vector<TVector3> ElePs;
0089
0090 if (!clusCands.empty()) {
0091 for (auto& clusCand : clusCands) {
0092 TVector3 positionVector(clusCand->superCluster()->position().x(),
0093 clusCand->superCluster()->position().y(),
0094 clusCand->superCluster()->position().z());
0095 ElePs.push_back(positionVector);
0096 }
0097 } else if (!eleCands.empty()) {
0098 for (auto& eleCand : eleCands) {
0099 TVector3 positionVector(eleCand->superCluster()->position().x(),
0100 eleCand->superCluster()->position().y(),
0101 eleCand->superCluster()->position().z());
0102 ElePs.push_back(positionVector);
0103 }
0104 } else if (!photonCands.empty()) {
0105 for (auto& photonCand : photonCands) {
0106 TVector3 positionVector(photonCand->superCluster()->position().x(),
0107 photonCand->superCluster()->position().y(),
0108 photonCand->superCluster()->position().z());
0109 ElePs.push_back(positionVector);
0110 }
0111 }
0112
0113 edm::Handle<TCollection> theJetCollectionHandle;
0114 iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
0115
0116 const TCollection& theJetCollection = *theJetCollectionHandle;
0117
0118
0119
0120 std::unique_ptr<TCollectionVector> allSelections(new TCollectionVector());
0121
0122
0123
0124 for (auto& EleP : ElePs) {
0125
0126
0127 TRefVector refVector;
0128
0129 for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0130 TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
0131 double DR = EleP.DeltaR(JetP);
0132
0133 if (DR > minDeltaR_)
0134 refVector.push_back(TRef(theJetCollectionHandle, j));
0135 }
0136 allSelections->push_back(refVector);
0137 }
0138
0139
0140 iEvent.put(std::move(allSelections));
0141
0142 return;
0143 }