File indexing completed on 2024-04-06 12:18:34
0001 #include <cmath>
0002 #include <vector>
0003
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0009 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0010 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0011 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0012 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0013 #include "HLTrigger/JetMET/interface/HLTJetCollForElePlusJets.h"
0014
0015 #include <TVector3.h>
0016
0017 template <typename T>
0018 HLTJetCollForElePlusJets<T>::HLTJetCollForElePlusJets(const edm::ParameterSet& iConfig)
0019 : hltElectronTag(iConfig.getParameter<edm::InputTag>("HltElectronTag")),
0020 sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
0021 minJetPt_(iConfig.getParameter<double>("MinJetPt")),
0022 maxAbsJetEta_(iConfig.getParameter<double>("MaxAbsJetEta")),
0023 minNJets_(iConfig.getParameter<unsigned int>("MinNJets")),
0024
0025 minDeltaR2_(iConfig.getParameter<double>("minDeltaR") * std::abs(iConfig.getParameter<double>("minDeltaR"))),
0026
0027 minSoftJetPt_(iConfig.getParameter<double>("MinSoftJetPt")),
0028 minDeltaEta_(iConfig.getParameter<double>("MinDeltaEta")) {
0029 typedef std::vector<T> TCollection;
0030 m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
0031 m_theJetToken = consumes<TCollection>(sourceJetTag);
0032 produces<TCollection>();
0033 }
0034
0035 template <typename T>
0036 void HLTJetCollForElePlusJets<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0037 edm::ParameterSetDescription desc;
0038 desc.add<edm::InputTag>("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
0039 desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
0040 desc.add<double>("MinJetPt", 30.);
0041 desc.add<double>("MaxAbsJetEta", 2.6);
0042 desc.add<unsigned int>("MinNJets", 1);
0043 desc.add<double>("minDeltaR", 0.5);
0044
0045 desc.add<double>("MinSoftJetPt", 25.);
0046 desc.add<double>("MinDeltaEta", -1.);
0047 descriptions.add(defaultModuleLabel<HLTJetCollForElePlusJets<T>>(), desc);
0048 }
0049
0050
0051
0052
0053
0054
0055 template <typename T>
0056 void HLTJetCollForElePlusJets<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0057 using namespace edm;
0058 using namespace std;
0059
0060 typedef vector<T> TCollection;
0061 typedef Ref<TCollection> TRef;
0062 typedef edm::RefVector<TCollection> TRefVector;
0063 typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
0064 TCollectionVector;
0065
0066 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0067 iEvent.getByToken(m_theElectronToken, PrevFilterOutput);
0068
0069
0070
0071 std::vector<edm::Ref<reco::RecoEcalCandidateCollection>> clusCands;
0072 PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
0073
0074 std::vector<edm::Ref<reco::ElectronCollection>> eleCands;
0075 PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
0076
0077 trigger::VRphoton photonCands;
0078 PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
0079
0080
0081 std::vector<TVector3> ElePs;
0082
0083 if (!clusCands.empty()) {
0084 for (auto& clusCand : clusCands) {
0085 TVector3 positionVector(clusCand->superCluster()->position().x(),
0086 clusCand->superCluster()->position().y(),
0087 clusCand->superCluster()->position().z());
0088 ElePs.push_back(positionVector);
0089 }
0090 } else if (!eleCands.empty()) {
0091 for (auto& eleCand : eleCands) {
0092 TVector3 positionVector(eleCand->superCluster()->position().x(),
0093 eleCand->superCluster()->position().y(),
0094 eleCand->superCluster()->position().z());
0095 ElePs.push_back(positionVector);
0096 }
0097 } else if (!photonCands.empty()) {
0098 for (auto& photonCand : photonCands) {
0099 TVector3 positionVector(photonCand->superCluster()->position().x(),
0100 photonCand->superCluster()->position().y(),
0101 photonCand->superCluster()->position().z());
0102 ElePs.push_back(positionVector);
0103 }
0104 }
0105
0106 edm::Handle<TCollection> theJetCollectionHandle;
0107 iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
0108
0109 const TCollection& theJetCollection = *theJetCollectionHandle;
0110
0111 std::unique_ptr<TCollection> theFilteredJetCollection(new TCollection);
0112
0113 std::unique_ptr<TCollectionVector> allSelections(new TCollectionVector());
0114
0115 bool foundSolution(false);
0116
0117 for (auto const& EleP : ElePs) {
0118 bool VBFJetPair = false;
0119 std::vector<int> store_jet;
0120 TRefVector refVector;
0121
0122 for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0123 TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
0124
0125 double const Deta = EleP.Eta() - JetP.Eta();
0126 double const Dphi = EleP.DeltaPhi(JetP);
0127 double const DR2 = Deta * Deta + Dphi * Dphi;
0128
0129 if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR2 > minDeltaR2_) {
0130 store_jet.push_back(j);
0131
0132 if (minDeltaEta_ > 0) {
0133 for (unsigned int k = j + 1; k < theJetCollection.size(); k++) {
0134 TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(), theJetCollection[k].pz());
0135
0136 if (std::abs(SoftJetP.Eta() - JetP.Eta()) <= minDeltaEta_) {
0137 continue;
0138 }
0139
0140 double const softDeta = EleP.Eta() - SoftJetP.Eta();
0141 double const softDphi = EleP.DeltaPhi(SoftJetP);
0142 double const softDR2 = softDeta * softDeta + softDphi * softDphi;
0143
0144 if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR2 > minDeltaR2_) {
0145 store_jet.push_back(k);
0146 VBFJetPair = true;
0147 }
0148 }
0149 }
0150 }
0151 }
0152
0153
0154 std::sort(store_jet.begin(), store_jet.end());
0155 store_jet.erase(unique(store_jet.begin(), store_jet.end()), store_jet.end());
0156
0157
0158 for (int& ijet : store_jet) {
0159
0160 refVector.push_back(TRef(theJetCollectionHandle, ijet));
0161
0162 if (!foundSolution)
0163 theFilteredJetCollection->push_back(theJetCollection[ijet]);
0164 }
0165
0166 allSelections->push_back(refVector);
0167
0168 if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
0169 foundSolution = true;
0170 else if (VBFJetPair && minDeltaEta_ > 0)
0171 foundSolution = true;
0172 else if (!foundSolution)
0173 theFilteredJetCollection->clear();
0174 }
0175
0176 iEvent.put(std::move(theFilteredJetCollection));
0177
0178 return;
0179 }