Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161
#ifndef DQMOFFLINE_TRIGGER_EGHLTTRIGTOOLS
#define DQMOFFLINE_TRIGGER_EGHLTTRIGTOOLS

#include "DQMOffline/Trigger/interface/EgHLTTrigCodes.h"
#include "DataFormats/HLTReco/interface/TriggerEvent.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Math/interface/deltaR.h"
#include <vector>

class HLTConfigProvider;

namespace egHLT {

  namespace trigTools {
    TrigCodes::TrigBitSet getFiltersPassed(const std::vector<std::pair<std::string, int> >& filters,
                                           const trigger::TriggerEvent* trigEvt,
                                           const std::string& hltTag,
                                           const TrigCodes& trigCodes);

    template <class T>
    void setFiltersObjPasses(std::vector<T>& objs,
                             const std::vector<std::string>& filters,
                             const std::vector<std::pair<std::string, std::string> >& l1PreAndSeedFilters,
                             const TrigCodes::TrigBitSet& evtTrigBits,
                             const TrigCodes& trigCodes,
                             const trigger::TriggerEvent* trigEvt,
                             const std::string& hltTag);

    template <class T, class U>
    void fillHLTposition(T& obj,
                         U& hltData,
                         const std::vector<std::string>& filters,
                         const trigger::TriggerEvent* trigEvt,
                         const std::string& hltTag);
    std::vector<int> getMinNrObjsRequiredByFilter(
        const std::vector<std::string>& filterName);  //slow function, call at begin job and cache results

    //reads hlt config and works out which are the active last filters stored in trigger summary, is sorted
    void getActiveFilters(const HLTConfigProvider& hltConfig,
                          std::vector<std::string>& activeFilters,
                          std::vector<std::string>& activeEleFilters,
                          std::vector<std::string>& activeEle2LegFilters,
                          std::vector<std::string>& activePhoFilters,
                          std::vector<std::string>& activePho2LegFilters);
    //---Morse test--------
    //void getPhoton30(const HLTConfigProvider& hltConfig,std::vector<std::string>& activeFilters);
    //------------------
    //filters a list of filternames removing any filters which are not in active filters, assumes active filters is sorted
    void filterInactiveTriggers(std::vector<std::string>& namesToFilter, std::vector<std::string>& activeFilters);
    //filters a list of filterName1:filterName2 removing any entry for which either filter is not in activeFilters, assumes active filters is sorted
    void filterInactiveTightLooseTriggers(std::vector<std::string>& namesToFilter,
                                          const std::vector<std::string>& activeFilters);

    void translateFiltersToPathNames(const HLTConfigProvider& hltConfig,
                                     const std::vector<std::string>& filters,
                                     std::vector<std::string>& paths);
    std::string getL1SeedFilterOfPath(const HLTConfigProvider& hltConfig, const std::string& path);

    //looks for string Et and then looks for a number after that (currently the standard of all E/g triggers)
    //returns 0 if unsuccessful
    float getEtThresFromName(const std::string& trigName);
    float getSecondEtThresFromName(const std::string& trigName);
  }  // namespace trigTools

