Histograms

LumiMonitor

MEbinning

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 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309
#include <string>

#include "DQM/TrackingMonitor/interface/GetLumi.h"
#include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"
#include "DataFormats/OnlineMetaData/interface/OnlineLuminosityRecord.h"
#include "DataFormats/Scalers/interface/LumiScalers.h"
#include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
#include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "FWCore/Utilities/interface/EDGetToken.h"
#include "Geometry/Records/interface/TrackerTopologyRcd.h"

namespace {
  typedef dqm::reco::DQMStore DQMStore;
  typedef dqm::reco::MonitorElement MonitorElement;
  struct MEbinning {
    int nbins;
    double xmin;
    double xmax;
  };

  struct Histograms {
    dqm::reco::MonitorElement* numberOfPixelClustersVsLS;
    dqm::reco::MonitorElement* numberOfPixelClustersVsLumi;
    dqm::reco::MonitorElement* lumiVsLS;
    dqm::reco::MonitorElement* puVsLS;
    dqm::reco::MonitorElement* pixelLumiVsLS;
    dqm::reco::MonitorElement* pixelLumiVsLumi;
  };
}  // namespace

//
// class declaration
//

class LumiMonitor : public DQMGlobalEDAnalyzer<Histograms> {
public:
  LumiMonitor(const edm::ParameterSet&);
  ~LumiMonitor() override = default;

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

  static void fillHistoPSetDescription(edm::ParameterSetDescription& pset,
                                       int const nbins,
                                       double const xmin,
                                       double const xmax);
  static void fillHistoLSPSetDescription(edm::ParameterSetDescription& pset, int const nbins);

private:
  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, Histograms&) const override;
  void dqmAnalyze(edm::Event const&, edm::EventSetup const&, Histograms const&) const override;

  static MEbinning getHistoPSet(const edm::ParameterSet& pset);
  static MEbinning getHistoLSPSet(const edm::ParameterSet& pset);

  std::string const folderName_;

  edm::EDGetTokenT<LumiScalersCollection> const lumiScalersToken_;
  edm::EDGetTokenT<OnlineLuminosityRecord> const onlineMetaDataDigisToken_;
  MEbinning const lumi_binning_;
  MEbinning const pu_binning_;
  MEbinning const ls_binning_;

  bool const doPixelLumi_;
  edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> const trkTopoToken_;
  edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> const pixelClustersToken_;
  bool const useBPixLayer1_;
  int const minNumberOfPixelsPerCluster_;
  float const minPixelClusterCharge_;
  MEbinning const pixelCluster_binning_;
  MEbinning const pixellumi_binning_;
  float const lumi_factor_per_bx_;
};

// -----------------------------
//  constructors and destructor
// -----------------------------

LumiMonitor::LumiMonitor(const edm::ParameterSet& config)
    : folderName_(config.getParameter<std::string>("folderName")),
      lumiScalersToken_(consumes(config.getParameter<edm::InputTag>("scalers"))),
      onlineMetaDataDigisToken_(consumes(config.getParameter<edm::InputTag>("onlineMetaDataDigis"))),
      lumi_binning_(getHistoPSet(
          config.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("lumiPSet"))),
      pu_binning_(
          getHistoPSet(config.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("puPSet"))),
      ls_binning_(getHistoLSPSet(
          config.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("lsPSet"))),
      doPixelLumi_(config.getParameter<bool>("doPixelLumi")),
      trkTopoToken_(doPixelLumi_ ? esConsumes<TrackerTopology, TrackerTopologyRcd>()
                                 : edm::ESGetToken<TrackerTopology, TrackerTopologyRcd>()),
      pixelClustersToken_(doPixelLumi_ ? consumes<edmNew::DetSetVector<SiPixelCluster>>(
                                             config.getParameter<edm::InputTag>("pixelClusters"))
                                       : edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>>()),
      useBPixLayer1_(doPixelLumi_ ? config.getParameter<bool>("useBPixLayer1") : false),
      minNumberOfPixelsPerCluster_(doPixelLumi_ ? config.getParameter<int>("minNumberOfPixelsPerCluster") : -1),
      minPixelClusterCharge_(doPixelLumi_ ? config.getParameter<double>("minPixelClusterCharge") : -1.),
      pixelCluster_binning_(doPixelLumi_ ? getHistoPSet(config.getParameter<edm::ParameterSet>("histoPSet")
                                                            .getParameter<edm::ParameterSet>("pixelClusterPSet"))
                                         : MEbinning{}),
      pixellumi_binning_(doPixelLumi_ ? getHistoPSet(config.getParameter<edm::ParameterSet>("histoPSet")
                                                         .getParameter<edm::ParameterSet>("pixellumiPSet"))
                                      : MEbinning{}),
      lumi_factor_per_bx_(useBPixLayer1_
                              ? GetLumi::FREQ_ORBIT * GetLumi::SECONDS_PER_LS / GetLumi::XSEC_PIXEL_CLUSTER
                              : GetLumi::FREQ_ORBIT * GetLumi::SECONDS_PER_LS / GetLumi::rXSEC_PIXEL_CLUSTER) {}

