CMS 3D CMS Logo

UtilFuncs.h
Go to the documentation of this file.
1 #ifndef DQMOffline_Trigger_UtilFuncs_h
2 #define DQMOffline_Trigger_UtilFuncs_h
3 
12 
13 namespace hltdqm {
14  bool passTrig(const float objEta,float objPhi, const trigger::TriggerEvent& trigEvt,
16  constexpr float kMaxDR2 = 0.1*0.1;
17 
18  edm::InputTag filterTag(filterName,"",processName);
19  trigger::size_type filterIndex = trigEvt.filterIndex(filterTag);
20  if(filterIndex<trigEvt.sizeFilters()){ //check that filter is in triggerEvent
21  const trigger::Keys& trigKeys = trigEvt.filterKeys(filterIndex);
22  const trigger::TriggerObjectCollection & trigObjColl(trigEvt.getObjects());
23  for(unsigned short trigKey : trigKeys){
24  const trigger::TriggerObject& trigObj = trigObjColl[trigKey];
25  if(reco::deltaR2(trigObj.eta(),trigObj.phi(),objEta,objPhi)<kMaxDR2) return true;
26  }
27  }
28  return false;
29  }
30 
31  //empty filters is auto pass
32  bool passTrig(const float objEta,float objPhi, const trigger::TriggerEvent& trigEvt,
33  const std::vector<std::string>& filterNames,bool orFilters,const std::string& processName){
34  if(orFilters){
35  if(filterNames.empty()) return true; //auto pass if empty filters
36  for(auto& filterName : filterNames){
37  if(passTrig(objEta,objPhi,trigEvt,filterName,processName)==true) return true;
38  }
39  return false;
40  }else{
41  for(auto& filterName : filterNames){
42  if(passTrig(objEta,objPhi,trigEvt,filterName,processName)==false) return false;
43  }
44  return true;
45  }
46  }
47 
48  //inspired by https://github.com/cms-sw/cmssw/blob/fc4f8bbe1258790e46e2d554aacea15c3e5d9afa/HLTrigger/HLTfilters/src/HLTHighLevel.cc#L124-L165
49  //triggers are ORed together
50  //empty pattern is auto pass
51  bool passTrig(const std::string& trigPattern,const edm::TriggerNames& trigNames,const edm::TriggerResults& trigResults){
52 
53  if(trigPattern.empty()) return true;
54 
55  std::vector<std::string> trigNamesToPass;
56  if(edm::is_glob(trigPattern)){
57 
58  //matches is vector of string iterators
59  const auto& matches = edm::regexMatch(trigNames.triggerNames(),trigPattern);
60  for(auto& name : matches) trigNamesToPass.push_back(*name);
61  }else{
62  trigNamesToPass.push_back(trigPattern); //not a pattern, much be a path
63  }
64  for(auto& trigName : trigNamesToPass){
65  size_t pathIndex = trigNames.triggerIndex(trigName);
66  if(pathIndex<trigResults.size() &&
67  trigResults.accept(pathIndex)) return true;
68  }
69 
70  return false;
71  }
72 
73 
74 }
75 
76 #endif
The single EDProduct to be saved for each event (AOD case)
Definition: TriggerEvent.h:25
trigger::size_type sizeFilters() const
Definition: TriggerEvent.h:135
float phi() const
Definition: TriggerObject.h:58
bool accept() const
Has at least one path accepted the event?
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:111
bool is_glob(std::string const &pattern)
Definition: RegexMatch.cc:18
trigger::size_type filterIndex(const edm::InputTag &filterTag) const
find index of filter in data-member vector from filter tag
Definition: TriggerEvent.h:123
float eta() const
Definition: TriggerObject.h:57
uint16_t size_type
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
#define constexpr
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:22
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.
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:81
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
std::vector< size_type > Keys
std::vector< std::vector< std::string >::const_iterator > regexMatch(std::vector< std::string > const &strings, std::regex const &regexp)
Definition: RegexMatch.cc:30
bool passTrig(const float objEta, float objPhi, const trigger::TriggerEvent &trigEvt, const std::string &filterName, const std::string &processName)
Definition: UtilFuncs.h:14