  //I have the horrible feeling that I'm converting into an intermediatry format and then coverting back again
  //Okay how this works
  //1) create a TrigBitSet for each particle set to 0 initally
  //2) loop over each filter, for each particle that passes the filter, set the appropriate bit in the TrigBitSet
  //3) after that, loop over each particle setting that its TrigBitSet which has been calculated
  //4) because L1 pre-scaled paths are special, we only set those if an event wide trigger has been set
  template <class T>
  void trigTools::setFiltersObjPasses(std::vector<T>& particles,
                                      const std::vector<std::string>& filters,
                                      const std::vector<std::pair<std::string, std::string> >& l1PreAndSeedFilters,
                                      const TrigCodes::TrigBitSet& evtTrigBits,
                                      const TrigCodes& trigCodes,
                                      const trigger::TriggerEvent* trigEvt,
                                      const std::string& hltTag) {
    std::vector<TrigCodes::TrigBitSet> partTrigBits(particles.size());
    const double maxDeltaR = 0.1;
    for (const auto& filter : filters) {
      size_t filterNrInEvt = trigEvt->filterIndex(edm::InputTag(filter, "", hltTag));
      const TrigCodes::TrigBitSet filterCode = trigCodes.getCode(filter.c_str());

      if (filterNrInEvt < trigEvt->sizeFilters()) {  //filter found in event, something passes it
        const trigger::Keys& trigKeys = trigEvt->filterKeys(
            filterNrInEvt);  //trigger::Keys is actually a vector<uint16_t> holding the position of trigger objects in the trigger collection passing the filter
        const trigger::TriggerObjectCollection& trigObjColl(trigEvt->getObjects());
        for (size_t partNr = 0; partNr < particles.size(); partNr++) {
          for (unsigned short trigKey : trigKeys) {
            float trigObjEta = trigObjColl[trigKey].eta();
            float trigObjPhi = trigObjColl[trigKey].phi();
            if (reco::deltaR(particles[partNr].eta(), particles[partNr].phi(), trigObjEta, trigObjPhi) < maxDeltaR) {
              partTrigBits[partNr] |= filterCode;
            }  //end dR<maxDeltaR trig obj match test
          }  //end loop over all objects passing filter
        }  //end loop over particles
      }  //end check if filter is present
    }  //end loop over all filters

    //okay the first element is the key, the second is the filter that exists in trigger event
    for (const auto& l1PreAndSeedFilter : l1PreAndSeedFilters) {
      const TrigCodes::TrigBitSet filterCode = trigCodes.getCode(l1PreAndSeedFilter.first.c_str());
      if ((filterCode & evtTrigBits) == filterCode) {  //check that filter has fired in the event

        size_t filterNrInEvt = trigEvt->filterIndex(edm::InputTag(l1PreAndSeedFilter.second, "", hltTag));

        if (filterNrInEvt < trigEvt->sizeFilters()) {  //filter found in event, something passes it
          const trigger::Keys& trigKeys = trigEvt->filterKeys(
              filterNrInEvt);  //trigger::Keys is actually a vector<uint16_t> holding the position of trigger objects in the trigger collection passing the filter
          const trigger::TriggerObjectCollection& trigObjColl(trigEvt->getObjects());
          for (size_t partNr = 0; partNr < particles.size(); partNr++) {
            for (unsigned short trigKey : trigKeys) {
              float trigObjEta = trigObjColl[trigKey].eta();
              float trigObjPhi = trigObjColl[trigKey].phi();
              if (reco::deltaR(particles[partNr].eta(), particles[partNr].phi(), trigObjEta, trigObjPhi) < maxDeltaR) {
                partTrigBits[partNr] |= filterCode;
              }  //end dR<maxDeltaR trig obj match test
            }  //end loop over all objects passing filter
          }  //end loop over particles
        }  //end check if filter is present
      }  //end check if path has fired in the event
    }  //end loop over all filters

    for (size_t partNr = 0; partNr < particles.size(); partNr++)
      particles[partNr].setTrigBits(partTrigBits[partNr]);
  }

  template <class T, class U>
  void trigTools::fillHLTposition(T& particle,
                                  U& hltData,
                                  const std::vector<std::string>& filters,
                                  const trigger::TriggerEvent* trigEvt,
                                  const std::string& hltTag) {
    std::vector<TrigCodes::TrigBitSet> partTrigBits(1);
    const double maxDeltaR = 0.1;
    for (const auto& filter : filters) {
      size_t filterNrInEvt = trigEvt->filterIndex(edm::InputTag(filter, "", hltTag));
      //const TrigCodes::TrigBitSet filterCode = trigCodes.getCode(filters[filterNrInVec].c_str());
      if (filterNrInEvt < trigEvt->sizeFilters()) {  //filter found in event, something passes it
        const trigger::Keys& trigKeys = trigEvt->filterKeys(
            filterNrInEvt);  //trigger::Keys is actually a vector<uint16_t> holding the position of trigger objects in the trigger collection passing the filter
        const trigger::TriggerObjectCollection& trigObjColl(trigEvt->getObjects());
        for (unsigned short trigKey : trigKeys) {
          float trigObjEta = trigObjColl[trigKey].eta();
          float trigObjPhi = trigObjColl[trigKey].phi();
          float trigObjE = trigObjColl[trigKey].energy();
          if (reco::deltaR(particle.superCluster()->eta(), particle.superCluster()->phi(), trigObjEta, trigObjPhi) <
              maxDeltaR) {
            hltData.HLTeta = trigObjEta;
            hltData.HLTphi = trigObjPhi;
            hltData.HLTenergy = trigObjE;
          }  //end dR<maxDeltaR trig obj match test
        }  //end loop over all objects passing filter`
      }  //end check if filter is present
    }  //end check if path has fired in the event
  }  //end loop over all filters

}  // namespace egHLT

#endif