File indexing completed on 2024-04-06 12:27:43
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012 #include <algorithm>
0013 #include <memory>
0014 #include <set>
0015 #include <utility>
0016
0017 #include "FWCore/Framework/interface/global/EDProducer.h"
0018 #include "FWCore/Framework/interface/Event.h"
0019 #include "FWCore/Framework/interface/EventSetup.h"
0020 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0023 #include "FWCore/Utilities/interface/InputTag.h"
0024 #include "FWCore/Utilities/interface/Exception.h"
0025 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0026 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0027 #include "DataFormats/TauReco/interface/PFTau.h"
0028 #include "DataFormats/Math/interface/deltaR.h"
0029 #include "CommonTools/Utils/interface/PtComparator.h"
0030
0031 class HLTDiPFJetPlusTausCandidatePFJetProducer : public edm::global::EDProducer<> {
0032 public:
0033 explicit HLTDiPFJetPlusTausCandidatePFJetProducer(const edm::ParameterSet&);
0034 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0035 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0036
0037 private:
0038 const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> tauSrc_;
0039 const edm::EDGetTokenT<reco::PFJetCollection> pfJetSrc_;
0040 const double extraTauPtCut_;
0041 const double mjjMin_, m2jjMin_;
0042 const double dRmin_, dRmin2_;
0043 GreaterByPt<reco::PFJet> pTComparator_;
0044 };
0045
0046 HLTDiPFJetPlusTausCandidatePFJetProducer::HLTDiPFJetPlusTausCandidatePFJetProducer(const edm::ParameterSet& iConfig)
0047 : tauSrc_(consumes(iConfig.getParameter<edm::InputTag>("tauSrc"))),
0048 pfJetSrc_(consumes(iConfig.getParameter<edm::InputTag>("pfJetSrc"))),
0049 extraTauPtCut_(iConfig.getParameter<double>("extraTauPtCut")),
0050 mjjMin_(iConfig.getParameter<double>("mjjMin")),
0051 m2jjMin_(mjjMin_ * mjjMin_),
0052 dRmin_(iConfig.getParameter<double>("dRmin")),
0053 dRmin2_(dRmin_ * dRmin_) {
0054 if (dRmin_ <= 0.) {
0055 throw cms::Exception("HLTDiPFJetPlusTausCandidatePFJetProducer")
0056 << "invalid value for parameter \"dRmin\" (must be > 0): " << dRmin_;
0057 }
0058 produces<reco::PFJetCollection>();
0059 }
0060
0061 void HLTDiPFJetPlusTausCandidatePFJetProducer::produce(edm::StreamID,
0062 edm::Event& iEvent,
0063 const edm::EventSetup&) const {
0064 auto const& pfJets = iEvent.get(pfJetSrc_);
0065
0066 auto cleanedPFJets = std::make_unique<reco::PFJetCollection>();
0067 cleanedPFJets->reserve(pfJets.size());
0068
0069 trigger::VRpftau taus;
0070 iEvent.get(tauSrc_).getObjects(trigger::TriggerTau, taus);
0071
0072 std::set<unsigned int> indices;
0073
0074 for (unsigned int iJet1 = 0; iJet1 < pfJets.size(); iJet1++) {
0075 for (unsigned int iJet2 = iJet1 + 1; iJet2 < pfJets.size(); iJet2++) {
0076 bool correctComb = false;
0077 const reco::PFJet& myPFJet1 = pfJets[iJet1];
0078 const reco::PFJet& myPFJet2 = pfJets[iJet2];
0079
0080 if (mjjMin_ >= 0. && (myPFJet1.p4() + myPFJet2.p4()).M2() < m2jjMin_)
0081 continue;
0082
0083 for (unsigned int iTau1 = 0; iTau1 < taus.size(); iTau1++) {
0084 if (reco::deltaR2(taus[iTau1]->p4(), myPFJet1.p4()) < dRmin2_)
0085 continue;
0086 if (reco::deltaR2(taus[iTau1]->p4(), myPFJet2.p4()) < dRmin2_)
0087 continue;
0088
0089 if (taus.size() == 1) {
0090 if (taus[iTau1]->pt() < extraTauPtCut_)
0091 continue;
0092
0093 correctComb = true;
0094 } else {
0095 for (unsigned int iTau2 = iTau1 + 1; iTau2 < taus.size(); iTau2++) {
0096 if (taus[iTau1]->pt() < extraTauPtCut_ && taus[iTau2]->pt() < extraTauPtCut_)
0097 continue;
0098
0099 if (reco::deltaR2(taus[iTau2]->p4(), myPFJet1.p4()) < dRmin2_)
0100 continue;
0101 if (reco::deltaR2(taus[iTau2]->p4(), myPFJet2.p4()) < dRmin2_)
0102 continue;
0103
0104 correctComb = true;
0105 break;
0106 }
0107 }
0108 if (correctComb)
0109 break;
0110 }
0111
0112 if (correctComb) {
0113 indices.insert(iJet1);
0114 indices.insert(iJet2);
0115 }
0116 }
0117
0118 for (const auto& i : indices)
0119 cleanedPFJets->emplace_back(pfJets[i]);
0120 }
0121
0122 std::sort(cleanedPFJets->begin(), cleanedPFJets->end(), pTComparator_);
0123 iEvent.put(std::move(cleanedPFJets));
0124 }
0125
0126 void HLTDiPFJetPlusTausCandidatePFJetProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0127 edm::ParameterSetDescription desc;
0128 desc.add<edm::InputTag>("pfJetSrc", edm::InputTag("hltAK4PFJetsCorrected"))->setComment("Input collection of PFJets");
0129 desc.add<edm::InputTag>("tauSrc", edm::InputTag("hltSinglePFTau20TrackPt1LooseChargedIsolationReg"))
0130 ->setComment("Input collection of PFTaus that have passed ID and isolation requirements");
0131 desc.add<double>("extraTauPtCut", 45)->setComment("In case of asymmetric tau pt cuts");
0132 desc.add<double>("mjjMin", 500)->setComment("VBF dijet mass condition");
0133 desc.add<double>("dRmin", 0.5)->setComment("Minimum dR between PFJets and filtered PFTaus");
0134 descriptions.setComment(
0135 "This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT "
0136 "filter.");
0137 descriptions.addWithDefaultLabel(desc);
0138 }
0139
0140 #include "FWCore/Framework/interface/MakerMacros.h"
0141 DEFINE_FWK_MODULE(HLTDiPFJetPlusTausCandidatePFJetProducer);