Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-31 02:19:39

0001 #ifndef HLTrigger_btau_HLTJetTag_h
0002 #define HLTrigger_btau_HLTJetTag_h
0003 
0004 /** \class HLTJetTag
0005  *
0006  *  This class is an HLTFilter (a specialized EDFilter) implementing
0007  *  tagged multi-jet trigger for b and tau.
0008  *  It should be run after the normal multi-jet trigger.
0009  *
0010  *
0011  *  \author Arnaud Gay, Ian Tomalin
0012  *  \maintainer Andrea Bocci
0013  *
0014  */
0015 
0016 #include "DataFormats/BTauReco/interface/JetTag.h"
0017 #include "DataFormats/Common/interface/Ref.h"
0018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0019 #include "DataFormats/Math/interface/deltaR.h"
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0023 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0024 #include "FWCore/Utilities/interface/Exception.h"
0025 #include "FWCore/Utilities/interface/InputTag.h"
0026 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0027 
0028 template <typename T>
0029 class HLTJetTag : public HLTFilter {
0030 public:
0031   explicit HLTJetTag(edm::ParameterSet const& config);
0032   ~HLTJetTag() override = default;
0033 
0034   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035 
0036   bool hltFilter(edm::Event& event,
0037                  edm::EventSetup const& setup,
0038                  trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
0039 
0040 private:
0041   auto logTrace() const {
0042     auto const moduleName = moduleDescription().moduleName();
0043     return LogTrace(moduleName) << "[" << moduleName << "] ";
0044   }
0045 
0046   bool findTagValueByMinDeltaR2(float& jetTagValue,
0047                                 T const& jet,
0048                                 reco::JetTagCollection const& jetTags,
0049                                 float maxDeltaR2) const;
0050 
0051   edm::InputTag const m_Jets;  // module label of input JetCollection
0052   edm::EDGetTokenT<std::vector<T>> const m_JetsToken;
0053   edm::InputTag const m_JetTags;  // module label of input JetTagCollection
0054   edm::EDGetTokenT<reco::JetTagCollection> const m_JetTagsToken;
0055   double const m_MinTag, m_MaxTag;               // tag discriminator cuts applied to each jet
0056   int const m_MinJets;                           // min. number of jets required to be tagged
0057   bool const m_MatchJetsByDeltaR;                // flag to enable jet matching by minimum Delta-R distance
0058   double const m_MaxJetDeltaR, m_MaxJetDeltaR2;  // max delta-R distance used in jet matching
0059   int const m_TriggerType;                       // type of trigger objects in output product
0060 };
0061 
0062 template <typename T>
0063 HLTJetTag<T>::HLTJetTag(edm::ParameterSet const& config)
0064     : HLTFilter(config),
0065       m_Jets(config.getParameter<edm::InputTag>("Jets")),
0066       m_JetsToken(consumes(m_Jets)),
0067       m_JetTags(config.getParameter<edm::InputTag>("JetTags")),
0068       m_JetTagsToken(consumes(m_JetTags)),
0069       m_MinTag(config.getParameter<double>("MinTag")),
0070       m_MaxTag(config.getParameter<double>("MaxTag")),
0071       m_MinJets(config.getParameter<int>("MinJets")),
0072       m_MatchJetsByDeltaR(config.getParameter<bool>("MatchJetsByDeltaR")),
0073       m_MaxJetDeltaR(config.getParameter<double>("MaxJetDeltaR")),
0074       m_MaxJetDeltaR2(m_MaxJetDeltaR * m_MaxJetDeltaR),
0075       m_TriggerType(config.getParameter<int>("TriggerType")) {
0076   if (m_MatchJetsByDeltaR and m_MaxJetDeltaR < 0) {
0077     throw cms::Exception("InvalidConfigurationParameter")
0078         << "invalid value for parameter \"MaxJetDeltaR\" (must be >= 0): " << m_MaxJetDeltaR;
0079   }
0080 
0081   logTrace() << "HLTFilter parameters:\n"
0082              << "\ttype of jets used: " << m_Jets.encode() << "\n"
0083              << "\ttype of tagged jets used: " << m_JetTags.encode() << "\n"
0084              << "\tmin/max tag value: [" << m_MinTag << ", " << m_MaxTag << "]\n"
0085              << "\tmin num. of tagged jets: " << m_MinJets << "\n"
0086              << "\tmatch jets by Delta-R: " << m_MatchJetsByDeltaR << "\n"
0087              << "\tmax Delta-R for jet matching: " << m_MaxJetDeltaR << "\n"
0088              << "\tTriggerType: " << m_TriggerType << "\n";
0089 }
0090 
0091 template <typename T>
0092 void HLTJetTag<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0093   edm::ParameterSetDescription desc;
0094   makeHLTFilterDescription(desc);
0095   desc.add<edm::InputTag>("Jets", edm::InputTag("hltJetCollection"));
0096   desc.add<edm::InputTag>("JetTags", edm::InputTag("hltJetTagCollection"));
0097   desc.add<double>("MinTag", 2.0);
0098   desc.add<double>("MaxTag", 999999.0);
0099   desc.add<int>("MinJets", 1);
0100   desc.add<bool>("MatchJetsByDeltaR", false);
0101   desc.add<double>("MaxJetDeltaR", 0.1);
0102   desc.add<int>("TriggerType", 0);
0103   descriptions.addWithDefaultLabel(desc);
0104 }
0105 
0106 template <typename T>
0107 bool HLTJetTag<T>::findTagValueByMinDeltaR2(float& jetTagValue,
0108                                             T const& jet,
0109                                             reco::JetTagCollection const& jetTags,
0110                                             float maxDeltaR2) const {
0111   bool ret = false;
0112   for (auto const& jetTag : jetTags) {
0113     auto const tmpDR2 = reco::deltaR2(jet, *(jetTag.first));
0114     if (tmpDR2 < maxDeltaR2) {
0115       maxDeltaR2 = tmpDR2;
0116       jetTagValue = jetTag.second;
0117       ret = true;
0118     }
0119   }
0120 
0121   return ret;
0122 }
0123 
0124 template <typename T>
0125 bool HLTJetTag<T>::hltFilter(edm::Event& event,
0126                              edm::EventSetup const& setup,
0127                              trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0128   auto const h_Jets = event.getHandle(m_JetsToken);
0129   if (saveTags())
0130     filterproduct.addCollectionTag(m_Jets);
0131 
0132   auto const h_JetTags = event.getHandle(m_JetTagsToken);
0133 
0134   // check if the product this one depends on is available
0135   auto const& handle = h_JetTags;
0136   auto const& dependent = handle->keyProduct();
0137   if (not dependent.isNull() and not dependent.hasCache()) {
0138     // only an empty AssociationVector can have an invalid dependent collection
0139     edm::StableProvenance const& dependent_provenance = event.getStableProvenance(dependent.id());
0140     if (dependent_provenance.productDescription().dropped())
0141       // FIXME the error message should be made prettier
0142       throw cms::Exception("ProductNotFound") << "Product " << handle.provenance()->branchName() << " requires product "
0143                                               << dependent_provenance.branchName() << ", which has been dropped";
0144   }
0145 
0146   int nJetTags = 0;
0147 
0148   for (size_t iJet = 0; iJet < h_Jets->size(); ++iJet) {
0149     auto const& jet = (*h_Jets)[iJet];
0150     auto const jetRef = edm::Ref<std::vector<T>>(h_Jets, iJet);
0151 
0152     float jetTag = -1.f;
0153     if (m_MatchJetsByDeltaR) {
0154       // if match-by-DR is used, update jetTag by reference and use it only if the matching is valid
0155       if (not findTagValueByMinDeltaR2(jetTag, jet, *h_JetTags, m_MaxJetDeltaR2)) {
0156         continue;
0157       }
0158     } else {
0159       // operator[] checks consistency between h_Jets and h_JetTags
0160       jetTag = (*h_JetTags)[reco::JetBaseRef(jetRef)];
0161     }
0162 
0163     // Check if jet is tagged
0164     bool const jetIsTagged = (m_MinTag <= jetTag) and (jetTag <= m_MaxTag);
0165 
0166     if (jetIsTagged) {
0167       // Store a reference to the jets which passed tagging cuts
0168       filterproduct.addObject(m_TriggerType, jetRef);
0169       ++nJetTags;
0170     }
0171 
0172     logTrace() << "Matched Jets -- Jet[" << iJet << "]"
0173                << ": pt=" << jet.pt() << ", eta=" << jet.eta() << ", phi=" << jet.phi() << ", tag=" << jetTag
0174                << ", isTagged=" << jetIsTagged;
0175   }
0176 
0177   // filter decision
0178   bool const accept = (nJetTags >= m_MinJets);
0179 
0180   logTrace() << " trigger accept = " << accept << ", nJetTags = " << nJetTags;
0181 
0182   return accept;
0183 }
0184 
0185 #endif  // HLTrigger_btau_HLTJetTag_h