CMS 3D CMS Logo

HLTSummaryFilter.cc
Go to the documentation of this file.
1 
11 #include "HLTSummaryFilter.h"
13 
15 
17 
18 //
19 // constructors and destructor
20 //
22  summaryTag_ (iConfig.getParameter<edm::InputTag>("summary")),
23  summaryToken_(consumes<trigger::TriggerEvent>(summaryTag_ )),
24  memberTag_ (iConfig.getParameter<edm::InputTag>("member" )),
25  cut_ (iConfig.getParameter<std::string> ("cut" )),
26  min_N_ (iConfig.getParameter<int> ("minN" )),
27  select_ (cut_ )
28 {
29  edm::LogInfo("HLTSummaryFilter")
30  << "Summary/member/cut/ncut : "
31  << summaryTag_.encode() << " "
32  << memberTag_.encode() << " "
33  << cut_<< " " << min_N_ ;
34 }
35 
37 
38 //
39 // member functions
40 //
41 void
45  // # trigger summary
46  desc.add<edm::InputTag>("summary",edm::InputTag("hltTriggerSummaryAOD","","HLT"));
47  // # filter or collection
48  desc.add<edm::InputTag>("member",edm::InputTag("hlt1jet30","","HLT"));
49  // # cut on trigger object
50  desc.add<std::string>("cut","pt>80");
51  // # min. # of passing objects needed
52  desc.add<int>("minN",1);
53  descriptions.add("hltSummaryFilter", desc);
54 
55 }
56 
57 // ------------ method called to produce the data ------------
58 bool
60 {
61  using namespace std;
62  using namespace edm;
63  using namespace reco;
64  using namespace trigger;
65 
67  iEvent.getByToken(summaryToken_,summary);
68 
69  if (!summary.isValid()) {
70  LogError("HLTSummaryFilter") << "Trigger summary product "
71  << summaryTag_.encode()
72  << " not found! Filter returns false always";
73  return false;
74  }
75 
76  size_type n(0);
77  size_type index(0);
78 
79  // check if we want to look at a filter and its passing physics objects
80  index=summary->filterIndex(memberTag_);
81  if (index<summary->sizeFilters()) {
82  const Keys& KEYS (summary->filterKeys(index));
83  const size_type n1(KEYS.size());
84  for (size_type i=0; i!=n1; ++i) {
85  const TriggerObject& TO( summary->getObjects().at(KEYS[i]) );
86  if (select_(TO)) n++;
87  }
88  const bool accept(n>=min_N_);
89  LogInfo("HLTSummaryFilter")
90  << " Filter objects: " << n << "/" << n1;
91  return accept;
92  }
93 
94  // check if we want to cut on all physics objects of a full "L3" collection
95  index=summary->collectionIndex(memberTag_);
96  if (index<summary->sizeCollections()) {
97  const Keys& KEYS (summary->collectionKeys());
98  const size_type n0 (index == 0? 0 : KEYS.at(index-1));
99  const size_type n1 (KEYS.at(index));
100  for (size_type i=n0; i!=n1; ++i) {
101  const TriggerObject& TO( summary->getObjects().at(i) );
102  if (select_(TO)) n++;
103  }
104  const bool accept(n>=min_N_);
105  LogInfo("HLTSummaryFilter")
106  << " Collection objects: " << n << "/" <<n1-n0;
107  return accept;
108  }
109 
110  // can't help you, bailing out!
111  const bool accept (false);
112  LogInfo("HLTSummaryFilter") << " Default decision: " << accept;
113  return accept;
114 
115 }
116 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< trigger::TriggerEvent > summaryToken_
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:111
trigger::size_type filterIndex(const edm::InputTag &filterTag) const
find index of filter in data-member vector from filter tag
Definition: TriggerEvent.h:123
StringCutObjectSelector< trigger::TriggerObject > select_
HLTSummaryFilter(const edm::ParameterSet &)
edm::InputTag memberTag_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
uint16_t size_type
std::string encode() const
Definition: InputTag.cc:159
const Keys & collectionKeys() const
Definition: TriggerEvent.h:97
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int n0
Definition: AMPTWrapper.h:34
const TriggerObjectCollection & getObjects() const
Definition: TriggerEvent.h:98
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
~HLTSummaryFilter() override
trigger::size_type collectionIndex(const edm::InputTag &collectionTag) const
find index of collection from collection tag
Definition: TriggerEvent.h:114
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
std::vector< size_type > Keys
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::InputTag summaryTag_