File indexing completed on 2025-04-02 23:19:42
0001 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0002 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0003 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidateIsolation.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 #include "FWCore/Utilities/interface/InputTag.h"
0010 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0011
0012 #include <vector>
0013 #include <memory>
0014
0015 namespace edm {
0016 class ConfigurationDescriptions;
0017 }
0018
0019 class HLTEgammaDoubleXGBoostCombFilter : public HLTFilter {
0020 public:
0021 explicit HLTEgammaDoubleXGBoostCombFilter(edm::ParameterSet const&);
0022 ~HLTEgammaDoubleXGBoostCombFilter() override {}
0023
0024 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0025
0026 private:
0027 bool hltFilter(edm::Event& event,
0028 const edm::EventSetup& setup,
0029 trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
0030
0031 const double highMassCut_;
0032 const std::vector<double> leadCutHighMass1_;
0033 const std::vector<double> subCutHighMass1_;
0034 const std::vector<double> leadCutHighMass2_;
0035 const std::vector<double> subCutHighMass2_;
0036 const std::vector<double> leadCutHighMass3_;
0037 const std::vector<double> subCutHighMass3_;
0038
0039 const edm::EDGetTokenT<reco::RecoEcalCandidateCollection> candToken_;
0040 const edm::EDGetTokenT<reco::RecoEcalCandidateIsolationMap> mvaToken_;
0041 };
0042
0043 HLTEgammaDoubleXGBoostCombFilter::HLTEgammaDoubleXGBoostCombFilter(edm::ParameterSet const& config)
0044 : HLTFilter(config),
0045 highMassCut_(config.getParameter<double>("highMassCut")),
0046 leadCutHighMass1_(config.getParameter<std::vector<double>>("leadCutHighMass1")),
0047 subCutHighMass1_(config.getParameter<std::vector<double>>("subCutHighMass1")),
0048 leadCutHighMass2_(config.getParameter<std::vector<double>>("leadCutHighMass2")),
0049 subCutHighMass2_(config.getParameter<std::vector<double>>("subCutHighMass2")),
0050 leadCutHighMass3_(config.getParameter<std::vector<double>>("leadCutHighMass3")),
0051 subCutHighMass3_(config.getParameter<std::vector<double>>("subCutHighMass3")),
0052 candToken_(consumes<reco::RecoEcalCandidateCollection>(config.getParameter<edm::InputTag>("candTag"))),
0053 mvaToken_(consumes<reco::RecoEcalCandidateIsolationMap>(config.getParameter<edm::InputTag>("mvaPhotonTag"))) {}
0054
0055 void HLTEgammaDoubleXGBoostCombFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0056 edm::ParameterSetDescription desc;
0057 makeHLTFilterDescription(desc);
0058
0059 desc.add<double>("highMassCut", 90.0);
0060 desc.add<std::vector<double>>("leadCutHighMass1", {0.92, 0.95});
0061 desc.add<std::vector<double>>("subCutHighMass1", {0.02, 0.04});
0062 desc.add<std::vector<double>>("leadCutHighMass2", {0.85, 0.85});
0063 desc.add<std::vector<double>>("subCutHighMass2", {0.04, 0.08});
0064 desc.add<std::vector<double>>("leadCutHighMass3", {0.30, 0.50});
0065 desc.add<std::vector<double>>("subCutHighMass3", {0.14, 0.20});
0066
0067 desc.add<edm::InputTag>("candTag", edm::InputTag("hltEgammaCandidatesUnseeded"));
0068 desc.add<edm::InputTag>("mvaPhotonTag", edm::InputTag("PhotonXGBoostProducer"));
0069
0070 descriptions.addWithDefaultLabel(desc);
0071 }
0072
0073 bool HLTEgammaDoubleXGBoostCombFilter::hltFilter(edm::Event& event,
0074 const edm::EventSetup& setup,
0075 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0076 const auto& recCollection = event.getHandle(candToken_);
0077 const auto& mvaMap = event.getHandle(mvaToken_);
0078
0079
0080 auto passesHighMassCuts = [&](float leadScore, float subScore, int leadEta, int subEta) {
0081 return (leadScore > leadCutHighMass1_[leadEta] && subScore > subCutHighMass1_[subEta]) ||
0082 (leadScore > leadCutHighMass2_[leadEta] && subScore > subCutHighMass2_[subEta]) ||
0083 (leadScore > leadCutHighMass3_[leadEta] && subScore > subCutHighMass3_[subEta]);
0084 };
0085
0086
0087 auto evaluatePair = [&](const edm::Ref<reco::RecoEcalCandidateCollection>& refi,
0088 const edm::Ref<reco::RecoEcalCandidateCollection>& refj) {
0089 float mvaScorei = (*mvaMap).find(refi)->val;
0090 float mvaScorej = (*mvaMap).find(refj)->val;
0091
0092 int etai = (std::abs(refi->eta()) < 1.5) ? 0 : 1;
0093 int etaj = (std::abs(refj->eta()) < 1.5) ? 0 : 1;
0094
0095 double mass = (refi->p4() + refj->p4()).M();
0096 if (mass < highMassCut_)
0097 return false;
0098
0099 if (mvaScorei >= mvaScorej) {
0100 return passesHighMassCuts(mvaScorei, mvaScorej, etai, etaj);
0101 } else {
0102 return passesHighMassCuts(mvaScorej, mvaScorei, etaj, etai);
0103 }
0104 };
0105
0106
0107 for (size_t i = 0; i < recCollection->size(); ++i) {
0108 edm::Ref<reco::RecoEcalCandidateCollection> refi(recCollection, i);
0109 for (size_t j = i + 1; j < recCollection->size(); ++j) {
0110 edm::Ref<reco::RecoEcalCandidateCollection> refj(recCollection, j);
0111 if (evaluatePair(refi, refj)) {
0112 return true;
0113 }
0114 }
0115 }
0116 return false;
0117 }
0118
0119 #include "FWCore/Framework/interface/MakerMacros.h"
0120 DEFINE_FWK_MODULE(HLTEgammaDoubleXGBoostCombFilter);