Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-03 00:42:13

0001 // TODO: change file name to SiPixelCompareRecHitsSoA.cc when CUDA code is removed
0002 
0003 #include "DQMServices/Core/interface/MonitorElement.h"
0004 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0005 #include "DQMServices/Core/interface/DQMStore.h"
0006 #include "DataFormats/Math/interface/approx_atan2.h"
0007 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0008 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0009 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0010 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0016 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0017 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0018 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0019 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0020 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0021 
0022 // TODO: change class name to SiPixelCompareRecHitsSoA when CUDA code is removed
0023 class SiPixelCompareRecHits : public DQMEDAnalyzer {
0024 public:
0025   using HitsSoA = reco::TrackingRecHitHost;
0026 
0027   explicit SiPixelCompareRecHits(const edm::ParameterSet&);
0028   ~SiPixelCompareRecHits() override = default;
0029   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0030   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0031   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0032   // analyzeSeparate is templated to accept distinct types of SoAs
0033   // The default use case is to use rechits from Alpaka reconstructed on CPU and GPU;
0034   template <typename U, typename V>
0035   void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event& iEvent);
0036   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0037 
0038 private:
0039   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0040   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoToken_;
0041   // these two are both on Host but originally they have been produced on Host or on Device
0042   const edm::EDGetTokenT<HitsSoA> tokenSoAHitsReference_;
0043   const edm::EDGetTokenT<HitsSoA> tokenSoAHitsTarget_;
0044   const std::string topFolderName_;
0045   const float mind2cut_;
0046   static constexpr uint32_t invalidHit_ = std::numeric_limits<uint32_t>::max();
0047   static constexpr float micron_ = 10000.;
0048   const TrackerGeometry* tkGeom_ = nullptr;
0049   const TrackerTopology* tTopo_ = nullptr;
0050   MonitorElement* hnHits_;
0051   MonitorElement* hBchargeL_[4];  // max 4 barrel hits
0052   MonitorElement* hBsizexL_[4];
0053   MonitorElement* hBsizeyL_[4];
0054   MonitorElement* hBposxL_[4];
0055   MonitorElement* hBposyL_[4];
0056   MonitorElement* hFchargeD_[2][12];  // max 12 endcap disks
0057   MonitorElement* hFsizexD_[2][12];
0058   MonitorElement* hFsizeyD_[2][12];
0059   MonitorElement* hFposxD_[2][12];
0060   MonitorElement* hFposyD_[2][12];
0061   //differences
0062   MonitorElement* hBchargeDiff_;
0063   MonitorElement* hFchargeDiff_;
0064   MonitorElement* hBsizeXDiff_;
0065   MonitorElement* hFsizeXDiff_;
0066   MonitorElement* hBsizeYDiff_;
0067   MonitorElement* hFsizeYDiff_;
0068   MonitorElement* hBposXDiff_;
0069   MonitorElement* hFposXDiff_;
0070   MonitorElement* hBposYDiff_;
0071   MonitorElement* hFposYDiff_;
0072 };
0073 
0074 //
0075 // constructors
0076 //
0077 
0078 SiPixelCompareRecHits::SiPixelCompareRecHits(const edm::ParameterSet& iConfig)
0079     : geomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0080       topoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()),
0081       tokenSoAHitsReference_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsReferenceSoA"))),
0082       tokenSoAHitsTarget_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsTargetSoA"))),
0083       topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
0084       mind2cut_(iConfig.getParameter<double>("minD2cut")) {}
0085 
0086 //
0087 // Begin Run
0088 //
0089 
0090 void SiPixelCompareRecHits::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0091   tkGeom_ = &iSetup.getData(geomToken_);
0092   tTopo_ = &iSetup.getData(topoToken_);
0093 }
0094 
0095 template <typename U, typename V>
0096 void SiPixelCompareRecHits::analyzeSeparate(U tokenRef, V tokenTar, const edm::Event& iEvent) {
0097   const auto& rhsoaHandleRef = iEvent.getHandle(tokenRef);
0098   const auto& rhsoaHandleTar = iEvent.getHandle(tokenTar);
0099 
0100   // Exit early if any handle is invalid
0101   if (!rhsoaHandleRef || !rhsoaHandleTar) {
0102     edm::LogWarning out("SiPixelCompareRecHits");
0103     if (!rhsoaHandleRef)
0104       out << "reference rechits not found; ";
0105     if (!rhsoaHandleTar)
0106       out << "target rechits not found; ";
0107     out << "the comparison will not run.";
0108     return;
0109   }
0110 
0111   const auto& rhsoaRef = *rhsoaHandleRef;
0112   const auto& rhsoaTar = *rhsoaHandleTar;
0113 
0114   auto const& soa2dRef = rhsoaRef.const_view();
0115   auto const& soa2dTar = rhsoaTar.const_view();
0116 
0117   uint32_t nHitsRef = soa2dRef.metadata().size();
0118   uint32_t nHitsTar = soa2dTar.metadata().size();
0119 
0120   hnHits_->Fill(nHitsRef, nHitsTar);
0121 
0122   // Map detector indices to target hits for quick access
0123   std::unordered_map<uint16_t, std::vector<size_t>> detectorIndexMap;
0124   detectorIndexMap.reserve(nHitsTar);
0125   for (size_t j = 0; j < nHitsTar; ++j) {
0126     detectorIndexMap[soa2dTar[j].detectorIndex()].push_back(j);
0127   }
0128 
0129   auto detIds = tkGeom_->detUnitIds();
0130 
0131   // Loop through reference hits
0132   for (uint32_t i = 0; i < nHitsRef; i++) {
0133     float minD = mind2cut_;
0134     uint32_t matchedHit = invalidHit_;
0135     uint16_t indRef = soa2dRef[i].detectorIndex();
0136     float xLocalRef = soa2dRef[i].xLocal();
0137     float yLocalRef = soa2dRef[i].yLocal();
0138 
0139     // Look up hits in target with matching detector index
0140     auto it = detectorIndexMap.find(indRef);
0141     if (it != detectorIndexMap.end()) {
0142       for (auto j : it->second) {
0143         float dx = xLocalRef - soa2dTar[j].xLocal();
0144         float dy = yLocalRef - soa2dTar[j].yLocal();
0145         float distance = dx * dx + dy * dy;
0146         if (distance < minD) {
0147           minD = distance;
0148           matchedHit = j;
0149         }
0150       }
0151     }
0152 
0153     // Gather reference hit properties
0154     DetId id = detIds[indRef];
0155     uint32_t chargeRef = soa2dRef[i].chargeAndStatus().charge;
0156     int16_t sizeXRef = (soa2dRef[i].clusterSizeX() + 7) / 8;
0157     int16_t sizeYRef = (soa2dRef[i].clusterSizeY() + 7) / 8;
0158 
0159     // Initialize target hit properties
0160     uint32_t chargeTar = 0;
0161     int16_t sizeXTar = -99;
0162     int16_t sizeYTar = -99;
0163     float xLocalTar = -999.;
0164     float yLocalTar = -999.;
0165 
0166     if (matchedHit != invalidHit_) {
0167       chargeTar = soa2dTar[matchedHit].chargeAndStatus().charge;
0168       sizeXTar = (soa2dTar[matchedHit].clusterSizeX() + 7) / 8;
0169       sizeYTar = (soa2dTar[matchedHit].clusterSizeY() + 7) / 8;
0170       xLocalTar = soa2dTar[matchedHit].xLocal();
0171       yLocalTar = soa2dTar[matchedHit].yLocal();
0172     }
0173 
0174     // Populate histograms based on subdetector type
0175     switch (id.subdetId()) {
0176       case PixelSubdetector::PixelBarrel:
0177         hBchargeL_[tTopo_->pxbLayer(id) - 1]->Fill(chargeRef, chargeTar);
0178         hBsizexL_[tTopo_->pxbLayer(id) - 1]->Fill(sizeXRef, sizeXTar);
0179         hBsizeyL_[tTopo_->pxbLayer(id) - 1]->Fill(sizeYRef, sizeYTar);
0180         hBposxL_[tTopo_->pxbLayer(id) - 1]->Fill(xLocalRef, xLocalTar);
0181         hBposyL_[tTopo_->pxbLayer(id) - 1]->Fill(yLocalRef, yLocalTar);
0182         hBchargeDiff_->Fill(chargeRef - chargeTar);
0183         hBsizeXDiff_->Fill(sizeXRef - sizeXTar);
0184         hBsizeYDiff_->Fill(sizeYRef - sizeYTar);
0185         hBposXDiff_->Fill(micron_ * (xLocalRef - xLocalTar));
0186         hBposYDiff_->Fill(micron_ * (yLocalRef - yLocalTar));
0187         break;
0188       case PixelSubdetector::PixelEndcap:
0189         hFchargeD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(chargeRef, chargeTar);
0190         hFsizexD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeXRef, sizeXTar);
0191         hFsizeyD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeYRef, sizeYTar);
0192         hFposxD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(xLocalRef, xLocalTar);
0193         hFposyD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(yLocalRef, yLocalTar);
0194         hFchargeDiff_->Fill(chargeRef - chargeTar);
0195         hFsizeXDiff_->Fill(sizeXRef - sizeXTar);
0196         hFsizeYDiff_->Fill(sizeYRef - sizeYTar);
0197         hFposXDiff_->Fill(micron_ * (xLocalRef - xLocalTar));
0198         hFposYDiff_->Fill(micron_ * (yLocalRef - yLocalTar));
0199         break;
0200     }
0201   }
0202 }
0203 
0204 //
0205 // -- Analyze
0206 //
0207 
0208 void SiPixelCompareRecHits::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0209   // The default use case is to use vertices from Alpaka reconstructed on CPU and GPU;
0210   // The function is left templated if any other cases need to be added
0211   analyzeSeparate(tokenSoAHitsReference_, tokenSoAHitsTarget_, iEvent);
0212 }
0213 
0214 //
0215 // -- Book Histograms
0216 //
0217 
0218 void SiPixelCompareRecHits::bookHistograms(DQMStore::IBooker& iBook,
0219                                            edm::Run const& iRun,
0220                                            edm::EventSetup const& iSetup) {
0221   iBook.cd();
0222   iBook.setCurrentFolder(topFolderName_);
0223 
0224   // clang-format off
0225   //Global
0226   hnHits_ = iBook.book2I("nHits", "ReferencevsTarget RecHits per event;#Reference RecHits;#Target RecHits", 200, 0, 5000,200, 0, 5000);
0227   //Barrel Layer
0228   for(unsigned int il=0;il<tkGeom_->numberOfLayers(PixelSubdetector::PixelBarrel);il++){
0229     hBchargeL_[il] = iBook.book2I(Form("recHitsBLay%dCharge",il+1), Form("ReferencevsTarget RecHits Charge Barrel Layer%d;Reference Charge;Target Charge",il+1), 250, 0, 100000, 250, 0, 100000);
0230     hBsizexL_[il] = iBook.book2I(Form("recHitsBLay%dSizex",il+1), Form("ReferencevsTarget RecHits SizeX Barrel Layer%d;Reference SizeX;Target SizeX",il+1), 30, 0, 30, 30, 0, 30);
0231     hBsizeyL_[il] = iBook.book2I(Form("recHitsBLay%dSizey",il+1), Form("ReferencevsTarget RecHits SizeY Barrel Layer%d;Reference SizeY;Target SizeY",il+1), 30, 0, 30, 30, 0, 30);
0232     hBposxL_[il] = iBook.book2D(Form("recHitsBLay%dPosx",il+1), Form("ReferencevsTarget RecHits x-pos in Barrel Layer%d;Reference pos x;Target pos x",il+1), 200, -5, 5, 200,-5,5);
0233     hBposyL_[il] = iBook.book2D(Form("recHitsBLay%dPosy",il+1), Form("ReferencevsTarget RecHits y-pos in Barrel Layer%d;Reference pos y;Target pos y",il+1), 200, -5, 5, 200,-5,5);
0234   }
0235   //Endcaps
0236   //Endcaps Disk
0237   for(int is=0;is<2;is++){
0238     int sign=is==0? -1:1;
0239     for(unsigned int id=0;id<tkGeom_->numberOfLayers(PixelSubdetector::PixelEndcap);id++){
0240       hFchargeD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dCharge",id*sign+sign), Form("ReferencevsTarget RecHits Charge Endcaps Disk%+d;Reference Charge;Target Charge",id*sign+sign), 250, 0, 100000, 250, 0, 100000);
0241       hFsizexD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dSizex",id*sign+sign), Form("ReferencevsTarget RecHits SizeX Endcaps Disk%+d;Reference SizeX;Target SizeX",id*sign+sign), 30, 0, 30, 30, 0, 30);
0242       hFsizeyD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dSizey",id*sign+sign), Form("ReferencevsTarget RecHits SizeY Endcaps Disk%+d;Reference SizeY;Target SizeY",id*sign+sign), 30, 0, 30, 30, 0, 30);
0243       hFposxD_[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosx",id*sign+sign), Form("ReferencevsTarget RecHits x-pos Endcaps Disk%+d;Reference pos x;Target pos x",id*sign+sign), 200, -5, 5, 200, -5, 5);
0244       hFposyD_[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosy",id*sign+sign), Form("ReferencevsTarget RecHits y-pos Endcaps Disk%+d;Reference pos y;Target pos y",id*sign+sign), 200, -5, 5, 200, -5, 5);
0245     }
0246   }
0247   //1D differences
0248   hBchargeDiff_ = iBook.book1D("rechitChargeDiffBpix","Charge differnce of rechits in BPix; rechit charge difference (Reference - Target)", 101, -50.5, 50.5);
0249   hFchargeDiff_ = iBook.book1D("rechitChargeDiffFpix","Charge differnce of rechits in FPix; rechit charge difference (Reference - Target)", 101, -50.5, 50.5);
0250   hBsizeXDiff_ = iBook.book1D("rechitsizeXDiffBpix","SizeX difference of rechits in BPix; rechit sizex difference (Reference - Target)", 21, -10.5, 10.5);
0251   hFsizeXDiff_ = iBook.book1D("rechitsizeXDiffFpix","SizeX difference of rechits in FPix; rechit sizex difference (Reference - Target)", 21, -10.5, 10.5);
0252   hBsizeYDiff_ = iBook.book1D("rechitsizeYDiffBpix","SizeY difference of rechits in BPix; rechit sizey difference (Reference - Target)", 21, -10.5, 10.5);
0253   hFsizeYDiff_ = iBook.book1D("rechitsizeYDiffFpix","SizeY difference of rechits in FPix; rechit sizey difference (Reference - Target)", 21, -10.5, 10.5);
0254   hBposXDiff_ = iBook.book1D("rechitsposXDiffBpix","x-position difference of rechits in BPix; rechit x-pos difference (Reference - Target)", 1000, -10, 10);
0255   hFposXDiff_ = iBook.book1D("rechitsposXDiffFpix","x-position difference of rechits in FPix; rechit x-pos difference (Reference - Target)", 1000, -10, 10);
0256   hBposYDiff_ = iBook.book1D("rechitsposYDiffBpix","y-position difference of rechits in BPix; rechit y-pos difference (Reference - Target)", 1000, -10, 10);
0257   hFposYDiff_ = iBook.book1D("rechitsposYDiffFpix","y-position difference of rechits in FPix; rechit y-pos difference (Reference - Target)", 1000, -10, 10);
0258 }
0259 
0260 
0261 void SiPixelCompareRecHits::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0262   // monitorpixelRecHitsSoAAlpaka
0263   edm::ParameterSetDescription desc;
0264   desc.add<edm::InputTag>("pixelHitsReferenceSoA", edm::InputTag("siPixelRecHitsPreSplittingAlpakaSerial"));
0265   desc.add<edm::InputTag>("pixelHitsTargetSoA", edm::InputTag("siPixelRecHitsPreSplittingAlpaka"));
0266   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelRecHitsCompareDeviceVSHost");
0267   desc.add<double>("minD2cut", 0.0001);
0268   descriptions.addWithDefaultLabel(desc);
0269 }
0270 
0271 using SiPixelCompareRecHits = SiPixelCompareRecHits;
0272 // keeping the old names to allow a smooth HLT migration
0273 using SiPixelPhase1CompareRecHits = SiPixelCompareRecHits;
0274 using SiPixelPhase2CompareRecHits = SiPixelCompareRecHits;
0275 using SiPixelHIonPhase1CompareRecHits = SiPixelCompareRecHits;
0276 
0277 #include "FWCore/Framework/interface/MakerMacros.h"
0278 // TODO: change module names to SiPixel*CompareRecHitsSoA when CUDA code is removed
0279 DEFINE_FWK_MODULE(SiPixelCompareRecHits);
0280 DEFINE_FWK_MODULE(SiPixelPhase1CompareRecHits);
0281 DEFINE_FWK_MODULE(SiPixelPhase2CompareRecHits);
0282 DEFINE_FWK_MODULE(SiPixelHIonPhase1CompareRecHits);
0283