HLTDQMFilterEffHists

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 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207
#ifndef DQMOnline_Trigger_HLTDQMHistColl_h
#define DQMOnline_Trigger_HLTDQMHistColl_h

//********************************************************************************
//
// Description:
//   This contains a collection of HLTDQMHists used to measure the efficiency of a
//   specified filter. It is resonsible for booking and filling the histograms
//   For every hist specified, it books two, one to record the
//   total objects passing sample selection passed to the class and then one for those objects
//   which then pass the HLT filter
//   The class contains a simple selection cuts (mostly intended for kinematic range cuts)
//   which are passed to the histograms as they are filled.
//   The cuts are passed to the histograms as some histograms will ignore certain cuts
//   for example, eff vs Et will ignore any global Et cuts applied
//
// Author : Sam Harper , RAL, May 2017
//
//***********************************************************************************

#include "DataFormats/HLTReco/interface/TriggerEvent.h"
#include "DQMOffline/Trigger/interface/HLTDQMHist.h"
#include "DQMOffline/Trigger/interface/VarRangeCutColl.h"
#include "DQMOffline/Trigger/interface/FunctionDefs.h"
#include "DQMOffline/Trigger/interface/UtilFuncs.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

template <typename ObjType>
class HLTDQMFilterEffHists {
public:
  typedef dqm::legacy::DQMStore DQMStore;

  explicit HLTDQMFilterEffHists(const edm::ParameterSet& config, std::string baseHistName, std::string hltProcess);

  static edm::ParameterSetDescription makePSetDescription();
  static edm::ParameterSetDescription makePSetDescriptionHistConfigs();

  void bookHists(DQMStore::IBooker& iBooker, const std::vector<edm::ParameterSet>& histConfigs);
  void fillHists(const ObjType& obj,
                 const edm::Event& event,
                 const edm::EventSetup& setup,
                 const trigger::TriggerEvent& trigEvt);

private:
  void book1D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig);
  void book2D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig);

private:
  std::vector<std::unique_ptr<HLTDQMHist<ObjType>>> histsPass_;
  std::vector<std::unique_ptr<HLTDQMHist<ObjType>>> histsTot_;
  VarRangeCutColl<ObjType> rangeCuts_;
  std::string filterName_;
  std::string histTitle_;
  std::string folderName_;
  std::string baseHistName_;
  std::string hltProcess_;
};

template <typename ObjType>
HLTDQMFilterEffHists<ObjType>::HLTDQMFilterEffHists(const edm::ParameterSet& config,
                                                    std::string baseHistName,
                                                    std::string hltProcess)
    : rangeCuts_(config.getParameter<std::vector<edm::ParameterSet>>("rangeCuts")),
      filterName_(config.getParameter<std::string>("filterName")),
      histTitle_(config.getParameter<std::string>("histTitle")),
      folderName_(config.getParameter<std::string>("folderName")),
      baseHistName_(std::move(baseHistName)),
      hltProcess_(std::move(hltProcess)) {}

template <typename ObjType>
edm::ParameterSetDescription HLTDQMFilterEffHists<ObjType>::makePSetDescription() {
  edm::ParameterSetDescription desc;
  desc.addVPSet("rangeCuts", VarRangeCut<ObjType>::makePSetDescription(), std::vector<edm::ParameterSet>());
  desc.add<std::string>("filterName", "");
  desc.add<std::string>("histTitle", "");
  desc.add<std::string>("folderName", "");
  return desc;
}

template <typename ObjType>
edm::ParameterSetDescription HLTDQMFilterEffHists<ObjType>::makePSetDescriptionHistConfigs() {
  edm::ParameterSetDescription desc;

  //what this is doing is trival and is left as an exercise to the reader
  auto histDescCases =
      "1D" >> (edm::ParameterDescription<std::vector<double>>("binLowEdges", std::vector<double>(), true) and
               edm::ParameterDescription<std::string>("nameSuffex", "", true) and
               edm::ParameterDescription<std::string>("vsVar", "", true)) or
      "2D" >> (edm::ParameterDescription<std::vector<double>>("xBinLowEdges", std::vector<double>(), true) and
               edm::ParameterDescription<std::vector<double>>("yBinLowEdges", std::vector<double>(), true) and
               edm::ParameterDescription<std::string>("nameSuffex", "", true) and
               edm::ParameterDescription<std::string>("xVar", "", true) and
               edm::ParameterDescription<std::string>("yVar", "", true));

  desc.ifValue(edm::ParameterDescription<std::string>("histType", "1D", true), std::move(histDescCases));
  desc.addVPSet("rangeCuts", VarRangeCut<ObjType>::makePSetDescription(), std::vector<edm::ParameterSet>());
  return desc;
}

template <typename ObjType>
void HLTDQMFilterEffHists<ObjType>::bookHists(DQMStore::IBooker& iBooker,
                                              const std::vector<edm::ParameterSet>& histConfigs) {
  iBooker.setCurrentFolder(folderName_);
  for (auto& histConfig : histConfigs) {
    std::string histType = histConfig.getParameter<std::string>("histType");
    if (histType == "1D") {
      book1D(iBooker, histConfig);
    } else if (histType == "2D") {
      book2D(iBooker, histConfig);
    } else {
      throw cms::Exception("ConfigError") << " histType " << histType << " not recognised" << std::endl;
    }
  }
}

