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.
6 
7 
9  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  if(doSim) {
27  for(int i = 0; i < 10; ++i)
28  {
29  char title1[200];
30  sprintf(title1, "CSCStripDigiResolution%d", i+1);
31  theResolutionPlots[i] = dbe_->book1D(title1, title1, 100, -5, 5);
32  }
33  }
34 
35 }
36 
37 
39 
40 // edm::LogInfo("CSCDigiValidation") << "RATIO for strips 4 to 5 : " << theRatio4to5Plot->getMean();
41 // edm::LogInfo("CSCDigiValidation") << "RATIO for strips 6 to 5 : " << theRatio6to5Plot->getMean();
42 // edm::LogInfo("CSCDigiValidation") << "NDIGIS per event : " << theNDigisPerEventPlot->getMean();
43 
44 }
45 
46 
48  const edm::EventSetup&)
49 {
51  e.getByLabel(theInputTag, strips);
52  if (!strips.isValid()) {
53  edm::LogError("CSCDigiValidation") << "Cannot get strips by label "
54  << theInputTag.encode();
55  }
56 
57  unsigned nDigisPerEvent = 0;
58 
59  for (CSCStripDigiCollection::DigiRangeIterator j=strips->begin(); j!=strips->end(); j++) {
60  std::vector<CSCStripDigi>::const_iterator digiItr = (*j).second.first;
61  std::vector<CSCStripDigi>::const_iterator last = (*j).second.second;
62  int nDigis = last-digiItr;
63  nDigisPerEvent += nDigis;
64  theNDigisPerLayerPlot->Fill(nDigis);
65 
66  double maxAmplitude = 0.;
67  int maxStrip = 0;
68 
69  for( ; digiItr != last; ++digiItr) {
70  // average up the pedestals
71  std::vector<int> adcCounts = digiItr->getADCCounts();
72  thePedestalSum += adcCounts[0];
73  thePedestalSum += adcCounts[1];
74  thePedestalCount += 2;
75  float pedestal = thePedestalSum/thePedestalCount;
76  if(adcCounts[4]-pedestal > maxAmplitude)
77  {
78  maxStrip = digiItr->getStrip();
79  maxAmplitude = adcCounts[4]-pedestal;
80  }
81 
82  // if we have enough pedestal statistics
83  if(thePedestalCount > 100)
84  {
85  fillPedestalPlots(*digiItr);
86 
87  // see if it's big enough to count as "signal"
88  if(adcCounts[5] > (thePedestalSum/thePedestalCount + 100))
89  {
90  fillSignalPlots(*digiItr);
91  }
92  }
93  }
94 
95 /*
96  int detId = (*j).first.rawId();
97 
98  edm::PSimHitContainer simHits = theSimHitMap->hits(detId);
99 
100  if(simHits.size() == 1)
101  {
102  const CSCLayer * layer = findLayer(detId);
103  int chamberType = layer->chamber()->specs()->chamberType();
104  plotResolution(simHits[0], maxStrip, layer, chamberType);
105  }
106 */
107  } // loop over digis
108 
109  theNDigisPerEventPlot->Fill(nDigisPerEvent);
110 }
111 
112 
114 {
115  std::vector<int> adcCounts = digi.getADCCounts();
116  thePedestalPlot->Fill(adcCounts[0]);
117  thePedestalPlot->Fill(adcCounts[1]);
118 }
119 
120 
121 
123 {
124  std::vector<int> adcCounts = digi.getADCCounts();
125  float pedestal = thePedestalSum/thePedestalCount;
126  theAmplitudePlot->Fill(adcCounts[4] - pedestal);
127  theRatio4to5Plot->Fill( (adcCounts[3]-pedestal) / (adcCounts[4]-pedestal) );
128  theRatio6to5Plot->Fill( (adcCounts[5]-pedestal) / (adcCounts[4]-pedestal) );
129 }
130 
131 
133  const CSCLayer * layer, int chamberType)
134 {
135  double hitX = hit.localPosition().x();
136  double hitY = hit.localPosition().y();
137  double digiX = layer->geometry()->xOfStrip(strip, hitY);
138  theResolutionPlots[chamberType-1]->Fill(digiX - hitX);
139 }
140 
141 
int i
Definition: DBlmapReader.cc:9
edm::InputTag theInputTag
MonitorElement * theNDigisPerEventPlot
void strip(std::string &input, const std::string &blanks=" \n\t")
Definition: stringTools.cc:16
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
MonitorElement * theRatio4to5Plot
std::vector< int > getADCCounts() const
Get ADC readings.
Definition: CSCStripDigi.cc:41
MonitorElement * theRatio6to5Plot
T y() const
Definition: PV3DBase.h:57
void analyze(const edm::Event &e, const edm::EventSetup &)
void fillSignalPlots(const CSCStripDigi &digi)
CSCStripDigiValidation(DQMStore *dbe, const edm::InputTag &inputTag, bool doSim)
std::string encode() const
Definition: InputTag.cc:72
void Fill(long long x)
float xOfStrip(int strip, float y=0.) const
Local3DPoint localPosition() const
Definition: PSimHit.h:44
MonitorElement * thePedestalPlot
MonitorElement * theNDigisPerLayerPlot
int j
Definition: DBlmapReader.cc:9
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
DQMStore * dbe_
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:56
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47