MEbinning LumiMonitor::getHistoPSet(const edm::ParameterSet& pset) {
  return MEbinning{
      pset.getParameter<int32_t>("nbins"),
      pset.getParameter<double>("xmin"),
      pset.getParameter<double>("xmax"),
  };
}

MEbinning LumiMonitor::getHistoLSPSet(const edm::ParameterSet& pset) {
  return MEbinning{pset.getParameter<int32_t>("nbins"), -0.5, pset.getParameter<int32_t>("nbins") - 0.5};
}

void LumiMonitor::bookHistograms(DQMStore::IBooker& booker,
                                 edm::Run const& run,
                                 edm::EventSetup const& setup,
                                 Histograms& histograms) const {
  booker.setCurrentFolder(folderName_);

  if (doPixelLumi_) {
    auto me = booker.book1D("numberOfPixelClustersVsLS",
                            "number of pixel clusters vs lumisection",
                            ls_binning_.nbins,
                            ls_binning_.xmin,
                            ls_binning_.xmax);
    me->setAxisTitle("lumisection", 1);
    me->setAxisTitle("number of pixel clusters", 2);
    histograms.numberOfPixelClustersVsLS = me;

    me = booker.bookProfile("numberOfPixelClustersVsLumi",
                            "number of pixel clusters vs online lumi",
                            lumi_binning_.nbins,
                            lumi_binning_.xmin,
                            lumi_binning_.xmax,
                            pixelCluster_binning_.xmin,
                            pixelCluster_binning_.xmax);
    me->setAxisTitle("online inst lumi E30 [Hz cm^{-2}]", 1);
    me->setAxisTitle("number of pixel clusters", 2);
    histograms.numberOfPixelClustersVsLumi = me;

    me = booker.bookProfile("pixelLumiVsLS",
                            "pixel-lumi vs lumisection",
                            ls_binning_.nbins,
                            ls_binning_.xmin,
                            ls_binning_.xmax,
                            pixellumi_binning_.xmin,
                            pixellumi_binning_.xmax);
    me->setAxisTitle("lumisection", 1);
    me->setAxisTitle("pixel-based inst lumi E30 [Hz cm^{-2}]", 2);
    histograms.pixelLumiVsLS = me;

    me = booker.bookProfile("pixelLumiVsLumi",
                            "pixel-lumi vs online lumi",
                            lumi_binning_.nbins,
                            lumi_binning_.xmin,
                            lumi_binning_.xmax,
                            pixellumi_binning_.xmin,
                            lumi_binning_.xmax);
    me->setAxisTitle("online inst lumi E30 [Hz cm^{-2}]", 1);
    me->setAxisTitle("pixel-based inst lumi E30 [Hz cm^{-2}]", 2);
    histograms.pixelLumiVsLumi = me;
  }

  auto me = booker.bookProfile("lumiVsLS",
                               "online lumi vs lumisection",
                               ls_binning_.nbins,
                               ls_binning_.xmin,
                               ls_binning_.xmax,
                               lumi_binning_.xmin,
                               lumi_binning_.xmax);
  me->setAxisTitle("lumisection", 1);
  me->setAxisTitle("online inst lumi E30 [Hz cm^{-2}]", 2);
  histograms.lumiVsLS = me;

  me = booker.bookProfile("puVsLS",
                          "online pileup vs lumisection",
                          ls_binning_.nbins,
                          ls_binning_.xmin,
                          ls_binning_.xmax,
                          pu_binning_.xmin,
                          pu_binning_.xmax);
  me->setAxisTitle("lumisection", 1);
  me->setAxisTitle("online pileup", 2);
  histograms.puVsLS = me;
}

