CMS 3D CMS Logo

CSCWireDigiValidation.cc
Go to the documentation of this file.
4 
8 
11  bool doSim):
12  CSCBaseValidation(inputTag),
13  doSim_(doSim),
14  theTimeBinPlots(),
15  theNDigisPerLayerPlots()
16 {
18 }
19 
21 {
22 }
23 
25 {
26  theNDigisPerEventPlot = iBooker.book1D("CSCWireDigisPerEvent", "CSC Wire Digis per event", 100, 0, 100);
27  for(int i = 0; i < 10; ++i)
28  {
29  char title1[200], title2[200], title3[200];
30  sprintf(title1, "CSCWireDigiTimeType%d", i+1);
31  sprintf(title2, "CSCWireDigisPerLayerType%d", i+1);
32  sprintf(title3, "CSCWireDigiResolution%d", i+1);
33  theTimeBinPlots[i] = iBooker.book1D(title1, title1, 9, 0, 8);
34  theNDigisPerLayerPlots[i] = iBooker.book1D(title2, title2, 100, 0, 20);
35  theResolutionPlots[i] = iBooker.book1D(title3, title3, 100, -10, 10);
36  }
37 }
38 
40 {
42 
43  e.getByToken(wires_Token_, wires);
44 
45  if (!wires.isValid()) {
46  edm::LogError("CSCDigiDump") << "Cannot get wires by label " << theInputTag.encode();
47  }
48 
49  unsigned nDigisPerEvent = 0;
50 
51  for (CSCWireDigiCollection::DigiRangeIterator j=wires->begin(); j!=wires->end(); j++) {
52  std::vector<CSCWireDigi>::const_iterator beginDigi = (*j).second.first;
53  std::vector<CSCWireDigi>::const_iterator endDigi = (*j).second.second;
54  int detId = (*j).first.rawId();
55 
56  const CSCLayer * layer = findLayer(detId);
57  int chamberType = layer->chamber()->specs()->chamberType();
58  int nDigis = endDigi-beginDigi;
59  nDigisPerEvent += nDigis;
60  theNDigisPerLayerPlots[chamberType-1]->Fill(nDigis);
61 
62  for( std::vector<CSCWireDigi>::const_iterator digiItr = beginDigi;
63  digiItr != endDigi; ++digiItr)
64  {
65  theTimeBinPlots[chamberType-1]->Fill(digiItr->getTimeBin());
66  }
67 
68  if(doSim_)
69  {
71  if(nDigis == 1 && simHits.size() == 1)
72  {
73  plotResolution(simHits[0], *beginDigi, layer, chamberType);
74  }
75  }
76  }
77 
78  theNDigisPerEventPlot->Fill(nDigisPerEvent);
79 }
80 
81 
83  const CSCWireDigi & digi,
84  const CSCLayer * layer,
85  int chamberType)
86 {
87  double hitX = hit.localPosition().x();
88  double hitY = hit.localPosition().y();
89  double digiY = layer->geometry()->yOfWireGroup(digi.getWireGroup(), hitX);
90  theResolutionPlots[chamberType-1]->Fill(digiY - hitY);
91 }
edm::InputTag theInputTag
MonitorElement * theNDigisPerEventPlot
MonitorElement * theResolutionPlots[10]
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
T y() const
Definition: PV3DBase.h:63
std::string encode() const
Definition: InputTag.cc:166
float yOfWireGroup(int wireGroup, float x=0.) const
void Fill(long long x)
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:42
Local3DPoint localPosition() const
Definition: PSimHit.h:44
const PSimHitMap * theSimHitMap
edm::EDGetTokenT< CSCWireDigiCollection > wires_Token_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:22
MonitorElement * theNDigisPerLayerPlots[10]
bool isValid() const
Definition: HandleBase.h:74
CSCWireDigiValidation(const edm::InputTag &inputTag, edm::ConsumesCollector &&iC, bool doSim)
void plotResolution(const PSimHit &hit, const CSCWireDigi &digi, const CSCLayer *layer, int chamberType)
MonitorElement * theTimeBinPlots[10]
int chamberType() const
int getWireGroup() const
default
Definition: CSCWireDigi.h:24
void bookHistograms(DQMStore::IBooker &)
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< PSimHit > PSimHitContainer
T x() const
Definition: PV3DBase.h:62
const CSCChamber * chamber() const
Definition: CSCLayer.h:52
const CSCLayerGeometry * geometry() const
Definition: CSCLayer.h:47
const CSCLayer * findLayer(int detId) const