Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:57

0001 #include "DQMServices/Core/interface/DQMStore.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "Validation/CSCRecHits/interface/CSCRecHit2DValidation.h"
0004 
0005 CSCRecHit2DValidation::CSCRecHit2DValidation(const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
0006     : CSCBaseValidation(ps), theNPerEventPlot(nullptr) {
0007   const auto &pset = ps.getParameterSet("cscRecHit");
0008   inputTag_ = pset.getParameter<edm::InputTag>("inputTag");
0009   rechits_Token_ = iC.consumes<CSCRecHit2DCollection>(inputTag_);
0010 }
0011 
0012 CSCRecHit2DValidation::~CSCRecHit2DValidation() {}
0013 
0014 void CSCRecHit2DValidation::bookHistograms(DQMStore::IBooker &iBooker) {
0015   theNPerEventPlot = iBooker.book1D("CSCRecHitsPerEvent", "Number of CSC Rec Hits per event", 100, 0, 500);
0016   // 10 chamber types, if you consider ME1/a and ME1/b separate
0017   for (int i = 1; i <= 10; ++i) {
0018     const std::string cn(CSCDetId::chamberName(i));
0019     const std::string t1("CSCRecHitResolution_" + cn);
0020     const std::string t2("CSCRecHitPull_" + cn);
0021     const std::string t3("CSCRecHitYResolution_" + cn);
0022 
0023     const std::string t4("CSCRecHitYPull_" + cn);
0024     const std::string t5("CSCRecHitPosInStrip_" + cn);
0025     const std::string t6("CSCSimHitPosInStrip_" + cn);
0026 
0027     const std::string t7("CSCRecHit_" + cn);
0028     const std::string t8("CSCSimHit_" + cn);
0029     const std::string t9("CSCTPeak_" + cn);
0030 
0031     theResolutionPlots[i - 1] = iBooker.book1D(t1, t1 + ";R*dPhi Resolution [cm];Entries", 100, -0.2, 0.2);
0032     thePullPlots[i - 1] = iBooker.book1D(t2, t2 + ";dPhi Pull;Entries", 100, -3, 3);
0033     theYResolutionPlots[i - 1] = iBooker.book1D(t3, t3 + ";Local Y Resolution [cm];Entries", 100, -5, 5);
0034     theYPullPlots[i - 1] = iBooker.book1D(t4, t4 + ";Local Y Pull;Entries", 100, -3, 3);
0035     theRecHitPosInStrip[i - 1] = iBooker.book1D(t5, t5 + ";Position in Strip;Entries", 100, -2, 2);
0036     theSimHitPosInStrip[i - 1] = iBooker.book1D(t6, t6 + ";Position in Strip;Entries", 100, -2, 2);
0037 
0038     theScatterPlots[i - 1] = iBooker.book2D(t7, t7 + ";Local Phi;Local Y [cm]", 200, -20, 20, 200, -250, 250);
0039     theSimHitScatterPlots[i - 1] = iBooker.book2D(t8, t8 + ";Local Phi;Local Y [cm]", 200, -20, 20, 200, -250, 250);
0040     theTPeaks[i - 1] = iBooker.book1D(t9, t9 + ";Peak Time [ns];Entries", 200, 0, 400);
0041   }
0042 }
0043 
0044 void CSCRecHit2DValidation::analyze(const edm::Event &e, const edm::EventSetup &eventSetup) {
0045   // get the collection of CSCRecHrecHitItrD
0046   edm::Handle<CSCRecHit2DCollection> hRecHits;
0047   e.getByToken(rechits_Token_, hRecHits);
0048   const CSCRecHit2DCollection *cscRecHits = hRecHits.product();
0049 
0050   unsigned nPerEvent = 0;
0051 
0052   for (auto recHitItr = cscRecHits->begin(); recHitItr != cscRecHits->end(); recHitItr++) {
0053     ++nPerEvent;
0054     int detId = (*recHitItr).cscDetId().rawId();
0055     edm::PSimHitContainer simHits = theSimHitMap->hits(detId);
0056     const CSCLayer *layer = findLayer(detId);
0057     int chamberType = layer->chamber()->specs()->chamberType();
0058     theTPeaks[chamberType - 1]->Fill(recHitItr->tpeak());
0059     if (simHits.size() == 1) {
0060       plotResolution(simHits[0], *recHitItr, layer, chamberType);
0061     }
0062     float localX = recHitItr->localPosition().x();
0063     float localY = recHitItr->localPosition().y();
0064     // find a local phi
0065     float globalR = layer->toGlobal(LocalPoint(0., 0., 0.)).perp();
0066     GlobalPoint axisThruChamber(globalR + localY, localX, 0.);
0067     float localPhi = axisThruChamber.phi().degrees();
0068     theScatterPlots[chamberType - 1]->Fill(localPhi, localY);
0069   }
0070   theNPerEventPlot->Fill(nPerEvent);
0071 
0072   if (doSim_) {
0073     // fill sim hits
0074     std::vector<int> layersWithSimHits = theSimHitMap->detsWithHits();
0075     for (unsigned i = 0; i < layersWithSimHits.size(); ++i) {
0076       edm::PSimHitContainer simHits = theSimHitMap->hits(layersWithSimHits[i]);
0077       for (auto hitItr = simHits.begin(); hitItr != simHits.end(); ++hitItr) {
0078         const CSCLayer *layer = findLayer(layersWithSimHits[i]);
0079         int chamberType = layer->chamber()->specs()->chamberType();
0080         float localX = hitItr->localPosition().x();
0081         float localY = hitItr->localPosition().y();
0082         // find a local phi
0083         float globalR = layer->toGlobal(LocalPoint(0., 0., 0.)).perp();
0084         GlobalPoint axisThruChamber(globalR + localY, localX, 0.);
0085         float localPhi = axisThruChamber.phi().degrees();
0086         theSimHitScatterPlots[chamberType - 1]->Fill(localPhi, localY);
0087       }
0088     }
0089   }
0090 }
0091 
0092 void CSCRecHit2DValidation::plotResolution(const PSimHit &simHit,
0093                                            const CSCRecHit2D &recHit,
0094                                            const CSCLayer *layer,
0095                                            int chamberType) {
0096   GlobalPoint simHitPos = layer->toGlobal(simHit.localPosition());
0097   GlobalPoint recHitPos = layer->toGlobal(recHit.localPosition());
0098 
0099   double dphi = recHitPos.phi() - simHitPos.phi();
0100   double rdphi = recHitPos.perp() * dphi;
0101   theResolutionPlots[chamberType - 1]->Fill(rdphi);
0102   thePullPlots[chamberType - 1]->Fill(rdphi / sqrt(recHit.localPositionError().xx()));
0103   double dy = recHit.localPosition().y() - simHit.localPosition().y();
0104   theYResolutionPlots[chamberType - 1]->Fill(dy);
0105   theYPullPlots[chamberType - 1]->Fill(dy / sqrt(recHit.localPositionError().yy()));
0106 
0107   const CSCLayerGeometry *layerGeometry = layer->geometry();
0108   float recStrip = layerGeometry->strip(recHit.localPosition());
0109   float simStrip = layerGeometry->strip(simHit.localPosition());
0110   theRecHitPosInStrip[chamberType - 1]->Fill(recStrip - int(recStrip));
0111   theSimHitPosInStrip[chamberType - 1]->Fill(simStrip - int(simStrip));
0112 }