Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:27

0001 // -*- C++ -*-
0002 //
0003 // Package:     SiPixelPhase1RecHits
0004 // Class:       SiPixelPhase1RecHits
0005 //
0006 
0007 // Original Author: Marcel Schneider
0008 
0009 #include "DQM/SiPixelPhase1Common/interface/SiPixelPhase1Base.h"
0010 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0011 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0012 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0015 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0016 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0017 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0018 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0019 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0020 #include "DataFormats/TrackReco/interface/Track.h"
0021 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0022 #include "DataFormats/VertexReco/interface/Vertex.h"
0023 
0024 namespace {
0025 
0026   class SiPixelPhase1RecHits final : public SiPixelPhase1Base {
0027     enum { NRECHITS, CLUST_X, CLUST_Y, ERROR_X, ERROR_Y, POS, CLUSTER_PROB, NONEDGE, NOTHERBAD };
0028 
0029   public:
0030     explicit SiPixelPhase1RecHits(const edm::ParameterSet& conf);
0031     void analyze(const edm::Event&, const edm::EventSetup&) override;
0032 
0033   private:
0034     edm::EDGetTokenT<reco::TrackCollection> srcToken_;
0035     edm::EDGetTokenT<reco::VertexCollection> offlinePrimaryVerticesToken_;
0036 
0037     edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeomToken_;
0038 
0039     bool onlyValid_;
0040     bool applyVertexCut_;
0041   };
0042 
0043   SiPixelPhase1RecHits::SiPixelPhase1RecHits(const edm::ParameterSet& iConfig) : SiPixelPhase1Base(iConfig) {
0044     srcToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("src"));
0045 
0046     offlinePrimaryVerticesToken_ = consumes<reco::VertexCollection>(std::string("offlinePrimaryVertices"));
0047 
0048     trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
0049 
0050     onlyValid_ = iConfig.getParameter<bool>("onlyValidHits");
0051 
0052     applyVertexCut_ = iConfig.getUntrackedParameter<bool>("VertexCut", true);
0053   }
0054 
0055   void SiPixelPhase1RecHits::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0056     if (!checktrigger(iEvent, iSetup, DCS))
0057       return;
0058 
0059     edm::ESHandle<TrackerGeometry> tracker = iSetup.getHandle(trackerGeomToken_);
0060     assert(tracker.isValid());
0061 
0062     edm::Handle<reco::TrackCollection> tracks;
0063     iEvent.getByToken(srcToken_, tracks);
0064     if (!tracks.isValid())
0065       return;
0066 
0067     edm::Handle<reco::VertexCollection> vertices;
0068     iEvent.getByToken(offlinePrimaryVerticesToken_, vertices);
0069 
0070     if (applyVertexCut_ && (!vertices.isValid() || vertices->empty()))
0071       return;
0072 
0073     for (auto const& track : *tracks) {
0074       if (applyVertexCut_ &&
0075           (track.pt() < 0.75 || std::abs(track.dxy(vertices->at(0).position())) > 5 * track.dxyError()))
0076         continue;
0077 
0078       bool isBpixtrack = false, isFpixtrack = false;
0079 
0080       auto const& trajParams = track.extra()->trajParams();
0081       auto hb = track.recHitsBegin();
0082       for (unsigned int h = 0; h < track.recHitsSize(); h++) {
0083         auto hit = *(hb + h);
0084         if (!trackerHitRTTI::isFromDet(*hit))
0085           continue;
0086 
0087         DetId id = hit->geographicalId();
0088         uint32_t subdetid = (id.subdetId());
0089 
0090         if (subdetid == PixelSubdetector::PixelBarrel)
0091           isBpixtrack = true;
0092         if (subdetid == PixelSubdetector::PixelEndcap)
0093           isFpixtrack = true;
0094       }
0095 
0096       if (!isBpixtrack && !isFpixtrack)
0097         continue;
0098 
0099       // then, look at each hit
0100       for (unsigned int h = 0; h < track.recHitsSize(); h++) {
0101         auto rechit = *(hb + h);
0102 
0103         if (!trackerHitRTTI::isFromDet(*rechit))
0104           continue;
0105 
0106         //continue if not a Pixel recHit
0107         DetId id = rechit->geographicalId();
0108         uint32_t subdetid = (id.subdetId());
0109 
0110         if (subdetid != PixelSubdetector::PixelBarrel && subdetid != PixelSubdetector::PixelEndcap)
0111           continue;
0112 
0113         bool isHitValid = rechit->getType() == TrackingRecHit::valid;
0114         if (onlyValid_ && !isHitValid)
0115           continue;  //useful to run on cosmics where the TrackEfficiency plugin is not used
0116 
0117         const SiPixelRecHit* prechit = dynamic_cast<const SiPixelRecHit*>(
0118             rechit);  //to be used to get the associated cluster and the cluster probability
0119 
0120         int sizeX = 0, sizeY = 0;
0121 
0122         if (isHitValid) {
0123           SiPixelRecHit::ClusterRef const& clust = prechit->cluster();
0124           sizeX = (*clust).sizeX();
0125           sizeY = (*clust).sizeY();
0126         }
0127 
0128         const PixelGeomDetUnit* geomdetunit = dynamic_cast<const PixelGeomDetUnit*>(tracker->idToDet(id));
0129         const PixelTopology& topol = geomdetunit->specificTopology();
0130 
0131         LocalPoint lp = trajParams[h].position();
0132         MeasurementPoint mp = topol.measurementPosition(lp);
0133 
0134         int row = (int)mp.x();
0135         int col = (int)mp.y();
0136 
0137         float rechit_x = lp.x();
0138         float rechit_y = lp.y();
0139 
0140         LocalError lerr = rechit->localPositionError();
0141         float lerr_x = sqrt(lerr.xx());
0142         float lerr_y = sqrt(lerr.yy());
0143 
0144         histo[NRECHITS].fill(id, &iEvent, col, row);  //in general a inclusive counter of missing/valid/inactive hits
0145         if (prechit->isOnEdge())
0146           histo[NONEDGE].fill(id, &iEvent, col, row);
0147         if (prechit->hasBadPixels())
0148           histo[NOTHERBAD].fill(id, &iEvent, col, row);
0149 
0150         if (isHitValid) {
0151           histo[CLUST_X].fill(sizeX, id, &iEvent, col, row);
0152           histo[CLUST_Y].fill(sizeY, id, &iEvent, col, row);
0153         }
0154 
0155         histo[ERROR_X].fill(lerr_x, id, &iEvent);
0156         histo[ERROR_Y].fill(lerr_y, id, &iEvent);
0157 
0158         histo[POS].fill(rechit_x, rechit_y, id, &iEvent);
0159 
0160         if (isHitValid) {
0161           double clusterProbability = prechit->clusterProbability(0);
0162           if (clusterProbability > 0)
0163             histo[CLUSTER_PROB].fill(log10(clusterProbability), id, &iEvent);
0164         }
0165       }
0166     }
0167 
0168     histo[NRECHITS].executePerEventHarvesting(&iEvent);
0169     histo[NONEDGE].executePerEventHarvesting(&iEvent);
0170     histo[NOTHERBAD].executePerEventHarvesting(&iEvent);
0171   }
0172 
0173 }  //namespace
0174 
0175 DEFINE_FWK_MODULE(SiPixelPhase1RecHits);