void LumiMonitor::dqmAnalyze(edm::Event const& event,
                             edm::EventSetup const& setup,
                             Histograms const& histograms) const {
  int const ls = event.id().luminosityBlock();

  float online_lumi = -1.f;
  float online_pu = -1.f;
  auto const lumiScalersHandle = event.getHandle(lumiScalersToken_);
  auto const onlineMetaDataDigisHandle = event.getHandle(onlineMetaDataDigisToken_);
  if (lumiScalersHandle.isValid() and not lumiScalersHandle->empty()) {
    auto const scalit = lumiScalersHandle->begin();
    online_lumi = scalit->instantLumi();
    online_pu = scalit->pileup();
  } else if (onlineMetaDataDigisHandle.isValid()) {
    online_lumi = onlineMetaDataDigisHandle->instLumi();
    online_pu = onlineMetaDataDigisHandle->avgPileUp();
  }
  histograms.lumiVsLS->Fill(ls, online_lumi);
  histograms.puVsLS->Fill(ls, online_pu);

  if (doPixelLumi_) {
    size_t pixel_clusters = 0;
    float pixel_lumi = -1.f;
    edm::Handle<edmNew::DetSetVector<SiPixelCluster>> pixelClusters;
    event.getByToken(pixelClustersToken_, pixelClusters);
    if (pixelClusters.isValid()) {
      auto const& tTopo = setup.getData(trkTopoToken_);

      // Count the number of clusters with at least a minimum
      // number of pixels per cluster and at least a minimum charge.
      for (auto pixCluDet = pixelClusters->begin(); pixCluDet != pixelClusters->end(); ++pixCluDet) {
        DetId detid = pixCluDet->detId();
        size_t subdetid = detid.subdetId();
        if (subdetid == (int)PixelSubdetector::PixelBarrel) {
          if (tTopo.layer(detid) == 1) {
            continue;
          }
        }

        for (auto pixClu = pixCluDet->begin(); pixClu != pixCluDet->end(); ++pixClu) {
          if ((pixClu->size() >= minNumberOfPixelsPerCluster_) and (pixClu->charge() >= minPixelClusterCharge_)) {
            ++pixel_clusters;
          }
        }
      }
      pixel_lumi = lumi_factor_per_bx_ * pixel_clusters / GetLumi::CM2_TO_NANOBARN;  // ?!?!
    } else {
      pixel_lumi = -1.;
    }

    histograms.numberOfPixelClustersVsLS->Fill(ls, pixel_clusters);
    histograms.numberOfPixelClustersVsLumi->Fill(online_lumi, pixel_clusters);
    histograms.pixelLumiVsLS->Fill(ls, pixel_lumi);
    histograms.pixelLumiVsLumi->Fill(online_lumi, pixel_lumi);
  }
}

void LumiMonitor::fillHistoPSetDescription(edm::ParameterSetDescription& pset,
                                           int const nbins,
                                           double const xmin,
                                           double const xmax) {
  pset.add<int>("nbins", nbins);
  pset.add<double>("xmin", xmin);
  pset.add<double>("xmax", xmax);
}

void LumiMonitor::fillHistoLSPSetDescription(edm::ParameterSetDescription& pset, int const nbins) {
  pset.add<int>("nbins", nbins);
}

void LumiMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("pixelClusters", edm::InputTag("hltSiPixelClusters"));
  desc.add<edm::InputTag>("scalers", edm::InputTag("hltScalersRawToDigi"));
  desc.add<edm::InputTag>("onlineMetaDataDigis", edm::InputTag("hltOnlineMetaDataDigis"));
  desc.add<std::string>("folderName", "HLT/LumiMonitoring");
  desc.add<bool>("doPixelLumi", false);
  desc.add<bool>("useBPixLayer1", false);
  desc.add<int>("minNumberOfPixelsPerCluster", 2);  // from DQM/PixelLumi/python/PixelLumiDQM_cfi.py
  desc.add<double>("minPixelClusterCharge", 15000.);

  edm::ParameterSetDescription histoPSet;

  edm::ParameterSetDescription lsPSet;
  fillHistoLSPSetDescription(lsPSet, 2500);
  histoPSet.add<edm::ParameterSetDescription>("lsPSet", lsPSet);

  edm::ParameterSetDescription puPSet;
  fillHistoPSetDescription(puPSet, 130, 0, 130);
  histoPSet.add<edm::ParameterSetDescription>("puPSet", puPSet);

  edm::ParameterSetDescription lumiPSet;
  fillHistoPSetDescription(lumiPSet, 5000, 0, 20000);
  histoPSet.add<edm::ParameterSetDescription>("lumiPSet", lumiPSet);

  edm::ParameterSetDescription pixellumiPSet;
  fillHistoPSetDescription(pixellumiPSet, 300, 0, 3);
  histoPSet.add<edm::ParameterSetDescription>("pixellumiPSet", pixellumiPSet);

  edm::ParameterSetDescription pixelClusterPSet;
  fillHistoPSetDescription(pixelClusterPSet, 200, -0.5, 19999.5);
  histoPSet.add("pixelClusterPSet", pixelClusterPSet);

  desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);

  descriptions.add("lumiMonitor", desc);
}

// Define this as a plug-in
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(LumiMonitor);