CMS 3D CMS Logo

SiPixelPhase1DeadFEDChannels.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1DeadFEDChannels
4 // Class: SiPixelPhase1DeadFEDChannels
5 //
6 
7 // Original Author: F.Fiori
8 
9 // C++ stuff
10 #include <iostream>
11 
12 // CMSSW stuff
22 
23 // DQM Stuff
25 
26 // Input data stuff
30 
31 // PixelDQM Framework
33 
34 namespace {
35 
36  class SiPixelPhase1DeadFEDChannels final : public SiPixelPhase1Base {
37  // List of quantities to be plotted.
38  enum {
39  DEADCHAN,
40  DEADCHANROC
41 
42  };
43 
44  public:
45  explicit SiPixelPhase1DeadFEDChannels(const edm::ParameterSet& conf);
46  void analyze(const edm::Event&, const edm::EventSetup&) override;
47 
48  private:
49  edm::EDGetTokenT<PixelFEDChannelCollection> pixelFEDChannelCollectionToken_;
52 
53  bool firstEvent_;
54  const TrackerGeometry* trackerGeometry_ = nullptr;
55  const SiPixelFedCabling* cablingMap = nullptr;
56  };
57 
58  SiPixelPhase1DeadFEDChannels::SiPixelPhase1DeadFEDChannels(const edm::ParameterSet& iConfig)
59  : SiPixelPhase1Base(iConfig) {
60  pixelFEDChannelCollectionToken_ = consumes<PixelFEDChannelCollection>(edm::InputTag("siPixelDigis"));
61  trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
62  cablingMapToken_ = esConsumes<SiPixelFedCablingMap, SiPixelFedCablingMapRcd>();
63  firstEvent_ = true;
64  };
65 
67  if (!checktrigger(iEvent, iSetup, DCS))
68  return;
69 
70  if (firstEvent_) {
71  edm::ESHandle<TrackerGeometry> tmpTkGeometry = iSetup.getHandle(trackerGeomToken_);
72  trackerGeometry_ = &(*tmpTkGeometry);
73 
74  edm::ESHandle<SiPixelFedCablingMap> pixelCabling = iSetup.getHandle(cablingMapToken_);
75  cablingMap = pixelCabling.product();
76 
77  firstEvent_ = false;
78  }
79 
81 
82  iEvent.getByToken(pixelFEDChannelCollectionToken_, input);
83  if (!input.isValid())
84  return;
85 
86  for (const auto& disabledOnDetId : *input) {
87  for (const auto& ch : disabledOnDetId) {
88  sipixelobjects::CablingPathToDetUnit path = {ch.fed, ch.link, 0};
89 
90  for (path.roc = 1; path.roc <= (ch.roc_last - ch.roc_first) + 1; path.roc++) {
91  const sipixelobjects::PixelROC* roc = cablingMap->findItem(path);
92  assert(roc != nullptr);
93  assert(roc->rawId() == disabledOnDetId.detId());
94 
95  const PixelGeomDetUnit* theGeomDet =
96  dynamic_cast<const PixelGeomDetUnit*>(trackerGeometry_->idToDet(roc->rawId()));
97  PixelTopology const* topology = &(theGeomDet->specificTopology());
99  topology->colsperroc() / 2}; //center of ROC
101  histo[DEADCHANROC].fill(disabledOnDetId.detId(), &iEvent, global.col, global.row);
102  }
103 
104  histo[DEADCHAN].fill(disabledOnDetId.detId(), &iEvent); // global count
105  }
106  }
107 
108  histo[DEADCHAN].executePerEventHarvesting(&iEvent);
109  }
110 } //namespace
111 
112 DEFINE_FWK_MODULE(SiPixelPhase1DeadFEDChannels);
virtual int rowsperroc() const =0
assert(be >=bs)
identify pixel inside single ROC
Definition: LocalPixel.h:7
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
static std::string const input
Definition: EdmProvDump.cc:50
global coordinates (row and column in DetUnit, as in PixelDigi)
Definition: GlobalPixel.h:6
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
virtual int colsperroc() const =0
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
row and collumn in ROC representation
Definition: LocalPixel.h:13
void analyze(edm::Event const &e, edm::EventSetup const &) override=0