File indexing completed on 2023-03-17 11:09:28
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "HLTSummaryFilter.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013
0014 #include "DataFormats/Common/interface/Handle.h"
0015
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017
0018
0019
0020
0021 HLTSummaryFilter::HLTSummaryFilter(const edm::ParameterSet& iConfig)
0022 : HLTFilter(iConfig),
0023 summaryTag_(iConfig.getParameter<edm::InputTag>("summary")),
0024 summaryToken_(consumes<trigger::TriggerEvent>(summaryTag_)),
0025 memberTag_(iConfig.getParameter<edm::InputTag>("member")),
0026 cut_(iConfig.getParameter<std::string>("cut")),
0027 min_N_(iConfig.getParameter<int>("minN")),
0028 select_(cut_) {
0029 edm::LogInfo("HLTSummaryFilter") << "Summary/member/cut/ncut : " << summaryTag_.encode() << " " << memberTag_.encode()
0030 << " " << cut_ << " " << min_N_;
0031 }
0032
0033 HLTSummaryFilter::~HLTSummaryFilter() = default;
0034
0035
0036
0037
0038 void HLTSummaryFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039 edm::ParameterSetDescription desc;
0040 makeHLTFilterDescription(desc);
0041
0042 desc.add<edm::InputTag>("summary", edm::InputTag("hltTriggerSummaryAOD", "", "HLT"));
0043
0044 desc.add<edm::InputTag>("member", edm::InputTag("hlt1jet30", "", "HLT"));
0045
0046 desc.add<std::string>("cut", "pt>80");
0047
0048 desc.add<int>("minN", 1);
0049 descriptions.add("hltSummaryFilter", desc);
0050 }
0051
0052
0053 bool HLTSummaryFilter::hltFilter(edm::Event& iEvent,
0054 const edm::EventSetup& iSetup,
0055 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0056 using namespace std;
0057 using namespace edm;
0058 using namespace reco;
0059 using namespace trigger;
0060
0061 Handle<TriggerEvent> summary;
0062 iEvent.getByToken(summaryToken_, summary);
0063
0064 if (!summary.isValid()) {
0065 LogError("HLTSummaryFilter") << "Trigger summary product " << summaryTag_.encode()
0066 << " not found! Filter returns false always";
0067 return false;
0068 }
0069
0070 size_type n(0);
0071 size_type index(0);
0072
0073
0074 index = summary->filterIndex(memberTag_);
0075 if (index < summary->sizeFilters()) {
0076 const Keys& KEYS(summary->filterKeys(index));
0077 const size_type n1(KEYS.size());
0078 for (size_type i = 0; i != n1; ++i) {
0079 const TriggerObject& TO(summary->getObjects().at(KEYS[i]));
0080 if (select_(TO))
0081 n++;
0082 }
0083 const bool accept(n >= min_N_);
0084 LogInfo("HLTSummaryFilter") << " Filter objects: " << n << "/" << n1;
0085 return accept;
0086 }
0087
0088
0089 index = summary->collectionIndex(memberTag_);
0090 if (index < summary->sizeCollections()) {
0091 const Keys& KEYS(summary->collectionKeys());
0092 const size_type n0(index == 0 ? 0 : KEYS.at(index - 1));
0093 const size_type n1(KEYS.at(index));
0094 for (size_type i = n0; i != n1; ++i) {
0095 const TriggerObject& TO(summary->getObjects().at(i));
0096 if (select_(TO))
0097 n++;
0098 }
0099 const bool accept(n >= min_N_);
0100 LogInfo("HLTSummaryFilter") << " Collection objects: " << n << "/" << n1 - n0;
0101 return accept;
0102 }
0103
0104
0105 const bool accept(false);
0106 LogInfo("HLTSummaryFilter") << " Default decision: " << accept;
0107 return accept;
0108 }
0109
0110 #include "FWCore/Framework/interface/MakerMacros.h"
0111 DEFINE_FWK_MODULE(HLTSummaryFilter);