Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:48

0001 #include "DQMServices/Core/interface/DQMStore.h"
0002 #include "DataFormats/Common/interface/Handle.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "Validation/MuonCSCDigis/interface/CSCComparatorDigiValidation.h"
0005 #include "DataFormats/CSCDigi/interface/CSCConstants.h"
0006 
0007 CSCComparatorDigiValidation::CSCComparatorDigiValidation(const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
0008     : CSCBaseValidation(ps), theTimeBinPlots(), theNDigisPerLayerPlots(), theStripDigiPlots(), the3StripPlots() {
0009   const auto &comps = ps.getParameterSet("cscComparatorDigi");
0010   inputTagComp_ = comps.getParameter<edm::InputTag>("inputTag");
0011   comparators_Token_ = iC.consumes<CSCComparatorDigiCollection>(inputTagComp_);
0012 
0013   const auto &strips = ps.getParameterSet("cscStripDigi");
0014   inputTagStrip_ = strips.getParameter<edm::InputTag>("inputTag");
0015   strips_Token_ = iC.consumes<CSCStripDigiCollection>(inputTagStrip_);
0016 }
0017 
0018 CSCComparatorDigiValidation::~CSCComparatorDigiValidation() {}
0019 
0020 void CSCComparatorDigiValidation::bookHistograms(DQMStore::IBooker &iBooker) {
0021   iBooker.setCurrentFolder("MuonCSCDigisV/CSCDigiTask/Comparator/Occupancy/");
0022 
0023   theNDigisPerEventPlot = iBooker.book1D("CSCComparatorDigisPerEvent",
0024                                          "CSC Comparator Digis per event;CSC Comparator Digis per event;Entries",
0025                                          100,
0026                                          0,
0027                                          100);
0028   // 10 chamber types, if you consider ME1/a and ME1/b separate
0029   for (int i = 1; i <= 10; ++i) {
0030     const std::string t1("CSCComparatorDigiTime_" + CSCDetId::chamberName(i));
0031     const std::string t2("CSCComparatorDigisPerLayer_" + CSCDetId::chamberName(i));
0032     const std::string t3("CSCComparatorStripAmplitude_" + CSCDetId::chamberName(i));
0033     const std::string t4("CSCComparator3StripAmplitude_" + CSCDetId::chamberName(i));
0034     theTimeBinPlots[i - 1] = iBooker.book1D(
0035         t1, "Comparator Time Bin " + CSCDetId::chamberName(i) + " ;Comparator Time Bin; Entries", 16, 0, 16);
0036     theNDigisPerLayerPlots[i - 1] = iBooker.book1D(
0037         t2,
0038         "Number of Comparator Digis " + CSCDetId::chamberName(i) + " ;Number of Comparator Digis; Entries",
0039         100,
0040         0,
0041         20);
0042     theStripDigiPlots[i - 1] = iBooker.book1D(
0043         t3, "Comparator Amplitude " + CSCDetId::chamberName(i) + " ;Comparator Amplitude; Entries", 100, 0, 1000);
0044     the3StripPlots[i - 1] = iBooker.book1D(
0045         t4,
0046         "Comparator-triplet Amplitude " + CSCDetId::chamberName(i) + " ;Comparator-triplet Amplitude; Entries",
0047         100,
0048         0,
0049         1000);
0050   }
0051   if (doSim_) {
0052     iBooker.setCurrentFolder("MuonCSCDigisV/CSCDigiTask/Comparator/Resolution/");
0053     for (int i = 1; i <= 10; ++i) {
0054       const std::string t1("CSCComparatorPosResolution_" + CSCDetId::chamberName(i));
0055       theResolutionPlots[i - 1] = iBooker.book1D(t1,
0056                                                  "Comparator X Position Resolution " + CSCDetId::chamberName(i) +
0057                                                      ";Comparator X Position Resolution [cm]; Entries",
0058                                                  100,
0059                                                  -5,
0060                                                  5);
0061     }
0062   }
0063 }
0064 
0065 void CSCComparatorDigiValidation::analyze(const edm::Event &e, const edm::EventSetup &) {
0066   edm::Handle<CSCComparatorDigiCollection> comparators;
0067   edm::Handle<CSCStripDigiCollection> stripDigis;
0068 
0069   e.getByToken(comparators_Token_, comparators);
0070   if (!comparators.isValid()) {
0071     edm::LogError("CSCComparatorDigiValidation") << "Cannot get comparators by label " << inputTagComp_.encode();
0072   }
0073   e.getByToken(strips_Token_, stripDigis);
0074   if (!stripDigis.isValid()) {
0075     edm::LogError("CSCComparatorDigiValidation") << "Cannot get strips by label " << inputTagComp_.encode();
0076   }
0077 
0078   unsigned nDigisPerEvent = 0;
0079 
0080   for (auto j = comparators->begin(); j != comparators->end(); j++) {
0081     auto digiItr = (*j).second.first;
0082     auto last = (*j).second.second;
0083 
0084     CSCDetId detId((*j).first);
0085     const CSCLayer *layer = findLayer(detId.rawId());
0086     int chamberType = layer->chamber()->specs()->chamberType();
0087     int nDigis = last - digiItr;
0088 
0089     theNDigisPerLayerPlots[chamberType - 1]->Fill(last - digiItr);
0090 
0091     // require exactly one comparator digi per layer
0092     if (doSim_) {
0093       const edm::PSimHitContainer &simHits = theSimHitMap->hits(detId);
0094       if (nDigis == 1 && simHits.size() == 1) {
0095         const int fullStrip(digiItr->getStrip() + digiItr->getCFEB() * CSCConstants::NUM_STRIPS_PER_CFEB);
0096         plotResolution(simHits[0], fullStrip, layer, chamberType);
0097       }
0098     }
0099 
0100     for (auto stripRange = stripDigis->get(detId); digiItr != last; ++digiItr) {
0101       ++nDigisPerEvent;
0102       theTimeBinPlots[chamberType - 1]->Fill(digiItr->getTimeBin());
0103 
0104       int strip = digiItr->getStrip();
0105       for (auto stripItr = stripRange.first; stripItr != stripRange.second; ++stripItr) {
0106         if (stripItr->getStrip() == strip) {
0107           std::vector<int> adc = stripItr->getADCCounts();
0108           float pedc = 0.5 * (adc[0] + adc[1]);
0109           float amp = adc[4] - pedc;
0110           theStripDigiPlots[chamberType - 1]->Fill(amp);
0111           // check neighbors
0112           if (stripItr != stripRange.first && stripItr != stripRange.second - 1) {
0113             std::vector<int> adcl = (stripItr - 1)->getADCCounts();
0114             std::vector<int> adcr = (stripItr + 1)->getADCCounts();
0115             float pedl = 0.5 * (adcl[0] + adcl[1]);
0116             float pedr = 0.5 * (adcr[0] + adcr[1]);
0117             float three = adcl[4] - pedl + adcr[4] - pedr + amp;
0118             the3StripPlots[chamberType - 1]->Fill(three);
0119           }
0120         }
0121       }
0122     }
0123   }
0124   theNDigisPerEventPlot->Fill(nDigisPerEvent);
0125 }
0126 
0127 void CSCComparatorDigiValidation::plotResolution(const PSimHit &hit, int strip, const CSCLayer *layer, int chamberType) {
0128   double hitX = hit.localPosition().x();
0129   double hitY = hit.localPosition().y();
0130   double digiX = layer->geometry()->xOfStrip(strip, hitY);
0131   theResolutionPlots[chamberType - 1]->Fill(digiX - hitX);
0132 }