Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:24

0001 #ifndef HLTrigger_HLTcore_HLTPrescaleProvider_h
0002 #define HLTrigger_HLTcore_HLTPrescaleProvider_h
0003 
0004 /** \class HLTPrescaleProvider
0005  *
0006  *  This class provides access routines to get hold of the HLT Configuration,
0007  *  as well as the prescales of Level-1 and High-Level triggers.
0008  *
0009  *  \author Martin Grunewald
0010  *
0011  *  Originally the functions in here were in HLTConfigProvider.
0012  *  The functions that use L1GtUtils and get products from the
0013  *  Event were moved into this class in 2015 when the consumes
0014  *  function calls were added. W. David Dagenhart
0015  */
0016 
0017 #include "HLTrigger/HLTcore/interface/FractionalPrescale.h"
0018 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0019 #include "L1Trigger/GlobalTriggerAnalyzer/interface/L1GtUtils.h"
0020 #include "L1Trigger/L1TGlobal/interface/L1TGlobalUtil.h"
0021 #include "DataFormats/L1TGlobal/interface/GlobalLogicParser.h"
0022 #include "DataFormats/Common/interface/TriggerResults.h"
0023 
0024 #include <memory>
0025 #include <string>
0026 #include <utility>
0027 #include <vector>
0028 #include <type_traits>
0029 
0030 namespace edm {
0031   class ConsumesCollector;
0032   class Event;
0033   class EventSetup;
0034   class ParameterSet;
0035   class Run;
0036   class ParameterSetDescription;
0037 }  // namespace edm
0038 
0039 class HLTPrescaleProvider {
0040 public:
0041   template <typename T>
0042   HLTPrescaleProvider(edm::ParameterSet const& pset, edm::ConsumesCollector&& iC, T& module);
0043 
0044   template <typename T>
0045   HLTPrescaleProvider(edm::ParameterSet const& pset, edm::ConsumesCollector& iC, T& module);
0046 
0047   /// Run-dependent initialisation (non-const method)
0048   ///   "init" return value indicates whether intitialisation has succeeded
0049   ///   "changed" parameter indicates whether the config has actually changed
0050   ///   This must be called at beginRun for most of the other functions in this class to succeed
0051   bool init(const edm::Run& iRun, const edm::EventSetup& iSetup, const std::string& processName, bool& changed);
0052 
0053   HLTConfigProvider const& hltConfigProvider() const { return hltConfigProvider_; }
0054   L1GtUtils const& l1GtUtils() const;
0055   l1t::L1TGlobalUtil const& l1tGlobalUtil() const;
0056 
0057   /// HLT prescale values via (L1) EventSetup
0058   /// current (default) prescale set index - to be taken from L1GtUtil via Event
0059   int prescaleSet(const edm::Event& iEvent, const edm::EventSetup& iSetup);
0060   // negative == error
0061 
0062   /// combining the two methods above
0063   template <typename T = unsigned int>
0064   T prescaleValue(const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& trigger) {
0065     const int set(prescaleSet(iEvent, iSetup));
0066     //there is a template specialisation for unsigned in which returns +1 which
0067     //emulates old behaviour
0068     return set < 0 ? -1 : hltConfigProvider_.prescaleValue<T>(static_cast<unsigned int>(set), trigger);
0069   }
0070 
0071   /// Combined L1T (pair.first) and HLT (pair.second) prescales per HLT path
0072   template <typename TL1 = int, typename THLT = TL1>
0073   std::pair<TL1, THLT> prescaleValues(const edm::Event& iEvent,
0074                                       const edm::EventSetup& iSetup,
0075                                       const std::string& trigger) {
0076     return {convertL1PS<TL1>(getL1PrescaleValue(iEvent, iSetup, trigger)),
0077             prescaleValue<THLT>(iEvent, iSetup, trigger)};
0078   }
0079   // any one negative => error in retrieving this (L1T or HLT) prescale
0080 
0081   // In case of a complex Boolean expression as L1 seed
0082   template <typename TL1 = int, typename THLT = TL1>
0083   std::pair<std::vector<std::pair<std::string, TL1>>, THLT> prescaleValuesInDetail(const edm::Event& iEvent,
0084                                                                                    const edm::EventSetup& iSetup,
0085                                                                                    const std::string& trigger) {
0086     std::pair<std::vector<std::pair<std::string, TL1>>, THLT> retval;
0087     for (auto& entry : getL1PrescaleValueInDetail(iEvent, iSetup, trigger)) {
0088       retval.first.emplace_back(std::move(entry.first), convertL1PS<TL1>(entry.second));
0089     }
0090     retval.second = prescaleValue<THLT>(iEvent, iSetup, trigger);
0091     return retval;
0092   }
0093   // Event rejected by HLTPrescaler on ith HLT path?
0094   bool rejectedByHLTPrescaler(const edm::TriggerResults& triggerResults, unsigned int i) const;
0095   static int l1PrescaleDenominator() { return kL1PrescaleDenominator_; }
0096 
0097   static void fillPSetDescription(edm::ParameterSetDescription& desc,
0098                                   unsigned int stageL1Trigger,
0099                                   edm::InputTag const& l1tAlgBlkInputTag,
0100                                   edm::InputTag const& l1tExtBlkInputTag,
0101                                   bool readPrescalesFromFile);
0102 
0103 private:
0104   static constexpr const char* l1tGlobalDecisionKeyword_ = "L1GlobalDecision";
0105 
0106   void checkL1GtUtils() const;
0107   void checkL1TGlobalUtil() const;
0108 
0109   template <typename T>
0110   T convertL1PS(double val) const {
0111     static_assert(std::is_same_v<T, double> or std::is_same_v<T, FractionalPrescale>,
0112                   "\n\n\tPlease use convertL1PS<double> or convertL1PS<FractionalPrescale>"
0113                   " (other types for L1T prescales are not supported anymore by HLTPrescaleProvider)"
0114                   "\n\tconvertL1PS is used inside prescaleValues and prescaleValuesInDetail,"
0115                   " so it might be necessary to specify template arguments for those calls,"
0116                   "\n\te.g. prescaleValues<double, FractionalPrescale>"
0117                   " (the 1st argument applies to L1T prescales, the 2nd to HLT prescales)\n");
0118     return T(val);
0119   }
0120 
0121   double getL1PrescaleValue(const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::string& trigger);
0122 
0123   std::vector<std::pair<std::string, double>> getL1PrescaleValueInDetail(const edm::Event& iEvent,
0124                                                                          const edm::EventSetup& iSetup,
0125                                                                          const std::string& trigger);
0126 
0127   static constexpr int kL1PrescaleDenominator_ = 100;
0128 
0129   HLTConfigProvider hltConfigProvider_;
0130 
0131   std::unique_ptr<L1GtUtils> l1GtUtils_;
0132   std::unique_ptr<l1t::L1TGlobalUtil> l1tGlobalUtil_;
0133 
0134   unsigned char count_[5] = {0, 0, 0, 0, 0};
0135   bool inited_ = false;
0136 };
0137 
0138 template <typename T>
0139 HLTPrescaleProvider::HLTPrescaleProvider(edm::ParameterSet const& pset, edm::ConsumesCollector&& iC, T& module)
0140     : HLTPrescaleProvider(pset, iC, module) {}
0141 
0142 template <typename T>
0143 HLTPrescaleProvider::HLTPrescaleProvider(edm::ParameterSet const& pset, edm::ConsumesCollector& iC, T& module) {
0144   unsigned int stageL1Trigger = pset.getParameter<unsigned int>("stageL1Trigger");
0145   if (stageL1Trigger <= 1) {
0146     l1GtUtils_ = std::make_unique<L1GtUtils>(pset, iC, false, module, L1GtUtils::UseEventSetupIn::Run);
0147   } else {
0148     l1tGlobalUtil_ = std::make_unique<l1t::L1TGlobalUtil>(pset, iC, module, l1t::UseEventSetupIn::Run);
0149   }
0150 }
0151 
0152 template <>
0153 FractionalPrescale HLTPrescaleProvider::convertL1PS(double val) const;
0154 
0155 #endif  // HLTrigger_HLTcore_HLTPrescaleProvider_h