Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQMServices/Core/interface/MonitorElement.h"
0002 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0003 #include "DQMServices/Core/interface/DQMStore.h"
0004 #include "DataFormats/Math/interface/approx_atan2.h"
0005 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0006 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0007 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0008 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0016 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0017 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0018 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0019 
0020 template <typename T>
0021 class SiPixelCompareRecHitsSoAAlpaka : public DQMEDAnalyzer {
0022 public:
0023   using HitsOnHost = TrackingRecHitHost<T>;
0024 
0025   explicit SiPixelCompareRecHitsSoAAlpaka(const edm::ParameterSet&);
0026   ~SiPixelCompareRecHitsSoAAlpaka() override = default;
0027   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0028   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0029   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0030   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0031 
0032 private:
0033   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0034   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoToken_;
0035   const edm::EDGetTokenT<HitsOnHost> tokenSoAHitsHost_;    //these two are both on Host but originally they have been
0036   const edm::EDGetTokenT<HitsOnHost> tokenSoAHitsDevice_;  //produced on Host or on Device
0037   const std::string topFolderName_;
0038   const float mind2cut_;
0039   static constexpr uint32_t invalidHit_ = std::numeric_limits<uint32_t>::max();
0040   static constexpr float micron_ = 10000.;
0041   const TrackerGeometry* tkGeom_ = nullptr;
0042   const TrackerTopology* tTopo_ = nullptr;
0043   MonitorElement* hnHits_;
0044   MonitorElement* hBchargeL_[4];  // max 4 barrel hits
0045   MonitorElement* hBsizexL_[4];
0046   MonitorElement* hBsizeyL_[4];
0047   MonitorElement* hBposxL_[4];
0048   MonitorElement* hBposyL_[4];
0049   MonitorElement* hFchargeD_[2][12];  // max 12 endcap disks
0050   MonitorElement* hFsizexD_[2][12];
0051   MonitorElement* hFsizeyD_[2][12];
0052   MonitorElement* hFposxD_[2][12];
0053   MonitorElement* hFposyD_[2][12];
0054   //differences
0055   MonitorElement* hBchargeDiff_;
0056   MonitorElement* hFchargeDiff_;
0057   MonitorElement* hBsizeXDiff_;
0058   MonitorElement* hFsizeXDiff_;
0059   MonitorElement* hBsizeYDiff_;
0060   MonitorElement* hFsizeYDiff_;
0061   MonitorElement* hBposXDiff_;
0062   MonitorElement* hFposXDiff_;
0063   MonitorElement* hBposYDiff_;
0064   MonitorElement* hFposYDiff_;
0065 };
0066 
0067 //
0068 // constructors
0069 //
0070 template <typename T>
0071 SiPixelCompareRecHitsSoAAlpaka<T>::SiPixelCompareRecHitsSoAAlpaka(const edm::ParameterSet& iConfig)
0072     : geomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0073       topoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()),
0074       tokenSoAHitsHost_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsSrcHost"))),
0075       tokenSoAHitsDevice_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsSrcDevice"))),
0076       topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
0077       mind2cut_(iConfig.getParameter<double>("minD2cut")) {}
0078 
0079 //
0080 // Begin Run
0081 //
0082 template <typename T>
0083 void SiPixelCompareRecHitsSoAAlpaka<T>::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0084   tkGeom_ = &iSetup.getData(geomToken_);
0085   tTopo_ = &iSetup.getData(topoToken_);
0086 }
0087 
0088 //
0089 // -- Analyze
0090 //
0091 template <typename T>
0092 void SiPixelCompareRecHitsSoAAlpaka<T>::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0093   const auto& rhsoaHandleHost = iEvent.getHandle(tokenSoAHitsHost_);
0094   const auto& rhsoaHandleDevice = iEvent.getHandle(tokenSoAHitsDevice_);
0095   if (not rhsoaHandleHost or not rhsoaHandleDevice) {
0096     edm::LogWarning out("SiPixelCompareRecHitsSoAAlpaka");
0097     if (not rhsoaHandleHost) {
0098       out << "reference (Host) rechits not found; ";
0099     }
0100     if (not rhsoaHandleDevice) {
0101       out << "target (Device) rechits not found; ";
0102     }
0103     out << "the comparison will not run.";
0104     return;
0105   }
0106 
0107   auto const& rhsoaHost = *rhsoaHandleHost;
0108   auto const& rhsoaDevice = *rhsoaHandleDevice;
0109 
0110   auto const& soa2dHost = rhsoaHost.const_view();
0111   auto const& soa2dDevice = rhsoaDevice.const_view();
0112 
0113   uint32_t nHitsHost = soa2dHost.metadata().size();
0114   uint32_t nHitsDevice = soa2dDevice.metadata().size();
0115 
0116   hnHits_->Fill(nHitsHost, nHitsDevice);
0117   auto detIds = tkGeom_->detUnitIds();
0118   for (uint32_t i = 0; i < nHitsHost; i++) {
0119     float minD = mind2cut_;
0120     uint32_t matchedHit = invalidHit_;
0121     uint16_t indHost = soa2dHost[i].detectorIndex();
0122     float xLocalHost = soa2dHost[i].xLocal();
0123     float yLocalHost = soa2dHost[i].yLocal();
0124     for (uint32_t j = 0; j < nHitsDevice; j++) {
0125       if (soa2dDevice.detectorIndex(j) == indHost) {
0126         float dx = xLocalHost - soa2dDevice[j].xLocal();
0127         float dy = yLocalHost - soa2dDevice[j].yLocal();
0128         float distance = dx * dx + dy * dy;
0129         if (distance < minD) {
0130           minD = distance;
0131           matchedHit = j;
0132         }
0133       }
0134     }
0135     DetId id = detIds[indHost];
0136     uint32_t chargeHost = soa2dHost[i].chargeAndStatus().charge;
0137     int16_t sizeXHost = std::ceil(float(std::abs(soa2dHost[i].clusterSizeX()) / 8.));
0138     int16_t sizeYHost = std::ceil(float(std::abs(soa2dHost[i].clusterSizeY()) / 8.));
0139     uint32_t chargeDevice = 0;
0140     int16_t sizeXDevice = -99;
0141     int16_t sizeYDevice = -99;
0142     float xLocalDevice = -999.;
0143     float yLocalDevice = -999.;
0144     if (matchedHit != invalidHit_) {
0145       chargeDevice = soa2dDevice[matchedHit].chargeAndStatus().charge;
0146       sizeXDevice = std::ceil(float(std::abs(soa2dDevice[matchedHit].clusterSizeX()) / 8.));
0147       sizeYDevice = std::ceil(float(std::abs(soa2dDevice[matchedHit].clusterSizeY()) / 8.));
0148       xLocalDevice = soa2dDevice[matchedHit].xLocal();
0149       yLocalDevice = soa2dDevice[matchedHit].yLocal();
0150     }
0151     switch (id.subdetId()) {
0152       case PixelSubdetector::PixelBarrel:
0153         hBchargeL_[tTopo_->pxbLayer(id) - 1]->Fill(chargeHost, chargeDevice);
0154         hBsizexL_[tTopo_->pxbLayer(id) - 1]->Fill(sizeXHost, sizeXDevice);
0155         hBsizeyL_[tTopo_->pxbLayer(id) - 1]->Fill(sizeYHost, sizeYDevice);
0156         hBposxL_[tTopo_->pxbLayer(id) - 1]->Fill(xLocalHost, xLocalDevice);
0157         hBposyL_[tTopo_->pxbLayer(id) - 1]->Fill(yLocalHost, yLocalDevice);
0158         hBchargeDiff_->Fill(chargeHost - chargeDevice);
0159         hBsizeXDiff_->Fill(sizeXHost - sizeXDevice);
0160         hBsizeYDiff_->Fill(sizeYHost - sizeYDevice);
0161         hBposXDiff_->Fill(micron_ * (xLocalHost - xLocalDevice));
0162         hBposYDiff_->Fill(micron_ * (yLocalHost - yLocalDevice));
0163         break;
0164       case PixelSubdetector::PixelEndcap:
0165         hFchargeD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(chargeHost, chargeDevice);
0166         hFsizexD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeXHost, sizeXDevice);
0167         hFsizeyD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeYHost, sizeYDevice);
0168         hFposxD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(xLocalHost, xLocalDevice);
0169         hFposyD_[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(yLocalHost, yLocalDevice);
0170         hFchargeDiff_->Fill(chargeHost - chargeDevice);
0171         hFsizeXDiff_->Fill(sizeXHost - sizeXDevice);
0172         hFsizeYDiff_->Fill(sizeYHost - sizeYDevice);
0173         hFposXDiff_->Fill(micron_ * (xLocalHost - xLocalDevice));
0174         hFposYDiff_->Fill(micron_ * (yLocalHost - yLocalDevice));
0175         break;
0176     }
0177   }
0178 }
0179 
0180 //
0181 // -- Book Histograms
0182 //
0183 template <typename T>
0184 void SiPixelCompareRecHitsSoAAlpaka<T>::bookHistograms(DQMStore::IBooker& iBook,
0185                                                        edm::Run const& iRun,
0186                                                        edm::EventSetup const& iSetup) {
0187   iBook.cd();
0188   iBook.setCurrentFolder(topFolderName_);
0189 
0190   // clang-format off
0191   //Global
0192   hnHits_ = iBook.book2I("nHits", "HostvsDevice RecHits per event;#Host RecHits;#Device RecHits", 200, 0, 5000,200, 0, 5000);
0193   //Barrel Layer
0194   for(unsigned int il=0;il<tkGeom_->numberOfLayers(PixelSubdetector::PixelBarrel);il++){
0195     hBchargeL_[il] = iBook.book2I(Form("recHitsBLay%dCharge",il+1), Form("HostvsDevice RecHits Charge Barrel Layer%d;Host Charge;Device Charge",il+1), 250, 0, 100000, 250, 0, 100000);
0196     hBsizexL_[il] = iBook.book2I(Form("recHitsBLay%dSizex",il+1), Form("HostvsDevice RecHits SizeX Barrel Layer%d;Host SizeX;Device SizeX",il+1), 30, 0, 30, 30, 0, 30);
0197     hBsizeyL_[il] = iBook.book2I(Form("recHitsBLay%dSizey",il+1), Form("HostvsDevice RecHits SizeY Barrel Layer%d;Host SizeY;Device SizeY",il+1), 30, 0, 30, 30, 0, 30);
0198     hBposxL_[il] = iBook.book2D(Form("recHitsBLay%dPosx",il+1), Form("HostvsDevice RecHits x-pos in Barrel Layer%d;Host pos x;Device pos x",il+1), 200, -5, 5, 200,-5,5);
0199     hBposyL_[il] = iBook.book2D(Form("recHitsBLay%dPosy",il+1), Form("HostvsDevice RecHits y-pos in Barrel Layer%d;Host pos y;Device pos y",il+1), 200, -5, 5, 200,-5,5);
0200   }
0201   //Endcaps
0202   //Endcaps Disk
0203   for(int is=0;is<2;is++){
0204     int sign=is==0? -1:1;
0205     for(unsigned int id=0;id<tkGeom_->numberOfLayers(PixelSubdetector::PixelEndcap);id++){
0206       hFchargeD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dCharge",id*sign+sign), Form("HostvsDevice RecHits Charge Endcaps Disk%+d;Host Charge;Device Charge",id*sign+sign), 250, 0, 100000, 250, 0, 100000);
0207       hFsizexD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dSizex",id*sign+sign), Form("HostvsDevice RecHits SizeX Endcaps Disk%+d;Host SizeX;Device SizeX",id*sign+sign), 30, 0, 30, 30, 0, 30);
0208       hFsizeyD_[is][id] = iBook.book2I(Form("recHitsFDisk%+dSizey",id*sign+sign), Form("HostvsDevice RecHits SizeY Endcaps Disk%+d;Host SizeY;Device SizeY",id*sign+sign), 30, 0, 30, 30, 0, 30);
0209       hFposxD_[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosx",id*sign+sign), Form("HostvsDevice RecHits x-pos Endcaps Disk%+d;Host pos x;Device pos x",id*sign+sign), 200, -5, 5, 200, -5, 5);
0210       hFposyD_[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosy",id*sign+sign), Form("HostvsDevice RecHits y-pos Endcaps Disk%+d;Host pos y;Device pos y",id*sign+sign), 200, -5, 5, 200, -5, 5);
0211     }
0212   }
0213   //1D differences
0214   hBchargeDiff_ = iBook.book1D("rechitChargeDiffBpix","Charge differnce of rechits in BPix; rechit charge difference (Host - Device)", 101, -50.5, 50.5);
0215   hFchargeDiff_ = iBook.book1D("rechitChargeDiffFpix","Charge differnce of rechits in FPix; rechit charge difference (Host - Device)", 101, -50.5, 50.5);
0216   hBsizeXDiff_ = iBook.book1D("rechitsizeXDiffBpix","SizeX difference of rechits in BPix; rechit sizex difference (Host - Device)", 21, -10.5, 10.5);
0217   hFsizeXDiff_ = iBook.book1D("rechitsizeXDiffFpix","SizeX difference of rechits in FPix; rechit sizex difference (Host - Device)", 21, -10.5, 10.5);
0218   hBsizeYDiff_ = iBook.book1D("rechitsizeYDiffBpix","SizeY difference of rechits in BPix; rechit sizey difference (Host - Device)", 21, -10.5, 10.5);
0219   hFsizeYDiff_ = iBook.book1D("rechitsizeYDiffFpix","SizeY difference of rechits in FPix; rechit sizey difference (Host - Device)", 21, -10.5, 10.5);
0220   hBposXDiff_ = iBook.book1D("rechitsposXDiffBpix","x-position difference of rechits in BPix; rechit x-pos difference (Host - Device)", 1000, -10, 10);
0221   hFposXDiff_ = iBook.book1D("rechitsposXDiffFpix","x-position difference of rechits in FPix; rechit x-pos difference (Host - Device)", 1000, -10, 10);
0222   hBposYDiff_ = iBook.book1D("rechitsposYDiffBpix","y-position difference of rechits in BPix; rechit y-pos difference (Host - Device)", 1000, -10, 10);
0223   hFposYDiff_ = iBook.book1D("rechitsposYDiffFpix","y-position difference of rechits in FPix; rechit y-pos difference (Host - Device)", 1000, -10, 10);
0224 }
0225 
0226 template<typename T>
0227 void SiPixelCompareRecHitsSoAAlpaka<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0228   // monitorpixelRecHitsSoAAlpaka
0229   edm::ParameterSetDescription desc;
0230   desc.add<edm::InputTag>("pixelHitsSrcHost", edm::InputTag("siPixelRecHitsPreSplittingAlpakaSerial"));
0231   desc.add<edm::InputTag>("pixelHitsSrcDevice", edm::InputTag("siPixelRecHitsPreSplittingAlpaka"));
0232   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelRecHitsCompareDeviceVSHost");
0233   desc.add<double>("minD2cut", 0.0001);
0234   descriptions.addWithDefaultLabel(desc);
0235 }
0236 
0237 using SiPixelPhase1CompareRecHitsSoAAlpaka = SiPixelCompareRecHitsSoAAlpaka<pixelTopology::Phase1>;
0238 using SiPixelPhase2CompareRecHitsSoAAlpaka = SiPixelCompareRecHitsSoAAlpaka<pixelTopology::Phase2>;
0239 using SiPixelHIonPhase1CompareRecHitsSoAAlpaka = SiPixelCompareRecHitsSoAAlpaka<pixelTopology::HIonPhase1>;
0240 
0241 #include "FWCore/Framework/interface/MakerMacros.h"
0242 DEFINE_FWK_MODULE(SiPixelPhase1CompareRecHitsSoAAlpaka);
0243 DEFINE_FWK_MODULE(SiPixelPhase2CompareRecHitsSoAAlpaka);
0244 DEFINE_FWK_MODULE(SiPixelHIonPhase1CompareRecHitsSoAAlpaka);