File indexing completed on 2024-04-06 12:15:51
0001 #include "HLTSumJetTag.h"
0002
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 #include "FWCore/Utilities/interface/Exception.h"
0005 #include "DataFormats/JetReco/interface/JetCollection.h"
0006 #include "DataFormats/Math/interface/deltaR.h"
0007 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0008
0009 #include <numeric>
0010
0011 template <typename T>
0012 HLTSumJetTag<T>::HLTSumJetTag(const edm::ParameterSet& config)
0013 : HLTFilter(config),
0014 m_Jets(config.getParameter<edm::InputTag>("Jets")),
0015 m_JetTags(config.getParameter<edm::InputTag>("JetTags")),
0016 m_JetsToken(consumes(m_Jets)),
0017 m_JetTagsToken(consumes(m_JetTags)),
0018 m_MinTag(config.getParameter<double>("MinTag")),
0019 m_MaxTag(config.getParameter<double>("MaxTag")),
0020 m_MinJetToSum(config.getParameter<int>("MinJetToSum")),
0021 m_MaxJetToSum(config.getParameter<int>("MaxJetToSum")),
0022 m_UseMeanValue(config.getParameter<bool>("UseMeanValue")),
0023 m_MatchByDeltaR(config.getParameter<bool>("MatchByDeltaR")),
0024 m_MaxDeltaR(config.getParameter<double>("MaxDeltaR")),
0025 m_TriggerType(config.getParameter<int>("TriggerType")) {
0026
0027 if (m_MatchByDeltaR and m_MaxDeltaR <= 0) {
0028 throw cms::Exception("HLTSumJetTag") << "invalid value for parameter \"MaxDeltaR\" (must be > 0): " << m_MaxDeltaR;
0029 }
0030
0031 if (m_MaxJetToSum >= 0 and m_MinJetToSum > m_MaxJetToSum) {
0032 throw cms::Exception("HLTSumJetTag")
0033 << "invalid value for min/max number of jets to sum: parameter \"MinJetToSum\" (" << m_MinJetToSum
0034 << ") must be <= than \"MaxJetToSum\" (" << m_MaxJetToSum << ") if \"MaxJetToSum\" >= 0";
0035 }
0036
0037 edm::LogInfo("") << " (HLTSumJetTag) trigger cuts: \n"
0038 << " \ttype of jets used: " << m_Jets.encode() << " \n"
0039 << " \ttype of tagged jets used: " << m_JetTags.encode() << " \n"
0040 << " \tmin/max tag value: [" << m_MinTag << ", " << m_MaxTag << "] \n"
0041 << " \tmin/max number of jets to sum: [" << m_MinJetToSum << ", " << m_MaxJetToSum << "] \n"
0042 << " \tuse mean value of jet tags: " << m_UseMeanValue << " \n"
0043 << " \tassign jet-tag values by Delta-R matching: " << m_MatchByDeltaR << "\n"
0044 << " \tmax Delta-R for jet-tag assignment by Delta-R matching: " << m_MaxDeltaR << "\n"
0045 << " \tTriggerType: " << m_TriggerType << " \n";
0046 }
0047
0048 template <typename T>
0049 void HLTSumJetTag<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050 edm::ParameterSetDescription desc;
0051 makeHLTFilterDescription(desc);
0052 desc.add<edm::InputTag>("Jets", edm::InputTag("hltJetCollection"));
0053 desc.add<edm::InputTag>("JetTags", edm::InputTag("hltJetTagCollection"));
0054 desc.add<double>("MinTag", 0.);
0055 desc.add<double>("MaxTag", 999999.0);
0056 desc.add<int>("MinJetToSum", 1);
0057 desc.add<int>("MaxJetToSum", -1);
0058 desc.add<bool>("UseMeanValue", true);
0059 desc.add<bool>("MatchByDeltaR", false);
0060 desc.add<double>("MaxDeltaR", 0.4);
0061 desc.add<int>("TriggerType", 86);
0062 descriptions.add(defaultModuleLabel<HLTSumJetTag<T>>(), desc);
0063 }
0064
0065 template <typename T>
0066 bool HLTSumJetTag<T>::hltFilter(edm::Event& event,
0067 const edm::EventSetup& setup,
0068 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0069 if (saveTags())
0070 filterproduct.addCollectionTag(m_Jets);
0071
0072 typedef edm::Ref<std::vector<T>> TRef;
0073
0074 auto const h_Jets = event.getHandle(m_JetsToken);
0075 auto const h_JetTags = event.getHandle(m_JetTagsToken);
0076
0077 std::vector<TRef> jetRefCollection;
0078 jetRefCollection.reserve(h_Jets->size());
0079
0080 std::vector<float> jetTagValues;
0081 jetTagValues.reserve(h_Jets->size());
0082
0083 if (m_MaxJetToSum == 0) {
0084
0085 LogDebug("HLTSumJetTag") << "Parameter \"MaxJetToSum\" set to zero --> Return False";
0086 return false;
0087 }
0088
0089
0090 auto const maxDeltaR2 = m_MaxDeltaR * m_MaxDeltaR;
0091 for (size_t iJet = 0; iJet < h_Jets->size(); ++iJet) {
0092 TRef const jetRef(h_Jets, iJet);
0093
0094 float jetTag = -1;
0095 if (m_MatchByDeltaR) {
0096
0097 if (not findTagValueByMinDeltaR2(jetTag, (*h_Jets)[iJet], *h_JetTags, maxDeltaR2)) {
0098 continue;
0099 }
0100 } else {
0101
0102 jetTag = (*h_JetTags)[reco::JetBaseRef(jetRef)];
0103 }
0104
0105 jetRefCollection.emplace_back(jetRef);
0106 jetTagValues.emplace_back(jetTag);
0107
0108 LogDebug("HLTSumJetTag") << "Input Jets -- Jet[" << iJet << "] (id = " << jetRefCollection.back().id() << ")"
0109 << ": tag=" << jetTag << ", pt=" << jetRefCollection.back()->pt()
0110 << ", eta=" << jetRefCollection.back()->eta()
0111 << ", phi=" << jetRefCollection.back()->phi();
0112 }
0113
0114
0115 std::vector<size_t> jetTagSortedIndices(jetTagValues.size());
0116 std::iota(jetTagSortedIndices.begin(), jetTagSortedIndices.end(), 0);
0117 std::sort(jetTagSortedIndices.begin(), jetTagSortedIndices.end(), [&](const size_t& i1, const size_t& i2) {
0118 return jetTagValues[i1] > jetTagValues[i2];
0119 });
0120
0121
0122 if (m_MaxJetToSum > 0 and int(jetTagSortedIndices.size()) > m_MaxJetToSum)
0123 jetTagSortedIndices.resize(m_MaxJetToSum);
0124
0125
0126 float sumJetTag = 0;
0127 for (auto const idx : jetTagSortedIndices) {
0128 sumJetTag += jetTagValues[idx];
0129 LogDebug("HLTSumJetTag") << "Selected Jets -- Jet[" << idx << "] (id = " << jetRefCollection[idx].id() << ")"
0130 << ": tag=" << jetTagValues[idx] << ", pt=" << jetRefCollection[idx]->pt()
0131 << ", eta=" << jetRefCollection[idx]->eta() << ", phi=" << jetRefCollection[idx]->phi();
0132 }
0133
0134 if (m_UseMeanValue and not jetTagSortedIndices.empty()) {
0135 sumJetTag /= jetTagSortedIndices.size();
0136 }
0137
0138
0139 auto const accept =
0140 (int(jetTagSortedIndices.size()) >= m_MinJetToSum and m_MinTag <= sumJetTag and sumJetTag <= m_MaxTag);
0141
0142 LogDebug("HLTSumJetTag") << "Filter Result = " << accept << " (SumJetTag = " << sumJetTag
0143 << ", Num. of Jets = " << jetTagSortedIndices.size()
0144 << ") [UseMeanValue = " << m_UseMeanValue << "]";
0145
0146
0147 if (accept) {
0148 for (auto const idx : jetTagSortedIndices) {
0149 filterproduct.addObject(m_TriggerType, jetRefCollection[idx]);
0150 }
0151 }
0152
0153 return accept;
0154 }
0155
0156 template <typename T>
0157 bool HLTSumJetTag<T>::findTagValueByMinDeltaR2(float& jetTagValue,
0158 const T& jet,
0159 const reco::JetTagCollection& jetTags,
0160 float maxDeltaR2) const {
0161 bool ret = false;
0162 for (const auto& jetTag : jetTags) {
0163 auto const tmpDR2 = reco::deltaR2(jet, *(jetTag.first));
0164 if (tmpDR2 < maxDeltaR2) {
0165 maxDeltaR2 = tmpDR2;
0166 jetTagValue = jetTag.second;
0167 ret = true;
0168 }
0169 }
0170
0171 return ret;
0172 }