CMS 3D CMS Logo

HLTEventAnalyzerAOD.cc
Go to the documentation of this file.
1 
16 
17 #include <cassert>
18 
19 //
20 // constructors and destructor
21 //
23  : processName_(ps.getParameter<std::string>("processName")),
24  triggerName_(ps.getParameter<std::string>("triggerName")),
25  triggerResultsTag_(ps.getParameter<edm::InputTag>("triggerResults")),
26  triggerResultsToken_(consumes<edm::TriggerResults>(triggerResultsTag_)),
27  triggerEventTag_(ps.getParameter<edm::InputTag>("triggerEvent")),
28  triggerEventToken_(consumes<trigger::TriggerEvent>(triggerEventTag_)),
29  hltPrescaleProvider_(ps, consumesCollector(), *this) {
30  using namespace std;
31  using namespace edm;
32 
33  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD configuration: " << endl
34  << " ProcessName = " << processName_ << endl
35  << " TriggerName = " << triggerName_ << endl
36  << " TriggerResultsTag = " << triggerResultsTag_.encode() << endl
37  << " TriggerEventTag = " << triggerEventTag_.encode() << endl;
38 }
39 
41 
42 //
43 // member functions
44 //
47  desc.add<std::string>("processName", "HLT");
48  desc.add<std::string>("triggerName", "@");
49  desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults", "", "HLT"));
50  desc.add<edm::InputTag>("triggerEvent", edm::InputTag("hltTriggerSummaryAOD", "", "HLT"));
51  descriptions.add("hltEventAnalyzerAOD", desc);
52 }
53 
54 void HLTEventAnalyzerAOD::endRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {}
55 
56 void HLTEventAnalyzerAOD::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
57  using namespace std;
58  using namespace edm;
59 
60  bool changed(true);
61  if (hltPrescaleProvider_.init(iRun, iSetup, processName_, changed)) {
63 
64  if (changed) {
65  // check if trigger name in (new) config
66  if (triggerName_ != "@") { // "@" means: analyze all triggers in config
67  const unsigned int n(hltConfig.size());
68  const unsigned int triggerIndex(hltConfig.triggerIndex(triggerName_));
69  if (triggerIndex >= n) {
70  LogVerbatim("HLTEventAnalyzerAOD")
71  << "HLTEventAnalyzerAOD::analyze:"
72  << " TriggerName " << triggerName_ << " not available in (new) config!" << endl;
73  LogVerbatim("HLTEventAnalyzerAOD") << "Available TriggerNames are: " << endl;
74  hltConfig.dump("Triggers");
75  }
76  }
77  hltConfig.dump("ProcessName");
78  hltConfig.dump("GlobalTag");
79  hltConfig.dump("TableName");
80  hltConfig.dump("Streams");
81  hltConfig.dump("Datasets");
82  hltConfig.dump("PrescaleTable");
83  hltConfig.dump("ProcessPSet");
84  }
85  } else {
86  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyze:"
87  << " config extraction failure with process name " << processName_ << endl;
88  }
89 }
90 
91 // ------------ method called to produce the data ------------
93  using namespace std;
94  using namespace edm;
95 
96  LogVerbatim("HLTEventAnalyzerAOD") << endl;
97 
98  // get event products
101  LogVerbatim("HLTEventAnalyzerAOD")
102  << "HLTEventAnalyzerAOD::analyze: Error in getting TriggerResults product from Event!" << endl;
103  return;
104  }
106  if (!triggerEventHandle_.isValid()) {
107  LogVerbatim("HLTEventAnalyzerAOD")
108  << "HLTEventAnalyzerAOD::analyze: Error in getting TriggerEvent product from Event!" << endl;
109  return;
110  }
111 
113 
114  // sanity check
115  assert(triggerResultsHandle_->size() == hltConfig.size());
116 
117  // analyze this event for the triggers requested
118  if (triggerName_ == "@") {
119  const unsigned int n(hltConfig.size());
120  for (unsigned int i = 0; i != n; ++i) {
121  analyzeTrigger(iEvent, iSetup, hltConfig.triggerName(i));
122  }
123  } else {
124  analyzeTrigger(iEvent, iSetup, triggerName_);
125  }
126 
127  return;
128 }
129 
131  const edm::EventSetup& iSetup,
132  const std::string& triggerName) {
133  using namespace std;
134  using namespace edm;
135  using namespace reco;
136  using namespace trigger;
137 
138  LogVerbatim("HLTEventAnalyzerAOD") << endl;
139 
141 
142  const unsigned int n(hltConfig.size());
143  const unsigned int triggerIndex(hltConfig.triggerIndex(triggerName));
144  assert(triggerIndex == iEvent.triggerNames(*triggerResultsHandle_).triggerIndex(triggerName));
145 
146  // abort on invalid trigger name
147  if (triggerIndex >= n) {
148  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyzeTrigger: path " << triggerName << " - not found!"
149  << endl;
150  return;
151  }
152 
153  const std::pair<int, int> prescales(hltPrescaleProvider_.prescaleValues(iEvent, iSetup, triggerName));
154  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyzeTrigger: path " << triggerName << " ["
155  << triggerIndex << "] "
156  << "prescales L1T,HLT: " << prescales.first << "," << prescales.second << endl;
157  const std::pair<std::vector<std::pair<std::string, int> >, int> prescalesInDetail(
158  hltPrescaleProvider_.prescaleValuesInDetail(iEvent, iSetup, triggerName));
159  std::ostringstream message;
160  for (unsigned int i = 0; i < prescalesInDetail.first.size(); ++i) {
161  message << " " << i << ":" << prescalesInDetail.first[i].first << "/" << prescalesInDetail.first[i].second;
162  }
163  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyzeTrigger: path " << triggerName << " ["
164  << triggerIndex << "] " << endl
165  << "prescales L1T: " << prescalesInDetail.first.size() << message.str() << endl
166  << " prescale HLT: " << prescalesInDetail.second << endl;
167 
168  // modules on this trigger path
169  const unsigned int m(hltConfig.size(triggerIndex));
170  const vector<string>& moduleLabels(hltConfig.moduleLabels(triggerIndex));
171 
172  // Results from TriggerResults product
173  LogVerbatim("HLTEventAnalyzerAOD") << " Trigger path status:"
174  << " WasRun=" << triggerResultsHandle_->wasrun(triggerIndex)
175  << " Accept=" << triggerResultsHandle_->accept(triggerIndex)
176  << " Error =" << triggerResultsHandle_->error(triggerIndex) << endl;
177  const unsigned int moduleIndex(triggerResultsHandle_->index(triggerIndex));
178  LogVerbatim("HLTEventAnalyzerAOD") << " Last active module - label/type: " << moduleLabels[moduleIndex] << "/"
179  << hltConfig.moduleType(moduleLabels[moduleIndex]) << " [" << moduleIndex
180  << " out of 0-" << (m - 1) << " on this path]" << endl;
181  assert(moduleIndex < m);
182 
183  // Results from TriggerEvent product - Attention: must look only for
184  // modules actually run in this path for this event!
185  for (unsigned int j = 0; j <= moduleIndex; ++j) {
186  const string& moduleLabel(moduleLabels[j]);
187  const string moduleType(hltConfig.moduleType(moduleLabel));
188  // check whether the module is packed up in TriggerEvent product
189  const unsigned int filterIndex(triggerEventHandle_->filterIndex(InputTag(moduleLabel, "", processName_)));
190  if (filterIndex < triggerEventHandle_->sizeFilters()) {
191  LogVerbatim("HLTEventAnalyzerAOD") << " 'L3' filter in slot " << j << " - label/type " << moduleLabel << "/"
192  << moduleType << endl;
193  const Vids& VIDS(triggerEventHandle_->filterIds(filterIndex));
194  const Keys& KEYS(triggerEventHandle_->filterKeys(filterIndex));
195  const size_type nI(VIDS.size());
196  const size_type nK(KEYS.size());
197  assert(nI == nK);
198  const size_type n(max(nI, nK));
199  LogVerbatim("HLTEventAnalyzerAOD") << " " << n << " accepted 'L3' objects found: " << endl;
201  for (size_type i = 0; i != n; ++i) {
202  const TriggerObject& TO(TOC[KEYS[i]]);
203  LogVerbatim("HLTEventAnalyzerAOD") << " " << i << " " << VIDS[i] << "/" << KEYS[i] << ": " << TO.id() << " "
204  << TO.pt() << " " << TO.eta() << " " << TO.phi() << " " << TO.mass() << endl;
205  }
206  }
207  }
208 
209  return;
210 }
unsigned int size() const
number of trigger paths in trigger table
void dump(const std::string &what) const
Dumping config info to cout.
bool wasrun() const
Was at least one path run?
virtual void analyzeTrigger(const edm::Event &, const edm::EventSetup &, const std::string &triggerName)
std::pair< std::vector< std::pair< std::string, int > >, int > prescaleValuesInDetail(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
const std::string moduleType(const std::string &module) const
C++ class name of module.
int id() const
getters
Definition: TriggerObject.h:51
const std::string & triggerName(unsigned int triggerIndex) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
float phi() const
Definition: TriggerObject.h:54
const edm::InputTag triggerResultsTag_
const edm::EDGetTokenT< edm::TriggerResults > triggerResultsToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool accept() const
Has at least one path accepted the event?
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:118
trigger::size_type filterIndex(const edm::InputTag &filterTag) const
find index of filter in data-member vector from filter tag
Definition: TriggerEvent.h:132
const std::string triggerName_
edm::Handle< trigger::TriggerEvent > triggerEventHandle_
float eta() const
Definition: TriggerObject.h:53
const std::string processName_
module config parameters
uint16_t size_type
std::string encode() const
Definition: InputTag.cc:159
HLTPrescaleProvider hltPrescaleProvider_
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
const Vids & filterIds(trigger::size_type index) const
Definition: TriggerEvent.h:117
int iEvent
Definition: GenABIO.cc:224
const edm::InputTag triggerEventTag_
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:24
const TriggerObjectCollection & getObjects() const
Definition: TriggerEvent.h:101
unsigned int size() const
Get number of paths stored.
unsigned int index(const unsigned int i) const
Get index (slot position) of module giving the decision of the ith path.
void beginRun(edm::Run const &, edm::EventSetup const &) override
bool error() const
Has any path encountered an error (exception)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void analyze(const edm::Event &, const edm::EventSetup &) override
const edm::EDGetTokenT< trigger::TriggerEvent > triggerEventToken_
bool isValid() const
Definition: HandleBase.h:70
const std::vector< std::string > & moduleLabels(unsigned int trigger) const
label(s) of module(s) on a trigger path
edm::Handle< edm::TriggerResults > triggerResultsHandle_
additional class data memebers
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:75
~HLTEventAnalyzerAOD() override
void endRun(edm::Run const &, edm::EventSetup const &) override
std::vector< size_type > Keys
HLTEventAnalyzerAOD(const edm::ParameterSet &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::pair< int, int > prescaleValues(const edm::Event &iEvent, const edm::EventSetup &iSetup, const std::string &trigger)
Combined L1T (pair.first) and HLT (pair.second) prescales per HLT path.
HLTConfigProvider const & hltConfigProvider() const
fixed size matrix
HLT enums.
float mass() const
Definition: TriggerObject.h:55
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:265
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
std::vector< int > Vids
Definition: Run.h:45