CMS 3D CMS Logo

TriggerBitChecker.cc
Go to the documentation of this file.
2 
4 #include <cassert>
5 #include <iostream>
6 
7 namespace heppy {
8 
9  TriggerBitChecker::TriggerBitChecker(const std::string &path) : paths_(1, returnPathStruct(path)) { rmstar(); }
10 
11  TriggerBitChecker::TriggerBitChecker(const std::vector<std::string> &paths) : paths_(paths.size()) {
12  for (size_t i = 0; i < paths.size(); ++i)
14  rmstar();
15  }
16 
18  pathStruct newPathStruct(path);
19  if (path[0] > 48 /*'0'*/ && path[0] <= 57 /*'9'*/) {
20  newPathStruct.first = atoi(path.substr(0, path.find('-')).c_str());
21  newPathStruct.last = atoi(path.substr(path.find('-') + 1, path.find(':') - path.find('-') - 1).c_str());
22  newPathStruct.pathName = path.substr(path.find(':') + 1);
23  }
24  return newPathStruct;
25  }
26 
28  if (result.parameterSetID() != lastID_) {
30  lastID_ = result.parameterSetID();
31  }
32  for (std::vector<unsigned int>::const_iterator it = indices_.begin(), ed = indices_.end(); it != ed; ++it) {
33  if (result.accept(*it))
34  return true;
35  }
36  return false;
37  }
38 
40  const edm::TriggerResults &result_tr,
41  const pat::PackedTriggerPrescales &result) const {
42  if (result_tr.parameterSetID() != lastID_) {
43  syncIndices(event, result_tr);
44  lastID_ = result_tr.parameterSetID();
45  }
46  bool outcome = true;
47  for (std::vector<unsigned int>::const_iterator it = indices_.begin(), ed = indices_.end(); it != ed; ++it) {
48  if (result.getPrescaleForIndex(*it) != 1) {
49  outcome = false;
50  break;
51  }
52  }
53  return outcome; // true only if all paths are unprescaled
54  }
55 
57  const edm::TriggerResults &result_tr,
58  const pat::PackedTriggerPrescales &result) const {
59  if (result_tr.parameterSetID() != lastID_) {
60  syncIndices(event, result_tr);
61  lastID_ = result_tr.parameterSetID();
62  }
63  if (indices_.empty()) {
64  // std::cout << " trying to check an inexistent trigger" << std::endl;
65  return -999;
66  }
67  if (indices_.size() > 1) {
68  std::cout << " trying to get prescale for multiple trigger objects at the same time" << std::endl;
69  assert(0);
70  }
71 
72  return result.getPrescaleForIndex(*(indices_.begin()));
73  }
74 
76  indices_.clear();
77  const edm::TriggerNames &names = event.triggerNames(result);
78  std::vector<pathStruct>::const_iterator itp, bgp = paths_.begin(), edp = paths_.end();
79  for (size_t i = 0, n = names.size(); i < n; ++i) {
80  const std::string &thispath = names.triggerName(i);
81  for (itp = bgp; itp != edp; ++itp) {
82  if (thispath.find(itp->pathName) == 0 && event.id().run() >= itp->first && event.id().run() <= itp->last)
83  indices_.push_back(i);
84  }
85  }
86  }
87 
89  std::vector<pathStruct>::iterator itp, bgp = paths_.begin(), edp = paths_.end();
90  for (itp = bgp; itp != edp; ++itp) {
91  std::string::size_type idx = itp->pathName.find('*');
92  if (idx != std::string::npos)
93  itp->pathName.erase(idx);
94  }
95  }
96 } // namespace heppy
void rmstar()
executes a &#39;rm -rf *&#39; in current directory
size
Write out results.
int getprescale(const edm::EventBase &event, const edm::TriggerResults &result_tr, const pat::PackedTriggerPrescales &result) const
bool check_unprescaled(const edm::EventBase &event, const edm::TriggerResults &result_tr, const pat::PackedTriggerPrescales &result) const
pathStruct returnPathStruct(const std::string &path) const
bool check(const edm::EventBase &event, const edm::TriggerResults &result) const
const ParameterSetID & parameterSetID() const
Get stored parameter set id.
assert(be >=bs)
uint16_t size_type
const std::string names[nVars_]
TAKEN FROM http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/ElectroWeakAnalysis/Utilities/src/PdfWeig...
Definition: AlphaT.h:16
std::vector< pathStruct > paths_
TriggerBitChecker(const std::string &path="DUMMY")
void syncIndices(const edm::EventBase &event, const edm::TriggerResults &result) const
sync indices with path names
std::vector< unsigned int > indices_
edm::ParameterSetID lastID_
Definition: event.py:1