CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCStripDigiValidation.cc
Go to the documentation of this file.
5 
6 
8  const edm::InputTag & inputTag,
10  bool doSim)
11 : CSCBaseValidation(dbe, inputTag),
12  thePedestalSum(0),
13  thePedestalCovarianceSum(0),
14  thePedestalCount(0),
15  theDoSimFlag(doSim),
16  thePedestalPlot( dbe_->book1D("CSCPedestal", "CSC Pedestal ", 400, 550, 650) ),
17  thePedestalTimeCorrelationPlot(0),
18  thePedestalNeighborCorrelationPlot(0),
19  theAmplitudePlot( dbe_->book1D("CSCStripAmplitude", "CSC Strip Amplitude", 200, 0, 2000) ),
20  theRatio4to5Plot( dbe_->book1D("CSCStrip4to5", "CSC Strip Ratio tbin 4 to tbin 5", 100, 0, 1) ),
21  theRatio6to5Plot( dbe_->book1D("CSCStrip6to5", "CSC Strip Ratio tbin 6 to tbin 5", 120, 0, 1.2) ),
22  theNDigisPerLayerPlot( dbe_->book1D("CSCStripDigisPerLayer", "Number of CSC Strip Digis per layer", 48, 0, 48) ),
23  theNDigisPerChamberPlot(0),
24  theNDigisPerEventPlot( dbe_->book1D("CSCStripDigisPerEvent", "Number of CSC Strip Digis per event", 100, 0, 500) )
25 {
26  strips_Token_ = iC.consumes<CSCStripDigiCollection>(inputTag);
27 
28  if(doSim) {
29  for(int i = 0; i < 10; ++i)
30  {
31  char title1[200];
32  sprintf(title1, "CSCStripDigiResolution%d", i+1);
33  theResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -5, 5);
34  }
35  }
36 
37 }
38 
39 
41 
42 // edm::LogInfo("CSCDigiValidation") << "RATIO for strips 4 to 5 : " << theRatio4to5Plot->getMean();
43 // edm::LogInfo("CSCDigiValidation") << "RATIO for strips 6 to 5 : " << theRatio6to5Plot->getMean();
44 // edm::LogInfo("CSCDigiValidation") << "NDIGIS per event : " << theNDigisPerEventPlot->getMean();
45 
46 }
47 
48 
50  const edm::EventSetup&)
51 {
53  e.getByToken(strips_Token_, strips);
54  if (!strips.isValid()) {
55  edm::LogError("CSCDigiValidation") << "Cannot get strips by label "
56  << theInputTag.encode();
57  }
58 
59  unsigned nDigisPerEvent = 0;
60 
61  for (CSCStripDigiCollection::DigiRangeIterator j=strips->begin(); j!=strips->end(); j++) {
62  std::vector<CSCStripDigi>::const_iterator digiItr = (*j).second.first;
63  std::vector<CSCStripDigi>::const_iterator last = (*j).second.second;
64  int nDigis = last-digiItr;
65  nDigisPerEvent += nDigis;
66  theNDigisPerLayerPlot->Fill(nDigis);
67 
68  double maxAmplitude = 0.;
69  // int maxStrip = 0;
70 
71  for( ; digiItr != last; ++digiItr) {
72  // average up the pedestals
73  std::vector<int> adcCounts = digiItr->getADCCounts();
74  thePedestalSum += adcCounts[0];
75  thePedestalSum += adcCounts[1];
76  thePedestalCount += 2;
78  if(adcCounts[4]-pedestal > maxAmplitude)
79  {
80  // maxStrip = digiItr->getStrip();
81  maxAmplitude = adcCounts[4]-pedestal;
82  }
83 
84  // if we have enough pedestal statistics
85  if(thePedestalCount > 100)
86  {
87  fillPedestalPlots(*digiItr);
88 
89  // see if it's big enough to count as "signal"
90  if(adcCounts[5] > (thePedestalSum/thePedestalCount + 100))
91  {
92  fillSignalPlots(*digiItr);
93  }
94  }
95  }
96 
97 /*
98  int detId = (*j).first.rawId();
99 
100  edm::PSimHitContainer simHits = theSimHitMap->hits(detId);
101 
102  if(simHits.size() == 1)
103  {
104  const CSCLayer * layer = findLayer(detId);
105  int chamberType = layer->chamber()->specs()->chamberType();
106  plotResolution(simHits[0], maxStrip, layer, chamberType);
107  }
108 */
109  } // loop over digis
110 
111  theNDigisPerEventPlot->Fill(nDigisPerEvent);
112 }
113 
114 
116 {
117  std::vector<int> adcCounts = digi.getADCCounts();
118  thePedestalPlot->Fill(adcCounts[0]);
119  thePedestalPlot->Fill(adcCounts[1]);
120 }
121 
122 
123 
125 {
126  std::vector<int> adcCounts = digi.getADCCounts();
128  theAmplitudePlot->Fill(adcCounts[4] - pedestal);
129  theRatio4to5Plot->Fill( (adcCounts[3]-pedestal) / (adcCounts[4]-pedestal) );
130  theRatio6to5Plot->Fill( (adcCounts[5]-pedestal) / (adcCounts[4]-pedestal) );
131 }
132 
133 
135  const CSCLayer * layer, int chamberType)
136 {
137  double hitX = hit.localPosition().x();
138  double hitY = hit.localPosition().y();
139  double digiX = layer->geometry()->xOfStrip(strip, hitY);
140  theResolutionPlots[chamberType-1]->Fill(digiX - hitX);
141 }
142 
143 
int i
Definition: DBlmapReader.cc:9
edm::InputTag theInputTag
MonitorElement * theNDigisPerEventPlot
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:954
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:434
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:164
void Fill(long long x)
float xOfStrip(int strip, float y=0.) const
Local3DPoint localPosition() const
Definition: PSimHit.h:44
edm::EDGetTokenT< CSCStripDigiCollection > strips_Token_
MonitorElement * thePedestalPlot
MonitorElement * theNDigisPerLayerPlot
int j
Definition: DBlmapReader.cc:9
bool isValid() const
Definition: HandleBase.h:76
DQMStore * dbe_
void plotResolution(const PSimHit &hit, int strip, const CSCLayer *layer, int chamberType)
MonitorElement * theAmplitudePlot
MonitorElement * theResolutionPlots[10]
CSCStripDigiValidation(DQMStore *dbe, const edm::InputTag &inputTag, edm::ConsumesCollector &&iC, bool doSim)
void fillPedestalPlots(const CSCStripDigi &digi)
T x() const
Definition: PV3DBase.h:62
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47