CMS 3D CMS Logo

SiStripBaselineValidator.cc
Go to the documentation of this file.
1 // Original Author: Ivan Amos Cali
2 // Created: Mon Jul 28 14:10:52 CEST 2008
3 //
4 //
5 
6 // system include files
7 #include <cassert>
8 #include <fstream>
9 #include <iostream>
10 #include <memory>
11 
12 // user include files
22 
23 //ROOT inclusion
24 #include "TH1F.h"
25 #include "TH2F.h"
26 #include "TString.h"
27 
28 class TFile;
29 
30 using namespace edm;
31 using namespace std;
32 
34  srcProcessedRawDigi_ = conf.getParameter<edm::InputTag>("srcProcessedRawDigi");
35  moduleRawDigiToken_ =
36  consumes<edm::DetSetVector<SiStripDigi> >(conf.getParameter<edm::InputTag>("srcProcessedRawDigi"));
37 }
38 
40  const edm::Run& run,
41  const edm::EventSetup& es) {
43  ibooker.setCurrentFolder("SiStrip/BaselineValidator");
44 
45  h1NumbadAPVsRes_ = ibooker.book1D("ResAPVs", ";#ResAPVs", 100, 1.0, 10001);
46 
47  h1ADC_vs_strip_ = ibooker.book2D("ADCvsAPVs", ";ADCvsAPVs", 768, -0.5, 767.5, 1023, -0.5, 1022.5);
48 
49  return;
50 }
51 
54  e.getByToken(moduleRawDigiToken_, moduleRawDigi);
55  edm::DetSetVector<SiStripDigi>::const_iterator itRawDigis = moduleRawDigi->begin();
56 
57  int NumResAPVs = 0;
58  for (; itRawDigis != moduleRawDigi->end(); ++itRawDigis) {
59 
60  edm::DetSet<SiStripDigi>::const_iterator itRaw = itRawDigis->begin();
61  int strip = 0, totStripAPV = 0, apv = 0, prevapv = itRaw->strip() / 128;
62 
63  for (; itRaw != itRawDigis->end(); ++itRaw) {
64 
65  strip = itRaw->strip();
66  apv = strip / 128;
67  float adc = itRaw->adc();
68  h1ADC_vs_strip_->Fill(strip, adc);
69 
70  if (prevapv != apv) {
71  if (totStripAPV > 64) {
72  NumResAPVs++;
73  }
74  prevapv = apv;
75  totStripAPV = 0;
76  }
77  if (adc > 0)
78  ++totStripAPV;
79 
80  }
81  if (totStripAPV > 64) {
82  NumResAPVs++;
83  }
84 
85  }
86 
88 
89  h1NumbadAPVsRes_->Fill(NumResAPVs);
90 
91 }
92 
93 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SiStripBaselineValidator(const edm::ParameterSet &)
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:325
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
HLT enums.
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:314
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
Definition: Run.h:45
uint16_t *__restrict__ uint16_t const *__restrict__ adc