template <typename ObjType>
void HLTDQMFilterEffHists<ObjType>::book1D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig) {
  auto binLowEdgesDouble = histConfig.getParameter<std::vector<double>>("binLowEdges");
  std::vector<float> binLowEdges;
  binLowEdges.reserve(binLowEdgesDouble.size());
  for (double lowEdge : binLowEdgesDouble)
    binLowEdges.push_back(lowEdge);
  auto nameSuffex = histConfig.getParameter<std::string>("nameSuffex");
  auto mePass = iBooker.book1D((baseHistName_ + filterName_ + nameSuffex + "_pass").c_str(),
                               (histTitle_ + nameSuffex + " Pass").c_str(),
                               binLowEdges.size() - 1,
                               &binLowEdges[0]);
  std::unique_ptr<HLTDQMHist<ObjType>> hist;
  auto vsVar = histConfig.getParameter<std::string>("vsVar");
  auto vsVarFunc = hltdqm::getUnaryFuncFloat<ObjType>(vsVar);
  if (!vsVarFunc) {
    throw cms::Exception("ConfigError") << " vsVar " << vsVar << " is giving null ptr (likely empty) in " << __FILE__
                                        << "," << __LINE__ << std::endl;
  }
  VarRangeCutColl<ObjType> rangeCuts(histConfig.getParameter<std::vector<edm::ParameterSet>>("rangeCuts"));
  hist = std::make_unique<HLTDQMHist1D<ObjType, float>>(mePass, vsVar, vsVarFunc, rangeCuts);
  histsPass_.emplace_back(std::move(hist));
  auto meTot = iBooker.book1D((baseHistName_ + filterName_ + nameSuffex + "_tot").c_str(),
                              (histTitle_ + nameSuffex + " Total").c_str(),
                              binLowEdges.size() - 1,
                              &binLowEdges[0]);
  hist = std::make_unique<HLTDQMHist1D<ObjType, float>>(meTot, vsVar, vsVarFunc, rangeCuts);
  histsTot_.emplace_back(std::move(hist));
}

template <typename ObjType>
void HLTDQMFilterEffHists<ObjType>::book2D(DQMStore::IBooker& iBooker, const edm::ParameterSet& histConfig) {
  auto xBinLowEdgesDouble = histConfig.getParameter<std::vector<double>>("xBinLowEdges");
  auto yBinLowEdgesDouble = histConfig.getParameter<std::vector<double>>("yBinLowEdges");
  std::vector<float> xBinLowEdges;
  std::vector<float> yBinLowEdges;
  xBinLowEdges.reserve(xBinLowEdgesDouble.size());
  for (double lowEdge : xBinLowEdgesDouble)
    xBinLowEdges.push_back(lowEdge);
  yBinLowEdges.reserve(yBinLowEdgesDouble.size());
  for (double lowEdge : yBinLowEdgesDouble)
    yBinLowEdges.push_back(lowEdge);
  auto nameSuffex = histConfig.getParameter<std::string>("nameSuffex");
  auto mePass = iBooker.book2D((baseHistName_ + filterName_ + nameSuffex + "_pass").c_str(),
                               (histTitle_ + nameSuffex + " Pass").c_str(),
                               xBinLowEdges.size() - 1,
                               &xBinLowEdges[0],
                               yBinLowEdges.size() - 1,
                               &yBinLowEdges[0]);
  std::unique_ptr<HLTDQMHist<ObjType>> hist;
  auto xVar = histConfig.getParameter<std::string>("xVar");
  auto yVar = histConfig.getParameter<std::string>("yVar");
  auto xVarFunc = hltdqm::getUnaryFuncFloat<ObjType>(xVar);
  auto yVarFunc = hltdqm::getUnaryFuncFloat<ObjType>(yVar);
  if (!xVarFunc || !yVarFunc) {
    throw cms::Exception("ConfigError") << " xVar " << xVar << " or yVar " << yVar
                                        << " is giving null ptr (likely empty str passed)" << std::endl;
  }
  VarRangeCutColl<ObjType> rangeCuts(histConfig.getParameter<std::vector<edm::ParameterSet>>("rangeCuts"));

  hist = std::make_unique<HLTDQMHist2D<ObjType, float>>(mePass, xVar, yVar, xVarFunc, yVarFunc, rangeCuts);
  histsPass_.emplace_back(std::move(hist));

  auto meTot = iBooker.book2D((baseHistName_ + filterName_ + nameSuffex + "_tot").c_str(),
                              (histTitle_ + nameSuffex + " Total").c_str(),
                              xBinLowEdges.size() - 1,
                              &xBinLowEdges[0],
                              yBinLowEdges.size() - 1,
                              &yBinLowEdges[0]);

  hist = std::make_unique<HLTDQMHist2D<ObjType, float>>(meTot, xVar, yVar, xVarFunc, yVarFunc, rangeCuts);
  histsTot_.emplace_back(std::move(hist));
}

template <typename ObjType>
void HLTDQMFilterEffHists<ObjType>::fillHists(const ObjType& obj,
                                              const edm::Event& event,
                                              const edm::EventSetup& setup,
                                              const trigger::TriggerEvent& trigEvt) {
  for (auto& hist : histsTot_) {
    hist->fill(obj, event, setup, rangeCuts_);
  }

  if (hltdqm::passTrig(obj.eta(), obj.phi(), trigEvt, filterName_, hltProcess_)) {
    for (auto& hist : histsPass_) {
      hist->fill(obj, event, setup, rangeCuts_);
    }
  }
}
#endif