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 {
31  using namespace std;
32  using namespace edm;
33 
34  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD configuration: " << endl
35  << " ProcessName = " << processName_ << endl
36  << " TriggerName = " << triggerName_ << endl
37  << " TriggerResultsTag = " << triggerResultsTag_.encode() << endl
38  << " TriggerEventTag = " << triggerEventTag_.encode() << endl;
39 
40 }
41 
43 {
44 }
45 
46 //
47 // member functions
48 //
49 void
52  desc.add<std::string>("processName","HLT");
53  desc.add<std::string>("triggerName","@");
54  desc.add<edm::InputTag>("triggerResults",edm::InputTag("TriggerResults","","HLT"));
55  desc.add<edm::InputTag>("triggerEvent",edm::InputTag("hltTriggerSummaryAOD","","HLT"));
56  descriptions.add("hltEventAnalyzerAOD", desc);
57 }
58 
59 void
60 HLTEventAnalyzerAOD::endRun(edm::Run const & iRun, edm::EventSetup const& iSetup) {}
61 
62 void
64 {
65  using namespace std;
66  using namespace edm;
67 
68  bool changed(true);
69  if (hltPrescaleProvider_.init(iRun,iSetup,processName_,changed)) {
70 
72 
73  if (changed) {
74  // check if trigger name in (new) config
75  if (triggerName_!="@") { // "@" means: analyze all triggers in config
76  const unsigned int n(hltConfig.size());
77  const unsigned int triggerIndex(hltConfig.triggerIndex(triggerName_));
78  if (triggerIndex>=n) {
79  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyze:"
80  << " TriggerName " << triggerName_
81  << " not available in (new) config!" << endl;
82  LogVerbatim("HLTEventAnalyzerAOD") << "Available TriggerNames are: " << endl;
83  hltConfig.dump("Triggers");
84  }
85  }
86  hltConfig.dump("ProcessName");
87  hltConfig.dump("GlobalTag");
88  hltConfig.dump("TableName");
89  hltConfig.dump("Streams");
90  hltConfig.dump("Datasets");
91  hltConfig.dump("PrescaleTable");
92  hltConfig.dump("ProcessPSet");
93  }
94  } else {
95  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyze:"
96  << " config extraction failure with process name "
97  << processName_ << endl;
98  }
99 }
100 
101 // ------------ method called to produce the data ------------
102 void
104 {
105  using namespace std;
106  using namespace edm;
107 
108  LogVerbatim("HLTEventAnalyzerAOD") << endl;
109 
110  // get event products
113  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyze: Error in getting TriggerResults product from Event!" << endl;
114  return;
115  }
117  if (!triggerEventHandle_.isValid()) {
118  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyze: Error in getting TriggerEvent product from Event!" << endl;
119  return;
120  }
121 
123 
124  // sanity check
125  assert(triggerResultsHandle_->size()==hltConfig.size());
126 
127  // analyze this event for the triggers requested
128  if (triggerName_=="@") {
129  const unsigned int n(hltConfig.size());
130  for (unsigned int i=0; i!=n; ++i) {
131  analyzeTrigger(iEvent,iSetup,hltConfig.triggerName(i));
132  }
133  } else {
134  analyzeTrigger(iEvent,iSetup,triggerName_);
135  }
136 
137  return;
138 
139 }
140 
141 void HLTEventAnalyzerAOD::analyzeTrigger(const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& triggerName) {
142 
143  using namespace std;
144  using namespace edm;
145  using namespace reco;
146  using namespace trigger;
147 
148  LogVerbatim("HLTEventAnalyzerAOD") << endl;
149 
151 
152  const unsigned int n(hltConfig.size());
153  const unsigned int triggerIndex(hltConfig.triggerIndex(triggerName));
154  assert(triggerIndex==iEvent.triggerNames(*triggerResultsHandle_).triggerIndex(triggerName));
155 
156  // abort on invalid trigger name
157  if (triggerIndex>=n) {
158  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyzeTrigger: path "
159  << triggerName << " - not found!" << endl;
160  return;
161  }
162 
163  const std::pair<int,int> prescales(hltPrescaleProvider_.prescaleValues(iEvent,iSetup,triggerName));
164  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyzeTrigger: path "
165  << triggerName << " [" << triggerIndex << "] "
166  << "prescales L1T,HLT: " << prescales.first << "," << prescales.second
167  << endl;
168  const std::pair<std::vector<std::pair<std::string,int> >,int> prescalesInDetail(hltPrescaleProvider_.prescaleValuesInDetail(iEvent,iSetup,triggerName));
169  std::ostringstream message;
170  for (unsigned int i=0; i<prescalesInDetail.first.size(); ++i) {
171  message << " " << i << ":" << prescalesInDetail.first[i].first << "/" << prescalesInDetail.first[i].second;
172  }
173  LogVerbatim("HLTEventAnalyzerAOD") << "HLTEventAnalyzerAOD::analyzeTrigger: path "
174  << triggerName << " [" << triggerIndex << "] "
175  << endl
176  << "prescales L1T: " << prescalesInDetail.first.size() << message.str()
177  << endl
178  << " prescale HLT: " << prescalesInDetail.second
179  << endl;
180 
181  // modules on this trigger path
182  const unsigned int m(hltConfig.size(triggerIndex));
183  const vector<string>& moduleLabels(hltConfig.moduleLabels(triggerIndex));
184 
185  // Results from TriggerResults product
186  LogVerbatim("HLTEventAnalyzerAOD") << " Trigger path status:"
187  << " WasRun=" << triggerResultsHandle_->wasrun(triggerIndex)
188  << " Accept=" << triggerResultsHandle_->accept(triggerIndex)
189  << " Error =" << triggerResultsHandle_->error(triggerIndex)
190  << endl;
191  const unsigned int moduleIndex(triggerResultsHandle_->index(triggerIndex));
192  LogVerbatim("HLTEventAnalyzerAOD") << " Last active module - label/type: "
193  << moduleLabels[moduleIndex] << "/" << hltConfig.moduleType(moduleLabels[moduleIndex])
194  << " [" << moduleIndex << " out of 0-" << (m-1) << " on this path]"
195  << endl;
196  assert (moduleIndex<m);
197 
198  // Results from TriggerEvent product - Attention: must look only for
199  // modules actually run in this path for this event!
200  for (unsigned int j=0; j<=moduleIndex; ++j) {
201  const string& moduleLabel(moduleLabels[j]);
202  const string moduleType(hltConfig.moduleType(moduleLabel));
203  // check whether the module is packed up in TriggerEvent product
204  const unsigned int filterIndex(triggerEventHandle_->filterIndex(InputTag(moduleLabel,"",processName_)));
205  if (filterIndex<triggerEventHandle_->sizeFilters()) {
206  LogVerbatim("HLTEventAnalyzerAOD") << " 'L3' filter in slot " << j << " - label/type " << moduleLabel << "/" << moduleType << endl;
207  const Vids& VIDS (triggerEventHandle_->filterIds(filterIndex));
208  const Keys& KEYS(triggerEventHandle_->filterKeys(filterIndex));
209  const size_type nI(VIDS.size());
210  const size_type nK(KEYS.size());
211  assert(nI==nK);
212  const size_type n(max(nI,nK));
213  LogVerbatim("HLTEventAnalyzerAOD") << " " << n << " accepted 'L3' objects found: " << endl;
215  for (size_type i=0; i!=n; ++i) {
216  const TriggerObject& TO(TOC[KEYS[i]]);
217  LogVerbatim("HLTEventAnalyzerAOD") << " " << i << " " << VIDS[i] << "/" << KEYS[i] << ": "
218  << TO.id() << " " << TO.pt() << " " << TO.eta() << " " << TO.phi() << " " << TO.mass()
219  << endl;
220  }
221  }
222  }
223 
224  return;
225 }
unsigned int size() const
number of trigger paths in trigger table
void dump(const std::string &what) const
Dumping config info to cout.
int i
Definition: DBlmapReader.cc:9
bool wasrun() const
Was at least one path run?
virtual void analyzeTrigger(const edm::Event &, const edm::EventSetup &, const std::string &triggerName)
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:234
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:55
const std::string & triggerName(unsigned int triggerIndex) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
float phi() const
Definition: TriggerObject.h:58
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: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
const std::string triggerName_
edm::Handle< trigger::TriggerEvent > triggerEventHandle_
float eta() const
Definition: TriggerObject.h:57
const std::string processName_
module config parameters
uint16_t size_type
std::string encode() const
Definition: InputTag.cc:165
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:110
int iEvent
Definition: GenABIO.cc:230
const edm::InputTag triggerEventTag_
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
const TriggerObjectCollection & getObjects() const
Definition: TriggerEvent.h:98
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.
virtual void beginRun(edm::Run const &, edm::EventSetup const &) override
bool error() const
Has any path encountered an error (exception)
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
const edm::EDGetTokenT< trigger::TriggerEvent > triggerEventToken_
bool isValid() const
Definition: HandleBase.h:75
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:81
virtual 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:59
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
std::vector< int > Vids
Definition: Run.h:42