CMS 3D CMS Logo

TriggerBitChecker.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_Heppy_TriggerBitChecker_h
2 #define PhysicsTools_Heppy_TriggerBitChecker_h
3 
4 #include <vector>
5 #include <string>
6 #include <iostream>
7 #include <cassert>
8 #include <type_traits>
9 
13 
14 namespace heppy {
15 
17  public:
18  struct pathStruct {
19  pathStruct(const std::string &s) : pathName(s), first(0), last(99999999) {}
20  pathStruct() : pathName(), first(0), last(99999999) {}
22  unsigned int first;
23  unsigned int last;
24  };
25 
26  TriggerBitChecker(const std::string &path = "DUMMY");
27  TriggerBitChecker(const std::vector<std::string> &paths);
29 
30  bool check(const edm::EventBase &event, const edm::TriggerResults &result) const;
31 
33  const edm::TriggerResults &result_tr,
34  const pat::PackedTriggerPrescales &result) const;
35 
36  // method templated to force correct choice of output type
37  // (as part of deprecating integer types for trigger prescales)
38  template <typename T = int>
40  const edm::TriggerResults &result_tr,
41  const pat::PackedTriggerPrescales &result) const;
42 
43  private:
44  // list of path name prefixes
45  std::vector<pathStruct> paths_;
46 
48  mutable std::vector<unsigned int> indices_;
49 
51  void syncIndices(const edm::EventBase &event, const edm::TriggerResults &result) const;
53 
55  void rmstar();
56  };
57 
58  template <typename T>
60  const edm::TriggerResults &result_tr,
61  const pat::PackedTriggerPrescales &result) const {
62  static_assert(std::is_same_v<T, double>,
63  "\n\n\tPlease use getprescale<double> "
64  "(other types for trigger prescales are not supported anymore by TriggerBitChecker)");
65  if (result_tr.parameterSetID() != lastID_) {
66  syncIndices(event, result_tr);
67  lastID_ = result_tr.parameterSetID();
68  }
69  if (indices_.empty()) {
70  return -999;
71  }
72  if (indices_.size() > 1) {
73  std::cout << " trying to get prescale for multiple trigger objects at the same time" << std::endl;
74  assert(0);
75  }
76 
77  return result.getPrescaleForIndex<T>(*(indices_.begin()));
78  }
79 
80 } // namespace heppy
81 
82 #endif
void rmstar()
executes a &#39;rm -rf *&#39; in current directory
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)
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_
T getprescale(const edm::EventBase &event, const edm::TriggerResults &result_tr, const pat::PackedTriggerPrescales &result) const
long double T
edm::ParameterSetID lastID_
Definition: event.py:1