CMS 3D CMS Logo

SiPixelPhase1Digis.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1Digis
4 // Class: SiPixelPhase1Digis
5 //
6 
7 // Original Author: Marcel Schneider
8 
9 // C++ stuff
10 #include <iostream>
11 
12 // CMSSW stuff
16 
17 // DQM Stuff
19 
20 // Input data stuff
23 
24 // PixelDQM Framework
26 
27 namespace {
28 
29  class SiPixelPhase1Digis final : public SiPixelPhase1Base {
30  // List of quantities to be plotted.
31  enum {
32  ADC, // digi ADC readouts
33  NDIGIS, // number of digis per event and module
34  NDIGISINCLUSIVE, //Total number of digis in BPix and FPix
35  NDIGIS_FED, // number of digis per event and FED
36  NDIGIS_FEDtrend, // number of digis per event and FED
37  EVENT, // event frequency
38  MAP, // digi hitmap per module
39  OCCUPANCY, // like map but coarser
40 
41  MAX_HIST // a sentinel that gives the number of quantities (not a plot).
42  };
43 
44  public:
45  explicit SiPixelPhase1Digis(const edm::ParameterSet& conf);
46 
47  void analyze(const edm::Event&, const edm::EventSetup&) override;
48 
49  private:
51  };
52 
53  SiPixelPhase1Digis::SiPixelPhase1Digis(const edm::ParameterSet& iConfig) : SiPixelPhase1Base(iConfig) {
54  srcToken_ = consumes<edm::DetSetVector<PixelDigi>>(iConfig.getParameter<edm::InputTag>("src"));
55  }
56 
57  void SiPixelPhase1Digis::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
58  if (!checktrigger(iEvent, iSetup, DCS))
59  return;
60 
62  iEvent.getByToken(srcToken_, input);
63  if (!input.isValid())
64  return;
65  bool hasDigis = false;
66 
68  for (it = input->begin(); it != input->end(); ++it) {
69  for (PixelDigi const& digi : *it) {
70  hasDigis = true;
71  histo[ADC].fill((double)digi.adc(), DetId(it->detId()), &iEvent, digi.column(), digi.row());
72  histo[MAP].fill(DetId(it->detId()), &iEvent, digi.column(), digi.row());
73  histo[OCCUPANCY].fill(DetId(it->detId()), &iEvent, digi.column(), digi.row());
74  histo[NDIGIS].fill(DetId(it->detId()), &iEvent); // count
75  histo[NDIGISINCLUSIVE].fill(DetId(it->detId()), &iEvent); // count
76  histo[NDIGIS_FED].fill(DetId(it->detId()), &iEvent);
77  histo[NDIGIS_FEDtrend].fill(DetId(it->detId()), &iEvent);
78  }
79  }
80  if (hasDigis)
81  histo[EVENT].fill(DetId(0), &iEvent);
82  histo[NDIGIS].executePerEventHarvesting(&iEvent);
83  histo[NDIGISINCLUSIVE].executePerEventHarvesting(&iEvent);
84  histo[NDIGIS_FED].executePerEventHarvesting(&iEvent);
85  histo[NDIGIS_FEDtrend].executePerEventHarvesting(&iEvent);
86  }
87 
88 } //namespace
89 
90 DEFINE_FWK_MODULE(SiPixelPhase1Digis);
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
static std::string const input
Definition: EdmProvDump.cc:47
int iEvent
Definition: GenABIO.cc:224
Definition: DetId.h:17
void analyze(edm::Event const &e, edm::EventSetup const &) override=0
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
struct ADC ADC