Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "RecoTauTag/HLTProducers/interface/PFTauL1TJetsMatching.h"
0002 #include "Math/GenVector/VectorUtil.h"
0003 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0004 #include "FWCore/Utilities/interface/EDMException.h"
0005 #include "DataFormats/TauReco/interface/PFTau.h"
0006 #include "DataFormats/Math/interface/deltaR.h"
0007 
0008 //
0009 // class declaration
0010 //
0011 PFTauL1TJetsMatching::PFTauL1TJetsMatching(const edm::ParameterSet& iConfig)
0012     : tauSrc_(consumes<reco::PFTauCollection>(iConfig.getParameter<edm::InputTag>("TauSrc"))),
0013       L1JetSrc_(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("L1JetSrc"))),
0014       matchingR2_(iConfig.getParameter<double>("MatchingdR") * iConfig.getParameter<double>("MatchingdR")),
0015       minTauPt_(iConfig.getParameter<double>("MinTauPt")),
0016       minL1TPt_(iConfig.getParameter<double>("MinL1TPt")) {
0017   produces<reco::PFTauCollection>();
0018 }
0019 PFTauL1TJetsMatching::~PFTauL1TJetsMatching() {}
0020 
0021 void PFTauL1TJetsMatching::produce(edm::StreamID iSId, edm::Event& iEvent, const edm::EventSetup& iES) const {
0022   std::unique_ptr<reco::PFTauCollection> L1TmatchedPFTau(new reco::PFTauCollection);
0023 
0024   edm::Handle<reco::PFTauCollection> taus;
0025   iEvent.getByToken(tauSrc_, taus);
0026 
0027   edm::Handle<trigger::TriggerFilterObjectWithRefs> L1Jets;
0028   iEvent.getByToken(L1JetSrc_, L1Jets);
0029 
0030   l1t::JetVectorRef jetCandRefVec;
0031   L1Jets->getObjects(trigger::TriggerL1Jet, jetCandRefVec);
0032 
0033   /* Loop over taus that must pass a certain minTauPt_ cut */
0034   /* then loop over L1T jets that must pass minL1TPt_ */
0035   /* and check whether they match, if yes -> include the taus in */
0036   /* the new L1T matched PFTau collection */
0037 
0038   for (unsigned int iTau = 0; iTau < taus->size(); iTau++) {
0039     bool isMatched = false;
0040     if ((*taus)[iTau].pt() > minTauPt_) {
0041       for (unsigned int iJet = 0; iJet < jetCandRefVec.size(); iJet++) {
0042         if (jetCandRefVec[iJet]->pt() > minL1TPt_) {
0043           if (reco::deltaR2((*taus)[iTau].p4(), jetCandRefVec[iJet]->p4()) < matchingR2_) {
0044             isMatched = true;
0045             break;
0046           }
0047         }
0048       }
0049     }
0050     if (isMatched)
0051       L1TmatchedPFTau->push_back((*taus)[iTau]);
0052   }
0053   iEvent.put(std::move(L1TmatchedPFTau));
0054 }
0055 
0056 void PFTauL1TJetsMatching::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0057   edm::ParameterSetDescription desc;
0058   desc.add<edm::InputTag>("L1JetSrc", edm::InputTag("hltL1VBFDiJetOR"))
0059       ->setComment("Input filter objects passing L1 seed");
0060   desc.add<edm::InputTag>("TauSrc", edm::InputTag("hltSelectedPFTausTrackFindingLooseChargedIsolationAgainstMuon"))
0061       ->setComment("Input collection of PFTaus");
0062   desc.add<double>("MatchingdR", 0.5)->setComment("Maximum dR for matching between PFTaus and L1 filter jets");
0063   desc.add<double>("MinTauPt", 20.0)->setComment("PFTaus above this pt will be considered");
0064   desc.add<double>("MinL1TPt", 115.0)->setComment("L1T Objects above this pt will be considered");
0065   descriptions.setComment(
0066       "This module produces a collection of PFTaus matched to the leading jet passing the L1 seed filter.");
0067   descriptions.add("PFTauL1TJetsMatching", desc);
0068 }
0069 
0070 #include "FWCore/Framework/interface/MakerMacros.h"
0071 DEFINE_FWK_MODULE(PFTauL1TJetsMatching);