Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:36:36

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   //producer collection (hltEgammaCandidates(Unseeded))
0077   const auto& recCollection = event.getHandle(candToken_);
0078 
0079   //get hold of photon MVA association map
0080   const auto& mvaMap = event.getHandle(mvaToken_);
0081 
0082   std::vector<math::XYZTLorentzVector> p4s(recCollection->size());
0083   std::vector<bool> isTight(recCollection->size());
0084 
0085   bool accept = false;
0086 
0087   for (size_t i = 0; i < recCollection->size(); i++) {
0088     edm::Ref<reco::RecoEcalCandidateCollection> refi(recCollection, i);
0089     float EtaSCi = refi->eta();
0090     int eta1 = (std::abs(EtaSCi) < 1.5) ? 0 : 1;
0091     float mvaScorei = (*mvaMap).find(refi)->val;
0092     math::XYZTLorentzVector p4i = refi->p4();
0093     for (size_t j = i + 1; j < recCollection->size(); j++) {
0094       edm::Ref<reco::RecoEcalCandidateCollection> refj(recCollection, j);
0095       float EtaSCj = refj->eta();
0096       int eta2 = (std::abs(EtaSCj) < 1.5) ? 0 : 1;
0097       float mvaScorej = (*mvaMap).find(refj)->val;
0098       math::XYZTLorentzVector p4j = refj->p4();
0099       math::XYZTLorentzVector pairP4 = p4i + p4j;
0100       double mass = pairP4.M();
0101       if (mass >= highMassCut_) {
0102         if (mvaScorei >= mvaScorej && ((mvaScorei > leadCutHighMass1_[eta1] && mvaScorej > subCutHighMass1_[eta2]) ||
0103                                        (mvaScorei > leadCutHighMass2_[eta1] && mvaScorej > subCutHighMass2_[eta2]) ||
0104                                        (mvaScorei > leadCutHighMass3_[eta1] && mvaScorej > subCutHighMass3_[eta2]))) {
0105           accept = true;
0106         }  //if scoreI > scoreJ
0107         else if (mvaScorej > mvaScorei &&
0108                  ((mvaScorej > leadCutHighMass1_[eta1] && mvaScorei > subCutHighMass1_[eta2]) ||
0109                   (mvaScorej > leadCutHighMass2_[eta1] && mvaScorei > subCutHighMass2_[eta2]) ||
0110                   (mvaScorej > leadCutHighMass3_[eta1] && mvaScorei > subCutHighMass3_[eta2]))) {
0111           accept = true;
0112         }  // if scoreJ > scoreI
0113       }  //If high mass
0114     }  //j loop
0115   }  //i loop
0116   return accept;
0117 }  //Definition
0118 
0119 #include "FWCore/Framework/interface/MakerMacros.h"
0120 DEFINE_FWK_MODULE(HLTEgammaDoubleXGBoostCombFilter);