CMS 3D CMS Logo

CSCStripDigiValidation.cc
Go to the documentation of this file.
5 
7  : CSCBaseValidation(inputTag),
8  thePedestalSum(0),
9  thePedestalCovarianceSum(0),
10  thePedestalCount(0),
11  thePedestalTimeCorrelationPlot(nullptr),
12  thePedestalNeighborCorrelationPlot(nullptr),
13  theNDigisPerChamberPlot(nullptr) {
15 }
16 
18 
20  thePedestalPlot = iBooker.book1D("CSCPedestal", "CSC Pedestal ", 400, 550, 650);
21  theAmplitudePlot = iBooker.book1D("CSCStripAmplitude", "CSC Strip Amplitude", 200, 0, 2000);
22  theRatio4to5Plot = iBooker.book1D("CSCStrip4to5", "CSC Strip Ratio tbin 4 to tbin 5", 100, 0, 1);
23  theRatio6to5Plot = iBooker.book1D("CSCStrip6to5", "CSC Strip Ratio tbin 6 to tbin 5", 120, 0, 1.2);
24  theNDigisPerLayerPlot = iBooker.book1D("CSCStripDigisPerLayer", "Number of CSC Strip Digis per layer", 48, 0, 48);
25  theNDigisPerEventPlot = iBooker.book1D("CSCStripDigisPerEvent", "Number of CSC Strip Digis per event", 100, 0, 500);
26  if (doSim) {
27  for (int i = 0; i < 10; ++i) {
28  char title1[200];
29  sprintf(title1, "CSCStripDigiResolution%d", i + 1);
30  theResolutionPlots[i] = iBooker.book1D(title1, title1, 100, -5, 5);
31  }
32  }
33 }
34 
37  e.getByToken(strips_Token_, strips);
38  if (!strips.isValid()) {
39  edm::LogError("CSCDigiValidation") << "Cannot get strips by label " << theInputTag.encode();
40  }
41 
42  unsigned nDigisPerEvent = 0;
43 
44  for (CSCStripDigiCollection::DigiRangeIterator j = strips->begin(); j != strips->end(); j++) {
45  std::vector<CSCStripDigi>::const_iterator digiItr = (*j).second.first;
46  std::vector<CSCStripDigi>::const_iterator last = (*j).second.second;
47  int nDigis = last - digiItr;
48  nDigisPerEvent += nDigis;
49  theNDigisPerLayerPlot->Fill(nDigis);
50 
51  double maxAmplitude = 0.;
52  // int maxStrip = 0;
53 
54  for (; digiItr != last; ++digiItr) {
55  // average up the pedestals
56  std::vector<int> adcCounts = digiItr->getADCCounts();
57  thePedestalSum += adcCounts[0];
58  thePedestalSum += adcCounts[1];
59  thePedestalCount += 2;
61  if (adcCounts[4] - pedestal > maxAmplitude) {
62  // maxStrip = digiItr->getStrip();
63  maxAmplitude = adcCounts[4] - pedestal;
64  }
65 
66  // if we have enough pedestal statistics
67  if (thePedestalCount > 100) {
68  fillPedestalPlots(*digiItr);
69 
70  // see if it's big enough to count as "signal"
71  if (adcCounts[5] > (thePedestalSum / thePedestalCount + 100)) {
72  fillSignalPlots(*digiItr);
73  }
74  }
75  }
76  } // loop over digis
77 
78  theNDigisPerEventPlot->Fill(nDigisPerEvent);
79 }
80 
82  std::vector<int> adcCounts = digi.getADCCounts();
83  thePedestalPlot->Fill(adcCounts[0]);
84  thePedestalPlot->Fill(adcCounts[1]);
85 }
86 
88  std::vector<int> adcCounts = digi.getADCCounts();
90  theAmplitudePlot->Fill(adcCounts[4] - pedestal);
91  theRatio4to5Plot->Fill((adcCounts[3] - pedestal) / (adcCounts[4] - pedestal));
92  theRatio6to5Plot->Fill((adcCounts[5] - pedestal) / (adcCounts[4] - pedestal));
93 }
94 
95 void CSCStripDigiValidation::plotResolution(const PSimHit &hit, int strip, const CSCLayer *layer, int chamberType) {
96  double hitX = hit.localPosition().x();
97  double hitY = hit.localPosition().y();
98  double digiX = layer->geometry()->xOfStrip(strip, hitY);
99  theResolutionPlots[chamberType - 1]->Fill(digiX - hitX);
100 }
edm::InputTag theInputTag
MonitorElement * theNDigisPerEventPlot
MonitorElement * theRatio4to5Plot
std::vector< int > const & getADCCounts() const
Get ADC readings.
Definition: CSCStripDigi.h:44
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void bookHistograms(DQMStore::IBooker &, bool doSim)
MonitorElement * theRatio6to5Plot
#define nullptr
T y() const
Definition: PV3DBase.h:63
void fillSignalPlots(const CSCStripDigi &digi)
std::string encode() const
Definition: InputTag.cc:159
void Fill(long long x)
float xOfStrip(int strip, float y=0.) const
Local3DPoint localPosition() const
Definition: PSimHit.h:52
CSCStripDigiValidation(const edm::InputTag &inputTag, edm::ConsumesCollector &&iC)
edm::EDGetTokenT< CSCStripDigiCollection > strips_Token_
MonitorElement * thePedestalPlot
MonitorElement * theNDigisPerLayerPlot
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
bool isValid() const
Definition: HandleBase.h:74
void plotResolution(const PSimHit &hit, int strip, const CSCLayer *layer, int chamberType)
MonitorElement * theAmplitudePlot
void analyze(const edm::Event &e, const edm::EventSetup &) override
MonitorElement * theResolutionPlots[10]
void fillPedestalPlots(const CSCStripDigi &digi)
T x() const
Definition: PV3DBase.h:62
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47