CMS 3D CMS Logo

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