File indexing completed on 2024-04-06 12:27:43
0001 #include "RecoTauTag/HLTProducers/interface/L1HLTTauMatching.h"
0002 #include "Math/GenVector/VectorUtil.h"
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
0005 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
0006 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0007 #include "FWCore/Utilities/interface/InputTag.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 #include "DataFormats/TauReco/interface/PFTau.h"
0010
0011
0012
0013
0014 using namespace reco;
0015 using namespace std;
0016 using namespace edm;
0017 using namespace l1extra;
0018
0019 L1HLTTauMatching::L1HLTTauMatching(const edm::ParameterSet& iConfig)
0020 : jetSrc(consumes<PFTauCollection>(iConfig.getParameter<InputTag>("JetSrc"))),
0021 tauTrigger(consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger"))),
0022 mEt_Min(iConfig.getParameter<double>("EtMin")) {
0023 produces<PFTauCollection>();
0024 }
0025
0026 void L1HLTTauMatching::produce(edm::StreamID iSId, edm::Event& iEvent, const edm::EventSetup& iES) const {
0027 using namespace edm;
0028 using namespace std;
0029 using namespace reco;
0030 using namespace trigger;
0031 using namespace l1extra;
0032
0033 unique_ptr<PFTauCollection> tauL2jets(new PFTauCollection);
0034
0035 constexpr double matchingR2 = 0.5 * 0.5;
0036
0037
0038 edm::Handle<PFTauCollection> tauJets;
0039 iEvent.getByToken(jetSrc, tauJets);
0040
0041 edm::Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus;
0042 iEvent.getByToken(tauTrigger, l1TriggeredTaus);
0043
0044 vector<l1extra::L1JetParticleRef> tauCandRefVec;
0045 vector<l1extra::L1JetParticleRef> jetCandRefVec;
0046 l1TriggeredTaus->getObjects(trigger::TriggerL1TauJet, tauCandRefVec);
0047 l1TriggeredTaus->getObjects(trigger::TriggerL1CenJet, jetCandRefVec);
0048
0049 math::XYZPoint a(0., 0., 0.);
0050
0051 for (unsigned int iL1Tau = 0; iL1Tau < tauCandRefVec.size(); iL1Tau++) {
0052 for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
0053
0054 const PFTau& myJet = (*tauJets)[iJet];
0055 double deltaR2 = ROOT::Math::VectorUtil::DeltaR2(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
0056 if (deltaR2 < matchingR2) {
0057
0058 if (myJet.leadChargedHadrCand().isNonnull()) {
0059 a = myJet.leadChargedHadrCand()->vertex();
0060 }
0061 PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
0062 if (myJet.pt() > mEt_Min) {
0063
0064 tauL2jets->push_back(myPFTau);
0065 }
0066 break;
0067 }
0068 }
0069 }
0070
0071 for (unsigned int iL1Tau = 0; iL1Tau < jetCandRefVec.size(); iL1Tau++) {
0072 for (unsigned int iJet = 0; iJet < tauJets->size(); iJet++) {
0073 const PFTau& myJet = (*tauJets)[iJet];
0074
0075 double deltaR2 = ROOT::Math::VectorUtil::DeltaR2(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
0076 if (deltaR2 < matchingR2) {
0077
0078 if (myJet.leadChargedHadrCand().isNonnull()) {
0079 a = myJet.leadChargedHadrCand()->vertex();
0080 }
0081
0082 PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
0083 if (myJet.pt() > mEt_Min) {
0084
0085 tauL2jets->push_back(myPFTau);
0086 }
0087 break;
0088 }
0089 }
0090 }
0091
0092 iEvent.put(std::move(tauL2jets));
0093
0094 }
0095
0096 void L1HLTTauMatching::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0097 edm::ParameterSetDescription desc;
0098 desc.add<edm::InputTag>("L1TauTrigger", edm::InputTag("hltL1sDoubleIsoTau40er"))
0099 ->setComment("Name of trigger filter");
0100 desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltSelectedPFTausTrackPt1MediumIsolationReg"))
0101 ->setComment("Input collection of PFTaus");
0102 desc.add<double>("EtMin", 0.0)->setComment("Minimal pT of PFTau to match");
0103 descriptions.setComment(
0104 "This module produces collection of PFTaus matched to L1ExtraTaus/Jets passing a HLT filter (Only p4 and vertex "
0105 "of returned PFTaus are set).");
0106 descriptions.add("L1HLTJetsMatching", desc);
0107 }
0108
0109 #include "FWCore/Framework/interface/MakerMacros.h"
0110 DEFINE_FWK_MODULE(L1HLTTauMatching);