Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-08-21 02:34:18

0001 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0002 #include "DataFormats/L1Trigger/interface/P2GTCandidate.h"
0003 #include "DataFormats/L1Trigger/interface/P2GTAlgoBlock.h"
0004 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 
0009 class HLTP2GTTauFilter : public HLTFilter {
0010 public:
0011   explicit HLTP2GTTauFilter(const edm::ParameterSet&);
0012   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0013   bool hltFilter(edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs&) const override;
0014 
0015 private:
0016   edm::InputTag m_l1GTAlgoBlockTag;
0017   edm::EDGetTokenT<l1t::P2GTAlgoBlockMap> m_algoBlockToken;
0018   std::vector<std::string> m_l1GTAlgoNames;
0019   double m_minPt;
0020   unsigned int m_minN;
0021   double m_maxAbsEta;
0022 };
0023 
0024 HLTP2GTTauFilter::HLTP2GTTauFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0025   m_l1GTAlgoBlockTag = iConfig.getParameter<edm::InputTag>("l1GTAlgoBlockTag");
0026   m_algoBlockToken = consumes<l1t::P2GTAlgoBlockMap>(m_l1GTAlgoBlockTag);
0027   m_l1GTAlgoNames = iConfig.getParameter<std::vector<std::string>>("l1GTAlgoNames");
0028   m_minPt = iConfig.getParameter<double>("minPt");
0029   m_minN = iConfig.getParameter<unsigned int>("minN");
0030   m_maxAbsEta = iConfig.getParameter<double>("maxAbsEta");
0031 }
0032 
0033 void HLTP2GTTauFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0034   edm::ParameterSetDescription desc;
0035   makeHLTFilterDescription(desc);
0036   desc.add<edm::InputTag>("l1GTAlgoBlockTag", edm::InputTag(""));
0037   desc.add<std::vector<std::string>>("l1GTAlgoNames", {});
0038   desc.add<double>("minPt", 24);
0039   desc.add<unsigned int>("minN", 1);
0040   desc.add<double>("maxAbsEta", 1e99);
0041   descriptions.add("HLTP2GTTauFilter", desc);
0042 }
0043 
0044 bool HLTP2GTTauFilter::hltFilter(edm::Event& iEvent,
0045                                  const edm::EventSetup& iSetup,
0046                                  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0047   std::vector<l1t::P2GTCandidateRef> vl1cands;
0048   if (saveTags()) {
0049     filterproduct.addCollectionTag(m_l1GTAlgoBlockTag);  // algorithm blocks
0050   }
0051   bool check_l1match = true;
0052   if (m_l1GTAlgoBlockTag.isUninitialized() || m_l1GTAlgoNames.empty())
0053     check_l1match = false;
0054   if (check_l1match) {
0055     const l1t::P2GTAlgoBlockMap& algos = iEvent.get(m_algoBlockToken);
0056     for (const auto& algoName : m_l1GTAlgoNames) {
0057       auto it = algos.find(algoName);
0058       if (it != algos.end() && it->second.decisionBeforeBxMaskAndPrescale()) {
0059         const l1t::P2GTCandidateVectorRef& objects = it->second.trigObjects();
0060         for (const l1t::P2GTCandidateRef& obj : objects) {
0061           if (obj->objectType() == l1t::P2GTCandidate::CL2Taus && obj->pt() >= m_minPt &&
0062               std::abs(obj->eta()) <= m_maxAbsEta) {
0063             vl1cands.push_back(obj);
0064           }
0065         }
0066       }
0067     }
0068   }
0069 
0070   if (saveTags()) {
0071     edm::InputTag tagOld;
0072     for (size_t i = 0; i < vl1cands.size(); ++i) {
0073       const auto& cand = vl1cands[i];
0074       const edm::ProductID pid(cand.id());
0075       const auto& prov = iEvent.getStableProvenance(pid);
0076       edm::InputTag tagNew(prov.moduleLabel(), prov.productInstanceName(), prov.processName());
0077       if (tagNew.encode() != tagOld.encode()) {
0078         filterproduct.addCollectionTag(tagNew);
0079         tagOld = tagNew;
0080       }
0081     }
0082   }
0083   for (const auto& vl1cand : vl1cands) {
0084     filterproduct.addObject(trigger::TriggerL1Tau, vl1cand);
0085   }
0086 
0087   return vl1cands.size() >= m_minN;
0088 }
0089 DEFINE_FWK_MODULE(HLTP2GTTauFilter);