CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCComparatorDigiValidation.cc
Go to the documentation of this file.
7 
8 
10  const edm::InputTag & inputTag, const edm::InputTag & stripDigiInputTag)
11 : CSCBaseValidation(dbe, inputTag),
12  theStripDigiInputTag(stripDigiInputTag),
13  theTimeBinPlots(),
14  theNDigisPerLayerPlots(),
15  theStripDigiPlots(),
16  the3StripPlots(),
17  theNDigisPerEventPlot( dbe_->book1D("CSCComparatorDigisPerEvent", "CSC Comparator Digis per event", 100, 0, 100) )
18 {
19  for(int i = 0; i < 10; ++i)
20  {
21  char title1[200], title2[200], title3[200], title4[200];
22  sprintf(title1, "CSCComparatorDigiTimeType%d", i+1);
23  sprintf(title2, "CSCComparatorDigisPerLayerType%d", i+1);
24  sprintf(title3, "CSCComparatorStripAmplitudeType%d", i+1);
25  sprintf(title4, "CSCComparator3StripAmplitudeType%d", i+1);
26  theTimeBinPlots[i] = dbe_->book1D(title1, title1, 9, 0, 8);
27  theNDigisPerLayerPlots[i] = dbe_->book1D(title2, title2, 100, 0, 20);
28  theStripDigiPlots[i] = dbe_->book1D(title3, title3, 100, 0, 1000);
29  the3StripPlots[i] = dbe_->book1D(title4, title4, 100, 0, 1000);
30 
31  }
32 }
33 
34 
35 
37 {
38 // for(int i = 0; i < 10; ++i)
39 // {
40 // edm::LogInfo("CSCDigiValidation") << "Mean of " << theTimeBinPlots[i]->getName()
41 // << " is " << theTimeBinPlots[i]->getMean()
42 // << " +/- " << theTimeBinPlots[i]->getRMS();
43 // edm::LogInfo("CSCDigiValidation") << "Mean charge of " << the3StripPlots[i]->getName()
44 // << " is " << the3StripPlots[i]->getMean();
45 // }
46 }
47 
48 
50 {
53 
54  e.getByLabel(theInputTag, comparators);
55  if (!comparators.isValid()) {
56  edm::LogError("CSCDigiDump") << "Cannot get comparators by label " << theInputTag.encode();
57  }
58  e.getByLabel(theStripDigiInputTag, stripDigis);
59  if (!stripDigis.isValid()) {
60  edm::LogError("CSCDigiDump") << "Cannot get comparators by label " << theInputTag.encode();
61  }
62 
63  unsigned nDigisPerEvent = 0;
64 
65  for (CSCComparatorDigiCollection::DigiRangeIterator j=comparators->begin(); j!=comparators->end(); j++) {
66  std::vector<CSCComparatorDigi>::const_iterator digiItr = (*j).second.first;
67  std::vector<CSCComparatorDigi>::const_iterator last = (*j).second.second;
68 
69  CSCDetId detId((*j).first);
70  const CSCLayer * layer = findLayer(detId.rawId());
71  int chamberType = layer->chamber()->specs()->chamberType();
72 
73  CSCStripDigiCollection::Range stripRange = stripDigis->get(detId);
74 
75  theNDigisPerLayerPlots[chamberType-1]->Fill(last-digiItr);
76 
77  for( ; digiItr != last; ++digiItr) {
78  ++nDigisPerEvent;
79  theTimeBinPlots[chamberType-1]->Fill(digiItr->getTimeBin());
80 
81  int strip = digiItr->getStrip();
82  for(std::vector<CSCStripDigi>::const_iterator stripItr = stripRange.first;
83  stripItr != stripRange.second; ++stripItr)
84  {
85  if(stripItr->getStrip() == strip)
86  {
87  std::vector<int> adc = stripItr->getADCCounts();
88  float pedc = 0.5*(adc[0]+adc[1]);
89  float amp = adc[4] - pedc;
90  theStripDigiPlots[chamberType-1]->Fill(amp);
91  // check neighbors
92  if(stripItr != stripRange.first && stripItr != stripRange.second-1)
93  {
94  std::vector<int> adcl = (stripItr-1)->getADCCounts();
95  std::vector<int> adcr = (stripItr+1)->getADCCounts();
96  float pedl = 0.5*(adcl[0]+adcl[1]);
97  float pedr = 0.5*(adcr[0]+adcr[1]);
98  float three = adcl[4]-pedl
99  + adcr[4]-pedr
100  + amp;
101  the3StripPlots[chamberType-1]->Fill(three);
102  }
103  }
104  }
105  }
106 
107  }
108 
109  theNDigisPerEventPlot->Fill(nDigisPerEvent);
110 }
111 
int adc(sample_type sample)
get the ADC sample (12 bits)
int i
Definition: DBlmapReader.cc:9
edm::InputTag theInputTag
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:717
CSCComparatorDigiValidation(DQMStore *dbe, const edm::InputTag &inputTag, const edm::InputTag &stripDigiInputTag)
void analyze(const edm::Event &, const edm::EventSetup &)
std::string encode() const
Definition: InputTag.cc:72
void Fill(long long x)
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:42
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:356
DQMStore * dbe_
int chamberType() const
std::pair< const_iterator, const_iterator > Range
const CSCChamber * chamber() const
Definition: CSCLayer.h:52
const CSCLayer * findLayer(int detId) const