CMS 3D CMS Logo

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