Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:52

0001 #ifndef DQMOffline_Trigger_UtilFuncs_h
0002 #define DQMOffline_Trigger_UtilFuncs_h
0003 
0004 #include "FWCore/Common/interface/TriggerNames.h"
0005 #include "FWCore/Utilities/interface/RegexMatch.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "DataFormats/Common/interface/TriggerResults.h"
0008 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0009 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0010 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0011 #include "DataFormats/Math/interface/deltaR.h"
0012 
0013 namespace hltdqm {
0014   inline bool passTrig(const float objEta,
0015                        float objPhi,
0016                        const trigger::TriggerEvent& trigEvt,
0017                        const std::string& filterName,
0018                        const std::string& processName) {
0019     constexpr float kMaxDR2 = 0.1 * 0.1;
0020 
0021     edm::InputTag filterTag(filterName, "", processName);
0022     trigger::size_type filterIndex = trigEvt.filterIndex(filterTag);
0023     if (filterIndex < trigEvt.sizeFilters()) {  //check that filter is in triggerEvent
0024       const trigger::Keys& trigKeys = trigEvt.filterKeys(filterIndex);
0025       const trigger::TriggerObjectCollection& trigObjColl(trigEvt.getObjects());
0026       for (unsigned short trigKey : trigKeys) {
0027         const trigger::TriggerObject& trigObj = trigObjColl[trigKey];
0028         if (reco::deltaR2(trigObj.eta(), trigObj.phi(), objEta, objPhi) < kMaxDR2)
0029           return true;
0030       }
0031     }
0032     return false;
0033   }
0034 
0035   //empty filters is auto pass
0036   inline bool passTrig(const float objEta,
0037                        float objPhi,
0038                        const trigger::TriggerEvent& trigEvt,
0039                        const std::vector<std::string>& filterNames,
0040                        bool orFilters,
0041                        const std::string& processName) {
0042     if (orFilters) {
0043       if (filterNames.empty())
0044         return true;  //auto pass if empty filters
0045       for (auto& filterName : filterNames) {
0046         if (passTrig(objEta, objPhi, trigEvt, filterName, processName) == true)
0047           return true;
0048       }
0049       return false;
0050     } else {
0051       for (auto& filterName : filterNames) {
0052         if (passTrig(objEta, objPhi, trigEvt, filterName, processName) == false)
0053           return false;
0054       }
0055       return true;
0056     }
0057   }
0058 
0059   //inspired by https://github.com/cms-sw/cmssw/blob/fc4f8bbe1258790e46e2d554aacea15c3e5d9afa/HLTrigger/HLTfilters/src/HLTHighLevel.cc#L124-L165
0060   //triggers are ORed together
0061   //empty pattern is auto pass
0062   inline bool passTrig(const std::string& trigPattern,
0063                        const edm::TriggerNames& trigNames,
0064                        const edm::TriggerResults& trigResults) {
0065     if (trigPattern.empty())
0066       return true;
0067 
0068     std::vector<std::string> trigNamesToPass;
0069     if (edm::is_glob(trigPattern)) {
0070       //matches is vector of string iterators
0071       const auto& matches = edm::regexMatch(trigNames.triggerNames(), trigPattern);
0072       for (auto& name : matches)
0073         trigNamesToPass.push_back(*name);
0074     } else {
0075       trigNamesToPass.push_back(trigPattern);  //not a pattern, much be a path
0076     }
0077     for (auto& trigName : trigNamesToPass) {
0078       size_t pathIndex = trigNames.triggerIndex(trigName);
0079       if (pathIndex < trigResults.size() && trigResults.accept(pathIndex))
0080         return true;
0081     }
0082 
0083     return false;
0084   }
0085 
0086 }  // namespace hltdqm
0087 
0088 #endif