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
// -*- C++ -*-
//
// Package:    SiStripMonitorCluster
// Class:      SiStripMonitorHLT
//
// class SiStripMonitorHLT SiStripMonitorHLT.cc
// DQM/SiStripMonitorCluster/src/SiStripMonitorHLT.cc

#include <vector>
#include <iostream>
#include <numeric>

#include "DQM/SiStripMonitorCluster/interface/SiStripMonitorHLT.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

SiStripMonitorHLT::SiStripMonitorHLT(const edm::ParameterSet& iConfig) {
  HLTDirectory = "HLTResults";
  conf_ = iConfig;

  filerDecisionToken_ = consumes<int>(conf_.getParameter<std::string>("HLTProducer"));
  sumOfClusterToken_ = consumes<uint>(conf_.getParameter<std::string>("HLTProducer"));
  clusterInSubComponentsToken_ =
      consumes<std::map<uint, std::vector<SiStripCluster> > >(conf_.getParameter<std::string>("HLTProducer"));
}

void SiStripMonitorHLT::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run& run, const edm::EventSetup& es) {
  ibooker.setCurrentFolder(HLTDirectory);
  std::string HLTProducer = conf_.getParameter<std::string>("HLTProducer");
  HLTDecision = ibooker.book1D(HLTProducer + "_HLTDecision", HLTProducer + "HLTDecision", 2, -0.5, 1.5);
  // all
  SumOfClusterCharges_all = ibooker.book1D("SumOfClusterCharges_all", "SumOfClusterCharges_all", 50, 0, 2000);
  ChargeOfEachClusterTIB_all =
      ibooker.book1D("ChargeOfEachClusterTIB_all", "ChargeOfEachClusterTIB_all", 400, -0.5, 400.5);
  ChargeOfEachClusterTOB_all =
      ibooker.book1D("ChargeOfEachClusterTOB_all", "ChargeOfEachClusterTOB_all", 400, -0.5, 400.5);
  ChargeOfEachClusterTEC_all =
      ibooker.book1D("ChargeOfEachClusterTEC_all", "ChargeOfEachClusterTEC_all", 400, -0.5, 400.5);
  NumberOfClustersAboveThreshold_all =
      ibooker.book1D("NumberOfClustersAboveThreshold_all", "NumberOfClustersAboveThreshold_all", 30, 30.5, 60.5);
  // 31 = TIB2, 32 = TIB2, 33 = TIB3, 51 = TOB1, 52=TOB2, 60 = TEC
  // accepted from HLT
  SumOfClusterCharges_hlt = ibooker.book1D("SumOfClusterCharges_hlt", "SumOfClusterCharges_hlt", 50, 0, 2000);
  ChargeOfEachClusterTIB_hlt =
      ibooker.book1D("ChargeOfEachClusterTIB_hlt", "ChargeOfEachClusterTIB_hlt", 400, -0.5, 400.5);
  ChargeOfEachClusterTOB_hlt =
      ibooker.book1D("ChargeOfEachClusterTOB_hlt", "ChargeOfEachClusterTOB_hlt", 400, -0.5, 400.5);
  ChargeOfEachClusterTEC_hlt =
      ibooker.book1D("ChargeOfEachClusterTEC_hlt", "ChargeOfEachClusterTEC_hlt", 400, -0.5, 400.5);
  NumberOfClustersAboveThreshold_hlt =
      ibooker.book1D("NumberOfClustersAboveThreshold_hlt", "NumberOfClustersAboveThreshold_hlt", 30, 30.5, 60.5);
}

void SiStripMonitorHLT::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
  // get from event
  std::string HLTProducer = conf_.getParameter<std::string>("HLTProducer");
  edm::Handle<int> filter_decision;
  iEvent.getByToken(filerDecisionToken_, filter_decision);  // filter decision
  edm::Handle<uint> sum_of_clustch;
  iEvent.getByToken(sumOfClusterToken_,
                    sum_of_clustch);  // sum of cluster charges
  // first element of pair: layer: TIB1, ...., TEC; second element: nr of
  // clusters above threshold
  edm::Handle<std::map<uint, std::vector<SiStripCluster> > > clusters_in_subcomponents;
  if (HLTProducer == "ClusterMTCCFilter")
    iEvent.getByToken(clusterInSubComponentsToken_, clusters_in_subcomponents);

  // trigger decision
  HLTDecision->Fill(*filter_decision);

  // sum of charges of clusters
  SumOfClusterCharges_all->Fill(*sum_of_clustch);
  if (*filter_decision)
    SumOfClusterCharges_hlt->Fill(*sum_of_clustch);

  // clusters in different layers
  if (HLTProducer == "ClusterMTCCFilter") {
    // loop over layers ("subcomponents")
    for (std::map<uint, std::vector<SiStripCluster> >::const_iterator it = clusters_in_subcomponents->begin();
         it != clusters_in_subcomponents->end();
         it++) {
      int generalized_layer = it->first;
      std::vector<SiStripCluster> theclusters = it->second;
      NumberOfClustersAboveThreshold_all->Fill(generalized_layer,
                                               theclusters.size());  // number of clusters in this generalized layer
      if (*filter_decision)
        NumberOfClustersAboveThreshold_hlt->Fill(generalized_layer, theclusters.size());
      // loop over clusters (and detids)
      for (std::vector<SiStripCluster>::const_iterator icluster = theclusters.begin(); icluster != theclusters.end();
           icluster++) {
        // calculate sum of amplitudes
        unsigned int amplclus = 0;
        for (auto ia = icluster->amplitudes().begin(); ia != icluster->amplitudes().end(); ia++) {
          if ((*ia) > 0)
            amplclus += (*ia);  // why should this be negative?
        }
        if (generalized_layer == 31 || generalized_layer == 32 ||
            generalized_layer == 33) {  // you can also ask the detid here whether is TIB
          ChargeOfEachClusterTIB_all->Fill(amplclus, 1.);
          if (*filter_decision)
            ChargeOfEachClusterTIB_hlt->Fill(amplclus, 1.);
        }
        if (generalized_layer == 51 || generalized_layer == 52) {
          ChargeOfEachClusterTOB_all->Fill(amplclus, 1.);
          if (*filter_decision)
            ChargeOfEachClusterTOB_hlt->Fill(amplclus, 1.);
        }
        if (generalized_layer == 60) {
          ChargeOfEachClusterTEC_all->Fill(amplclus, 1.);
          if (*filter_decision)
            ChargeOfEachClusterTEC_hlt->Fill(amplclus, 1.);
        }
      }
    }
  }
}