File indexing completed on 2024-04-06 12:27:44
0001 #include "RecoTauTag/HLTProducers/interface/PFJetsTauOverlapRemoval.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
0010
0011 PFJetsTauOverlapRemoval::PFJetsTauOverlapRemoval(const edm::ParameterSet& iConfig)
0012 : tauSrc_(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("TauSrc"))),
0013 pfJetSrc_(consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("PFJetSrc"))),
0014 matchingR2_(iConfig.getParameter<double>("Min_dR") * iConfig.getParameter<double>("Min_dR")) {
0015 produces<reco::PFJetCollection>();
0016 }
0017 PFJetsTauOverlapRemoval::~PFJetsTauOverlapRemoval() {}
0018
0019 void PFJetsTauOverlapRemoval::produce(edm::StreamID iSId, edm::Event& iEvent, const edm::EventSetup& iES) const {
0020 std::unique_ptr<reco::PFJetCollection> cleanedPFJets(new reco::PFJetCollection);
0021
0022 edm::Handle<trigger::TriggerFilterObjectWithRefs> tauJets;
0023 iEvent.getByToken(tauSrc_, tauJets);
0024
0025 edm::Handle<reco::PFJetCollection> PFJets;
0026 iEvent.getByToken(pfJetSrc_, PFJets);
0027
0028 trigger::VRpftau taus;
0029 tauJets->getObjects(trigger::TriggerTau, taus);
0030
0031 if (PFJets->size() > 1) {
0032 for (unsigned int iJet = 0; iJet < PFJets->size(); iJet++) {
0033 bool isMatched = false;
0034 const reco::PFJet& myPFJet = (*PFJets)[iJet];
0035 for (unsigned int iTau = 0; iTau < taus.size(); iTau++) {
0036 if (reco::deltaR2(taus[iTau]->p4(), myPFJet.p4()) < matchingR2_) {
0037 isMatched = true;
0038 break;
0039 }
0040 }
0041 if (isMatched == false)
0042 cleanedPFJets->push_back(myPFJet);
0043 }
0044 }
0045 iEvent.put(std::move(cleanedPFJets));
0046 }
0047
0048 void PFJetsTauOverlapRemoval::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0049 edm::ParameterSetDescription desc;
0050 desc.add<edm::InputTag>("PFJetSrc", edm::InputTag("hltAK4PFJetsCorrected"))->setComment("Input collection of PFJets");
0051 desc.add<edm::InputTag>("TauSrc", edm::InputTag("hltSinglePFTau20TrackPt1LooseChargedIsolationReg"))
0052 ->setComment("Input collection of PFTaus that have passed ID and isolation requirements");
0053 desc.add<double>("Min_dR", 0.5)->setComment("Minimum dR outside of which PFJets will be saved");
0054 descriptions.setComment(
0055 "This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT "
0056 "filter.");
0057 descriptions.add("PFJetsTauOverlapRemoval", desc);
0058 }
0059
0060 #include "FWCore/Framework/interface/MakerMacros.h"
0061 DEFINE_FWK_MODULE(PFJetsTauOverlapRemoval);