Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:44

0001 /*
0002  * \class TauTagFilter
0003  *
0004  * Filter tau candidates based on tagger scores.
0005  *
0006  * \author Konstantin Androsov, EPFL and ETHZ
0007  */
0008 
0009 #include "DataFormats/BTauReco/interface/JetTag.h"
0010 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Utilities/interface/Exception.h"
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "FWCore/Utilities/interface/InputTag.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0017 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0018 #include "RecoTauTag/RecoTau/interface/TauWPThreshold.h"
0019 
0020 class TauTagFilter : public HLTFilter {
0021 public:
0022   using TauCollection = reco::PFJetCollection;
0023   using TauTagCollection = reco::JetTagCollection;
0024   using TauRef = reco::PFJetRef;
0025   using Selector = tau::TauWPThreshold;
0026   using LorentzVectorM = math::PtEtaPhiMLorentzVector;
0027 
0028   explicit TauTagFilter(const edm::ParameterSet& cfg)
0029       : HLTFilter(cfg),
0030         nExpected_(cfg.getParameter<int>("nExpected")),
0031         tausSrc_(cfg.getParameter<edm::InputTag>("taus")),
0032         tausToken_(consumes<TauCollection>(tausSrc_)),
0033         tauTagsToken_(consumes<TauTagCollection>(cfg.getParameter<edm::InputTag>("tauTags"))),
0034         tauPtCorrToken_(mayConsume<TauTagCollection>(cfg.getParameter<edm::InputTag>("tauPtCorr"))),
0035         seedsSrc_(mayConsume<trigger::TriggerFilterObjectWithRefs>(cfg.getParameter<edm::InputTag>("seeds"))),
0036         seedTypes_(cfg.getParameter<std::vector<int>>("seedTypes")),
0037         selector_(cfg.getParameter<std::string>("selection")),
0038         minPt_(cfg.getParameter<double>("minPt")),
0039         maxEta_(cfg.getParameter<double>("maxEta")),
0040         usePtCorr_(cfg.getParameter<bool>("usePtCorr")),
0041         matchWithSeeds_(cfg.getParameter<bool>("matchWithSeeds") && cfg.getParameter<double>("matchingdR") >= 0),
0042         matchingdR2_(std::pow(cfg.getParameter<double>("matchingdR"), 2)) {
0043     if (cfg.getParameter<bool>("matchWithSeeds") && cfg.getParameter<double>("matchingdR") < 0)
0044       edm::LogWarning("TauTagFilter") << "Matching with seeds is disabled because matchingdR < 0";
0045 
0046     extractMomenta();  // checking that all seed types are supported
0047   }
0048 
0049   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050     edm::ParameterSetDescription desc;
0051     makeHLTFilterDescription(desc);
0052     desc.add<int>("nExpected", 2)->setComment("number of expected taus per event");
0053     desc.add<edm::InputTag>("taus", edm::InputTag(""))->setComment("input collection of taus");
0054     desc.add<edm::InputTag>("tauTags", edm::InputTag(""))->setComment("input collection of tau tagger scores");
0055     desc.add<edm::InputTag>("tauPtCorr", edm::InputTag(""))
0056         ->setComment("input collection of multiplicative tau pt corrections");
0057     desc.add<edm::InputTag>("seeds", edm::InputTag(""))->setComment("input collection of seeds");
0058     desc.add<std::vector<int>>("seedTypes",
0059                                {trigger::TriggerL1Tau, trigger::TriggerL1Jet, trigger::TriggerTau, trigger::TriggerJet})
0060         ->setComment("list of seed object types");
0061     desc.add<std::string>("selection", "0")->setComment("selection formula");
0062     desc.add<double>("minPt", 20)->setComment("minimal tau pt");
0063     desc.add<double>("maxEta", 2.5)->setComment("maximal tau abs(eta)");
0064     desc.add<bool>("usePtCorr", false)->setComment("use multiplicative tau pt corrections");
0065     desc.add<bool>("matchWithSeeds", false)->setComment("apply match with seeds");
0066     desc.add<double>("matchingdR", 0.5)->setComment("deltaR for matching with seeds");
0067     descriptions.addWithDefaultLabel(desc);
0068   }
0069 
0070   bool hltFilter(edm::Event& event,
0071                  const edm::EventSetup& eventsetup,
0072                  trigger::TriggerFilterObjectWithRefs& filterproduct) const override {
0073     if (saveTags())
0074       filterproduct.addCollectionTag(tausSrc_);
0075 
0076     int nTauPassed = 0;
0077 
0078     const auto tausHandle = event.getHandle(tausToken_);
0079     const auto& taus = *tausHandle;
0080 
0081     const std::vector<LorentzVectorM> seed_p4s = extractMomenta(&event);
0082     auto hasMatch = [&](const LorentzVectorM& p4) {
0083       for (const auto& seed_p4 : seed_p4s) {
0084         if (reco::deltaR2(p4, seed_p4) < matchingdR2_)
0085           return true;
0086       }
0087       return false;
0088     };
0089 
0090     const auto& tauTags = event.get(tauTagsToken_);
0091     const TauTagCollection* tauPtCorr = nullptr;
0092     if (usePtCorr_)
0093       tauPtCorr = &event.get(tauPtCorrToken_);
0094 
0095     if (taus.size() != tauTags.size())
0096       throw cms::Exception("Inconsistent Data", "TauTagFilter::hltFilter") << "taus.size() != tauTags.size()";
0097     if (usePtCorr_ && taus.size() != tauPtCorr->size())
0098       throw cms::Exception("Inconsistent Data", "TauTagFilter::hltFilter") << "taus.size() != tauPtCorr.size()";
0099 
0100     for (size_t tau_idx = 0; tau_idx < taus.size(); ++tau_idx) {
0101       const auto& tau = taus[tau_idx];
0102       double pt = tau.pt();
0103       if (usePtCorr_)
0104         pt *= (*tauPtCorr)[tau_idx].second;
0105       const double eta = std::abs(tau.eta());
0106       if (pt > minPt_ && eta < maxEta_ && (!matchWithSeeds_ || hasMatch(tau.polarP4()))) {
0107         const double tag = tauTags[tau_idx].second;
0108         const double tag_thr = selector_(tau);
0109         if (tag > tag_thr) {
0110           filterproduct.addObject(trigger::TriggerTau, TauRef(tausHandle, tau_idx));
0111           nTauPassed++;
0112         }
0113       }
0114     }
0115 
0116     return nTauPassed >= nExpected_;
0117   }
0118 
0119 private:
0120   std::vector<LorentzVectorM> extractMomenta(const edm::Event* event = nullptr) const {
0121     std::vector<LorentzVectorM> seed_p4s;
0122     if (matchWithSeeds_) {
0123       const trigger::TriggerFilterObjectWithRefs* seeds = nullptr;
0124       if (event)
0125         seeds = &event->get(seedsSrc_);
0126       for (const int seedType : seedTypes_) {
0127         if (seedType == trigger::TriggerL1Tau) {
0128           extractMomenta<l1t::TauVectorRef>(seeds, seedType, seed_p4s);
0129         } else if (seedType == trigger::TriggerL1Jet) {
0130           extractMomenta<l1t::JetVectorRef>(seeds, seedType, seed_p4s);
0131         } else if (seedType == trigger::TriggerTau) {
0132           extractMomenta<std::vector<reco::PFTauRef>>(seeds, seedType, seed_p4s);
0133         } else if (seedType == trigger::TriggerJet) {
0134           extractMomenta<std::vector<reco::PFJetRef>>(seeds, seedType, seed_p4s);
0135         } else
0136           throw cms::Exception("Invalid seed type", "TauTagFilter::extractMomenta")
0137               << "Unsupported seed type: " << seedType;
0138       }
0139     }
0140     return seed_p4s;
0141   }
0142 
0143   template <typename Collection>
0144   static void extractMomenta(const trigger::TriggerRefsCollections* triggerObjects,
0145                              int objType,
0146                              std::vector<LorentzVectorM>& p4s) {
0147     if (triggerObjects) {
0148       Collection objects;
0149       triggerObjects->getObjects(objType, objects);
0150       for (const auto& obj : objects)
0151         p4s.push_back(obj->polarP4());
0152     }
0153   }
0154 
0155 private:
0156   const int nExpected_;
0157   const edm::InputTag tausSrc_;
0158   const edm::EDGetTokenT<TauCollection> tausToken_;
0159   const edm::EDGetTokenT<TauTagCollection> tauTagsToken_, tauPtCorrToken_;
0160   const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> seedsSrc_;
0161   const std::vector<int> seedTypes_;
0162   const Selector selector_;
0163   const double minPt_, maxEta_;
0164   const bool usePtCorr_;
0165   const bool matchWithSeeds_;
0166   const double matchingdR2_;
0167 };
0168 
0169 //define this as a plug-in
0170 #include "FWCore/Framework/interface/MakerMacros.h"
0171 DEFINE_FWK_MODULE(TauTagFilter);