PixelMEs

PixelVTXMonitor

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
/*
 * \file PixelVTXMonitor.cc
 * \author S. Dutta
 * Last Update:
 *
 * Description: Pixel Vertex Monitoring for different HLT paths
 *
*/

// system includes
#include <string>
#include <vector>
#include <map>

// user includes
#include "DQMServices/Core/interface/DQMEDAnalyzer.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/TriggerResults.h"
#include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
#include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/LuminosityBlock.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"

// ROOT includes
#include "TPRegexp.h"

//
// class declaration
//

class PixelVTXMonitor : public DQMEDAnalyzer {
public:
  typedef dqm::legacy::MonitorElement MonitorElement;
  typedef dqm::legacy::DQMStore DQMStore;
  PixelVTXMonitor(const edm::ParameterSet&);
  ~PixelVTXMonitor() override = default;

protected:
  void bookHistograms(DQMStore::IBooker& iBooker, const edm::Run& iRun, const edm::EventSetup& iSetup) override;

private:
  void dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) override;
  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;

  edm::ParameterSet parameters_;

  const std::string moduleName_;
  const std::string folderName_;
  const edm::InputTag pixelClusterInputTag_;
  const edm::InputTag pixelVertexInputTag_;
  const edm::InputTag hltInputTag_;
  const edm::EDGetTokenT<SiPixelClusterCollectionNew> pixelClusterInputTagToken_;
  const edm::EDGetTokenT<reco::VertexCollection> pixelVertexInputTagToken_;
  const edm::EDGetTokenT<edm::TriggerResults> hltInputTagToken_;
  const float minVtxDoF_;

  HLTConfigProvider hltConfig_;

  struct PixelMEs {
    MonitorElement* clusME;
    MonitorElement* vtxME;
  };

  std::map<std::string, PixelMEs> histoMap_;
};

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

PixelVTXMonitor::PixelVTXMonitor(const edm::ParameterSet& ps)
    : parameters_(ps),
      moduleName_(parameters_.getParameter<std::string>("ModuleName")),
      folderName_(parameters_.getParameter<std::string>("FolderName")),
      pixelClusterInputTag_(parameters_.getUntrackedParameter<edm::InputTag>("PixelClusterInputTag")),
      pixelVertexInputTag_(parameters_.getUntrackedParameter<edm::InputTag>("PixelVertexInputTag")),
      hltInputTag_(parameters_.getUntrackedParameter<edm::InputTag>("HLTInputTag")),
      pixelClusterInputTagToken_(consumes<SiPixelClusterCollectionNew>(pixelClusterInputTag_)),
      pixelVertexInputTagToken_(consumes<reco::VertexCollection>(pixelVertexInputTag_)),
      hltInputTagToken_(consumes<edm::TriggerResults>(hltInputTag_)),
      minVtxDoF_(parameters_.getParameter<double>("MinVtxDoF")) {}

void PixelVTXMonitor::bookHistograms(DQMStore::IBooker& iBooker, const edm::Run&, const edm::EventSetup&) {
  std::vector<std::string> hltPathsOfInterest =
      parameters_.getParameter<std::vector<std::string> >("HLTPathsOfInterest");
  if (hltPathsOfInterest.empty())
    return;

  const std::vector<std::string>& pathList = hltConfig_.triggerNames();
  std::vector<std::string> selectedPaths;
  for (const auto& it : pathList) {
    int nmatch = 0;
    for (const auto& kt : hltPathsOfInterest) {
      nmatch += TPRegexp(kt).Match(it);
    }
    if (!nmatch)
      continue;
    else
      selectedPaths.push_back(it);
  }

  edm::ParameterSet ClusHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1ClusPar");
  edm::ParameterSet VtxHistoPar = parameters_.getParameter<edm::ParameterSet>("TH1VtxPar");

  std::string currentFolder = moduleName_ + "/" + folderName_;
  iBooker.setCurrentFolder(currentFolder);

  PixelMEs local_MEs;
  for (const auto& tag : selectedPaths) {
    std::map<std::string, PixelMEs>::iterator iPos = histoMap_.find(tag);
    if (iPos == histoMap_.end()) {
      std::string hname, htitle;

      hname = "nPxlClus_";
      hname += tag;
      htitle = "# of Pixel Clusters (";
      htitle += tag + ")";
      local_MEs.clusME = iBooker.book1D(hname,
                                        htitle,
                                        ClusHistoPar.getParameter<int32_t>("Xbins"),
                                        ClusHistoPar.getParameter<double>("Xmin"),
                                        ClusHistoPar.getParameter<double>("Xmax"));

      hname = "nPxlVtx_";
      hname += tag;
      htitle = "# of Pixel Vertices (";
      htitle += tag + ")";
      local_MEs.vtxME = iBooker.book1D(hname,
                                       htitle,
                                       VtxHistoPar.getParameter<int32_t>("Xbins"),
                                       VtxHistoPar.getParameter<double>("Xmin"),
                                       VtxHistoPar.getParameter<double>("Xmax"));

      histoMap_.insert(std::make_pair(tag, local_MEs));
    }
  }
}

void PixelVTXMonitor::dqmBeginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
  bool changed = true;
  if (hltConfig_.init(iRun, iSetup, hltInputTag_.process(), changed)) {
    // if init returns TRUE, initialisation has succeeded!
    edm::LogInfo("PixelVTXMonitor") << "HLT config with process name " << hltInputTag_.process()
                                    << " successfully extracted";
  } else {
    // if init returns FALSE, initialisation has NOT succeeded, which indicates a problem
    // with the file and/or code and needs to be investigated!
    edm::LogError("PixelVTXMonotor") << "Error! HLT config extraction with process name " << hltInputTag_.process()
                                     << " failed";
    // In this case, all access methods will return empty values!
  }
}
void PixelVTXMonitor::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
  if (histoMap_.empty())
    return;

  //Access Pixel Clusters
  edm::Handle<SiPixelClusterCollectionNew> siPixelClusters = iEvent.getHandle(pixelClusterInputTagToken_);

  if (!siPixelClusters.isValid()) {
    edm::LogError("PixelVTXMonotor") << "Could not find Cluster Collection " << pixelClusterInputTag_;
    return;
  }
  unsigned nClusters = siPixelClusters->size();

  //Access Pixel Verteces
  edm::Handle<reco::VertexCollection> pixelVertices = iEvent.getHandle(pixelVertexInputTagToken_);
  if (!pixelVertices.isValid()) {
    edm::LogError("PixelVTXMonotor") << "Could not find Vertex Collection " << pixelVertexInputTag_;
    return;
  }

  int nVtx = 0;
  for (const auto& ivtx : *pixelVertices) {
    if (minVtxDoF_ == -1)
      nVtx++;
    else {
      if ((ivtx.isValid() == true) && (ivtx.isFake() == false) && (ivtx.ndof() >= minVtxDoF_) &&
          (ivtx.tracksSize() != 0))
        nVtx++;
    }
  }
  // Access Trigger Results
  edm::Handle<edm::TriggerResults> triggerResults = iEvent.getHandle(hltInputTagToken_);
  if (!triggerResults.isValid())
    return;

  for (const auto& it : histoMap_) {
    std::string path = it.first;
    MonitorElement* me_clus = it.second.clusME;
    MonitorElement* me_vtx = it.second.vtxME;
    unsigned int index = hltConfig_.triggerIndex(path);
    if (index < triggerResults->size() && triggerResults->accept(index)) {
      if (me_vtx)
        me_vtx->Fill(nVtx);
      if (me_clus)
        me_clus->Fill(nClusters);
    }
  }
}

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