File indexing completed on 2024-04-06 12:18:41
0001 #include "HLTMultipletFilter.h"
0002
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/Common/interface/RefToBase.h"
0005
0006 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0007 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
0008
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0012
0013 #include <cmath>
0014
0015 HLTMultipletFilter::HLTMultipletFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0016 hltEGammaSeedLabel_ = iConfig.getParameter<edm::InputTag>("L1EGammaInputTag");
0017 hltEtSumSeedLabel_ = iConfig.getParameter<edm::InputTag>("L1EtSumInputTag");
0018 hltJetSeedLabel_ = iConfig.getParameter<edm::InputTag>("L1JetInputTag");
0019 hltMuonSeedLabel_ = iConfig.getParameter<edm::InputTag>("L1MuonInputTag");
0020 hltTauSeedLabel_ = iConfig.getParameter<edm::InputTag>("L1TauInputTag");
0021 minN_ = iConfig.getParameter<int>("MinN");
0022 ibxMin_ = iConfig.getParameter<int>("IBxMin");
0023 ibxMax_ = iConfig.getParameter<int>("IBxMax");
0024 minEta_ = iConfig.getParameter<double>("MinEta");
0025 maxEta_ = iConfig.getParameter<double>("MaxEta");
0026 minPhi_ = iConfig.getParameter<double>("MinPhi");
0027 maxPhi_ = iConfig.getParameter<double>("MaxPhi");
0028 minPt_ = iConfig.getParameter<double>("MinPt");
0029
0030 if (hltEGammaSeedLabel_ == edm::InputTag()) {
0031 flag_[EGamma] = false;
0032 } else {
0033 flag_[EGamma] = true;
0034 hltEGammaToken_ = consumes<l1t::EGammaBxCollection>(hltEGammaSeedLabel_);
0035 }
0036 if (hltEtSumSeedLabel_ == edm::InputTag()) {
0037 flag_[EtSum] = false;
0038 } else {
0039 flag_[EtSum] = true;
0040 hltEtSumToken_ = consumes<l1t::EtSumBxCollection>(hltEtSumSeedLabel_);
0041 }
0042 if (hltJetSeedLabel_ == edm::InputTag()) {
0043 flag_[Jet] = false;
0044 } else {
0045 flag_[Jet] = true;
0046 hltJetToken_ = consumes<l1t::JetBxCollection>(hltJetSeedLabel_);
0047 }
0048 if (hltMuonSeedLabel_ == edm::InputTag()) {
0049 flag_[Muon] = false;
0050 } else {
0051 flag_[Muon] = true;
0052 hltMuonToken_ = consumes<l1t::MuonBxCollection>(hltMuonSeedLabel_);
0053 }
0054 if (hltTauSeedLabel_ == edm::InputTag()) {
0055 flag_[Tau] = false;
0056 } else {
0057 flag_[Tau] = true;
0058 hltTauToken_ = consumes<l1t::TauBxCollection>(hltTauSeedLabel_);
0059 }
0060 edm::LogVerbatim("Report") << "Input Parameters:: minN = " << minN_ << " Bx Range = " << ibxMin_ << ":" << ibxMax_
0061 << " minPt = " << minPt_ << " Eta " << minEta_ << ":" << maxEta_ << " Phi " << minPhi_
0062 << ":" << maxPhi_ << " GT Seed for EGamma " << hltEGammaSeedLabel_ << ", EtSum "
0063 << hltEtSumSeedLabel_ << ", Jet " << hltJetSeedLabel_ << ", Muon " << hltMuonSeedLabel_
0064 << ", and Tau " << hltTauSeedLabel_ << std::endl;
0065 }
0066
0067 HLTMultipletFilter::~HLTMultipletFilter() = default;
0068
0069 void HLTMultipletFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0070 edm::ParameterSetDescription desc;
0071 makeHLTFilterDescription(desc);
0072 desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag());
0073 desc.add<edm::InputTag>("L1EtSumInputTag", edm::InputTag());
0074 desc.add<edm::InputTag>("L1JetInputTag", edm::InputTag("hltCaloStage2Digis:Jet"));
0075 desc.add<edm::InputTag>("L1MuonInputTag", edm::InputTag());
0076 desc.add<edm::InputTag>("L1TauInputTag", edm::InputTag("hltCaloStage2Digis:Tau"));
0077 desc.add<int>("MinN", 1);
0078 desc.add<int>("IBxMin", 0);
0079 desc.add<int>("IBxMax", 0);
0080 desc.add<double>("MinEta", 1.305);
0081 desc.add<double>("MaxEta", 3.000);
0082 desc.add<double>("MinPhi", 5.4105);
0083 desc.add<double>("MaxPhi", 5.5796);
0084 desc.add<double>("MinPt", 20.0);
0085 descriptions.add("hltMultipletFilter", desc);
0086 }
0087
0088 bool HLTMultipletFilter::hltFilter(edm::Event& iEvent,
0089 const edm::EventSetup& iSetup,
0090 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0091
0092 if (saveTags()) {
0093 if (flag_[EGamma])
0094 filterproduct.addCollectionTag(hltEGammaSeedLabel_);
0095 if (flag_[EtSum])
0096 filterproduct.addCollectionTag(hltEtSumSeedLabel_);
0097 if (flag_[Jet])
0098 filterproduct.addCollectionTag(hltJetSeedLabel_);
0099 if (flag_[Muon])
0100 filterproduct.addCollectionTag(hltMuonSeedLabel_);
0101 if (flag_[Tau])
0102 filterproduct.addCollectionTag(hltTauSeedLabel_);
0103 }
0104
0105 int nobj(0);
0106
0107 if (flag_[EGamma]) {
0108 nobj += objects(iEvent, hltEGammaToken_, hltEGammaSeedLabel_, EGamma);
0109 if (nobj >= minN_)
0110 return true;
0111 }
0112 if (flag_[EtSum]) {
0113 nobj += objects(iEvent, hltEtSumToken_, hltEtSumSeedLabel_, EtSum);
0114 if (nobj >= minN_)
0115 return true;
0116 }
0117 if (flag_[Jet]) {
0118 nobj += objects(iEvent, hltJetToken_, hltJetSeedLabel_, Jet);
0119 if (nobj >= minN_)
0120 return true;
0121 }
0122 if (flag_[Muon]) {
0123 nobj += objects(iEvent, hltMuonToken_, hltMuonSeedLabel_, Muon);
0124 if (nobj >= minN_)
0125 return true;
0126 }
0127 if (flag_[Tau]) {
0128 nobj += objects(iEvent, hltTauToken_, hltTauSeedLabel_, Tau);
0129 if (nobj >= minN_)
0130 return true;
0131 }
0132 return false;
0133 }
0134
0135 template <typename T1>
0136 int HLTMultipletFilter::objects(edm::Event& iEvent,
0137 edm::EDGetTokenT<T1> const& hltToken,
0138 edm::InputTag const& hltSeedLabel,
0139 HLTMultipletFilter::Types type) const {
0140 int nobj(0);
0141 edm::Handle<T1> objs;
0142 iEvent.getByToken(hltToken, objs);
0143 if (!objs.isValid()) {
0144 edm::LogWarning("Report") << "Collection with input tag " << hltSeedLabel
0145 << " requested, but not found in the event.";
0146 } else {
0147 edm::LogVerbatim("Report") << "Collection for type " << type << " has " << objs->size() << " in " << ibxMin_ << ":"
0148 << ibxMax_ << " BX's";
0149 for (int ibx = ibxMin_; ibx <= ibxMax_; ++ibx) {
0150 for (auto p = objs->begin(ibx); p != objs->end(ibx); ++p) {
0151 if (p->pt() > minPt_) {
0152 if (p->eta() > minEta_ && p->eta() < maxEta_) {
0153 double phi = p->phi();
0154 if (phi < 0)
0155 phi += 2 * M_PI;
0156 if (phi > minPhi_ && phi < maxPhi_)
0157 ++nobj;
0158 }
0159 }
0160 }
0161 }
0162 }
0163 return nobj;
0164 }
0165
0166
0167 #include "FWCore/Framework/interface/MakerMacros.h"
0168 DEFINE_FWK_MODULE(HLTMultipletFilter);