CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CSCComparatorDigiValidation.cc
Go to the documentation of this file.
6 
8  : CSCBaseValidation(ps), theTimeBinPlots(), theNDigisPerLayerPlots(), theStripDigiPlots(), the3StripPlots() {
9  const auto &comps = ps.getParameterSet("cscComparatorDigi");
10  inputTagComp_ = comps.getParameter<edm::InputTag>("inputTag");
12 
13  const auto &strips = ps.getParameterSet("cscStripDigi");
14  inputTagStrip_ = strips.getParameter<edm::InputTag>("inputTag");
16 }
17 
19 
21  iBooker.setCurrentFolder("MuonCSCDigisV/CSCDigiTask/Comparator/Occupancy/");
22 
23  theNDigisPerEventPlot = iBooker.book1D("CSCComparatorDigisPerEvent",
24  "CSC Comparator Digis per event;CSC Comparator Digis per event;Entries",
25  100,
26  0,
27  100);
28  // 10 chamber types, if you consider ME1/a and ME1/b separate
29  for (int i = 1; i <= 10; ++i) {
30  const std::string t1("CSCComparatorDigiTime_" + CSCDetId::chamberName(i));
31  const std::string t2("CSCComparatorDigisPerLayer_" + CSCDetId::chamberName(i));
32  const std::string t3("CSCComparatorStripAmplitude_" + CSCDetId::chamberName(i));
33  const std::string t4("CSCComparator3StripAmplitude_" + CSCDetId::chamberName(i));
34  theTimeBinPlots[i - 1] = iBooker.book1D(
35  t1, "Comparator Time Bin " + CSCDetId::chamberName(i) + " ;Comparator Time Bin; Entries", 16, 0, 16);
36  theNDigisPerLayerPlots[i - 1] = iBooker.book1D(
37  t2,
38  "Number of Comparator Digis " + CSCDetId::chamberName(i) + " ;Number of Comparator Digis; Entries",
39  100,
40  0,
41  20);
42  theStripDigiPlots[i - 1] = iBooker.book1D(
43  t3, "Comparator Amplitude " + CSCDetId::chamberName(i) + " ;Comparator Amplitude; Entries", 100, 0, 1000);
44  the3StripPlots[i - 1] = iBooker.book1D(
45  t4,
46  "Comparator-triplet Amplitude " + CSCDetId::chamberName(i) + " ;Comparator-triplet Amplitude; Entries",
47  100,
48  0,
49  1000);
50  }
51  if (doSim_) {
52  iBooker.setCurrentFolder("MuonCSCDigisV/CSCDigiTask/Comparator/Resolution/");
53  for (int i = 1; i <= 10; ++i) {
54  const std::string t1("CSCComparatorPosResolution_" + CSCDetId::chamberName(i));
55  theResolutionPlots[i - 1] = iBooker.book1D(t1,
56  "Comparator X Position Resolution " + CSCDetId::chamberName(i) +
57  ";Comparator X Position Resolution [cm]; Entries",
58  100,
59  -5,
60  5);
61  }
62  }
63 }
64 
68 
69  e.getByToken(comparators_Token_, comparators);
70  if (!comparators.isValid()) {
71  edm::LogError("CSCComparatorDigiValidation") << "Cannot get comparators by label " << inputTagComp_.encode();
72  }
73  e.getByToken(strips_Token_, stripDigis);
74  if (!stripDigis.isValid()) {
75  edm::LogError("CSCComparatorDigiValidation") << "Cannot get strips by label " << inputTagComp_.encode();
76  }
77 
78  unsigned nDigisPerEvent = 0;
79 
80  for (auto j = comparators->begin(); j != comparators->end(); j++) {
81  auto digiItr = (*j).second.first;
82  auto last = (*j).second.second;
83 
84  CSCDetId detId((*j).first);
85  const CSCLayer *layer = findLayer(detId.rawId());
86  int chamberType = layer->chamber()->specs()->chamberType();
87  int nDigis = last - digiItr;
88 
89  theNDigisPerLayerPlots[chamberType - 1]->Fill(last - digiItr);
90 
91  // require exactly one comparator digi per layer
92  if (doSim_) {
94  if (nDigis == 1 && simHits.size() == 1) {
95  const int fullStrip(digiItr->getStrip() + digiItr->getCFEB() * CSCConstants::NUM_STRIPS_PER_CFEB);
96  plotResolution(simHits[0], fullStrip, layer, chamberType);
97  }
98  }
99 
100  for (auto stripRange = stripDigis->get(detId); digiItr != last; ++digiItr) {
101  ++nDigisPerEvent;
102  theTimeBinPlots[chamberType - 1]->Fill(digiItr->getTimeBin());
103 
104  int strip = digiItr->getStrip();
105  for (auto stripItr = stripRange.first; stripItr != stripRange.second; ++stripItr) {
106  if (stripItr->getStrip() == strip) {
107  std::vector<int> adc = stripItr->getADCCounts();
108  float pedc = 0.5 * (adc[0] + adc[1]);
109  float amp = adc[4] - pedc;
110  theStripDigiPlots[chamberType - 1]->Fill(amp);
111  // check neighbors
112  if (stripItr != stripRange.first && stripItr != stripRange.second - 1) {
113  std::vector<int> adcl = (stripItr - 1)->getADCCounts();
114  std::vector<int> adcr = (stripItr + 1)->getADCCounts();
115  float pedl = 0.5 * (adcl[0] + adcl[1]);
116  float pedr = 0.5 * (adcr[0] + adcr[1]);
117  float three = adcl[4] - pedl + adcr[4] - pedr + amp;
118  the3StripPlots[chamberType - 1]->Fill(three);
119  }
120  }
121  }
122  }
123  }
124  theNDigisPerEventPlot->Fill(nDigisPerEvent);
125 }
126 
127 void CSCComparatorDigiValidation::plotResolution(const PSimHit &hit, int strip, const CSCLayer *layer, int chamberType) {
128  double hitX = hit.localPosition().x();
129  double hitY = hit.localPosition().y();
130  double digiX = layer->geometry()->xOfStrip(strip, hitY);
131  theResolutionPlots[chamberType - 1]->Fill(digiX - hitX);
132 }
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T y() const
Definition: PV3DBase.h:60
Log< level::Error, false > LogError
std::string encode() const
Definition: InputTag.cc:159
constexpr std::array< uint8_t, layerIndexSize > layer
float xOfStrip(int strip, float y=0.) const
void Fill(long long x)
void analyze(const edm::Event &, const edm::EventSetup &) override
void plotResolution(const PSimHit &hit, int strip, const CSCLayer *layer, int chamberType)
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:39
edm::EDGetTokenT< CSCComparatorDigiCollection > comparators_Token_
Local3DPoint localPosition() const
Definition: PSimHit.h:52
const PSimHitMap * theSimHitMap
tuple strips
#turn off noise in all subdetectors simHcalUnsuppressedDigis.doNoise = False mix.digitizers.hcal.doNoise = False simEcalUnsuppressedDigis.doNoise = False mix.digitizers.ecal.doNoise = False simEcalUnsuppressedDigis.doESNoise = False simSiPixelDigis.AddNoise = False mix.digitizers.pixel.AddNoise = False simSiStripDigis.Noise = False mix.digitizers.strip.AddNoise = False
Definition: DigiDM_cff.py:32
const edm::PSimHitContainer & hits(int detId) const
Definition: PSimHitMap.cc:20
void bookHistograms(DQMStore::IBooker &)
bool isValid() const
Definition: HandleBase.h:70
std::string chamberName() const
Definition: CSCDetId.cc:92
int chamberType() const
ParameterSet const & getParameterSet(std::string const &) const
tuple simHits
Definition: trackerHits.py:16
std::vector< PSimHit > PSimHitContainer
CSCComparatorDigiValidation(const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
tuple last
Definition: dqmdumpme.py:56
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
edm::EDGetTokenT< CSCStripDigiCollection > strips_Token_
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
uint16_t *__restrict__ uint16_t const *__restrict__ adc