File indexing completed on 2024-04-06 12:18:33
0001 #include <cmath>
0002
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0005 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0006 #include "DataFormats/Math/interface/deltaR.h"
0007 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0008 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0009 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
0010 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0015 #include "HLTrigger/JetMET/interface/HLTJetCollectionsForLeptonPlusJets.h"
0016
0017 template <typename jetType>
0018 HLTJetCollectionsForLeptonPlusJets<jetType>::HLTJetCollectionsForLeptonPlusJets(const edm::ParameterSet& iConfig)
0019 : hltLeptonTag(iConfig.getParameter<edm::InputTag>("HltLeptonTag")),
0020 sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
0021
0022 minDeltaR2_(iConfig.getParameter<double>("minDeltaR") * std::abs(iConfig.getParameter<double>("minDeltaR"))) {
0023 using namespace edm;
0024 using namespace std;
0025 typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
0026 JetCollectionVector;
0027 m_theLeptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltLeptonTag);
0028 m_theJetToken = consumes<std::vector<jetType>>(sourceJetTag);
0029 produces<JetCollectionVector>();
0030 }
0031
0032 template <typename jetType>
0033 void HLTJetCollectionsForLeptonPlusJets<jetType>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0034 edm::ParameterSetDescription desc;
0035 desc.add<edm::InputTag>("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
0036 desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("caloJetCollection"));
0037 desc.add<double>("minDeltaR", 0.5);
0038 descriptions.add(defaultModuleLabel<HLTJetCollectionsForLeptonPlusJets<jetType>>(), desc);
0039 }
0040
0041
0042
0043
0044
0045
0046
0047 template <typename jetType>
0048 void HLTJetCollectionsForLeptonPlusJets<jetType>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0049 using namespace edm;
0050 using namespace std;
0051
0052 typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
0053 JetCollectionVector;
0054 typedef vector<jetType> JetCollection;
0055 typedef edm::RefVector<JetCollection> JetRefVector;
0056 typedef edm::Ref<JetCollection> JetRef;
0057
0058 Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0059 iEvent.getByToken(m_theLeptonToken, PrevFilterOutput);
0060
0061
0062
0063 vector<Ref<reco::RecoEcalCandidateCollection>> clusCands;
0064 PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
0065
0066 vector<Ref<reco::ElectronCollection>> eleCands;
0067 PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
0068
0069 trigger::VRphoton photonCands;
0070 PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
0071
0072 vector<reco::RecoChargedCandidateRef> muonCands;
0073 PrevFilterOutput->getObjects(trigger::TriggerMuon, muonCands);
0074
0075 Handle<JetCollection> theJetCollectionHandle;
0076 iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
0077
0078 const JetCollection& theJetCollection = *theJetCollectionHandle;
0079
0080 unique_ptr<JetCollectionVector> allSelections(new JetCollectionVector());
0081
0082 if (!clusCands.empty()) {
0083 for (auto const& clusCand : clusCands) {
0084 JetRefVector refVector;
0085 for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0086 if (reco::deltaR2(clusCand->superCluster()->position(), theJetCollection[j]) > minDeltaR2_)
0087 refVector.push_back(JetRef(theJetCollectionHandle, j));
0088 }
0089 allSelections->push_back(refVector);
0090 }
0091 }
0092
0093 if (!eleCands.empty()) {
0094 for (auto const& eleCand : eleCands) {
0095 JetRefVector refVector;
0096 for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0097 if (reco::deltaR2(eleCand->superCluster()->position(), theJetCollection[j]) > minDeltaR2_)
0098 refVector.push_back(JetRef(theJetCollectionHandle, j));
0099 }
0100 allSelections->push_back(refVector);
0101 }
0102 }
0103
0104 if (!photonCands.empty()) {
0105 for (auto const& photonCand : photonCands) {
0106 JetRefVector refVector;
0107 for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0108 if (reco::deltaR2(photonCand->superCluster()->position(), theJetCollection[j]) > minDeltaR2_)
0109 refVector.push_back(JetRef(theJetCollectionHandle, j));
0110 }
0111 allSelections->push_back(refVector);
0112 }
0113 }
0114
0115 if (!muonCands.empty()) {
0116 for (auto const& muonCand : muonCands) {
0117 JetRefVector refVector;
0118 for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0119 if (reco::deltaR2(muonCand->p4(), theJetCollection[j]) > minDeltaR2_)
0120 refVector.push_back(JetRef(theJetCollectionHandle, j));
0121 }
0122 allSelections->push_back(refVector);
0123 }
0124 }
0125
0126 iEvent.put(std::move(allSelections));
0127
0128 return;
0129 }