Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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 "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0015 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0016 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0017 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0018 
0019 class SiPixelMonitorRecHitsSoAAlpaka : public DQMEDAnalyzer {
0020 public:
0021   using HitsOnHost = reco::TrackingRecHitHost;
0022 
0023   explicit SiPixelMonitorRecHitsSoAAlpaka(const edm::ParameterSet&);
0024   ~SiPixelMonitorRecHitsSoAAlpaka() override = default;
0025   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0026   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0027   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0028   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029 
0030 private:
0031   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0032   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> topoToken_;
0033   const edm::EDGetTokenT<HitsOnHost> tokenSoAHits_;
0034   const std::string topFolderName_;
0035   const TrackerGeometry* tkGeom_ = nullptr;
0036   const TrackerTopology* tTopo_ = nullptr;
0037   MonitorElement* hnHits;
0038   MonitorElement* hBFposZP;
0039   MonitorElement* hBFposZR;
0040   MonitorElement* hBposXY;
0041   MonitorElement* hBposZP;
0042   MonitorElement* hBcharge;
0043   MonitorElement* hBsizex;
0044   MonitorElement* hBsizey;
0045   MonitorElement* hBposZPL[4];  // max 4 barrel hits
0046   MonitorElement* hBchargeL[4];
0047   MonitorElement* hBsizexL[4];
0048   MonitorElement* hBsizeyL[4];
0049   MonitorElement* hFposXY;
0050   MonitorElement* hFposZP;
0051   MonitorElement* hFcharge;
0052   MonitorElement* hFsizex;
0053   MonitorElement* hFsizey;
0054   MonitorElement* hFposXYD[2][12];  // max 12 endcap disks
0055   MonitorElement* hFchargeD[2][12];
0056   MonitorElement* hFsizexD[2][12];
0057   MonitorElement* hFsizeyD[2][12];
0058 };
0059 
0060 //
0061 // constructors
0062 //
0063 
0064 SiPixelMonitorRecHitsSoAAlpaka::SiPixelMonitorRecHitsSoAAlpaka(const edm::ParameterSet& iConfig)
0065     : geomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0066       topoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>()),
0067       tokenSoAHits_(consumes(iConfig.getParameter<edm::InputTag>("pixelHitsSrc"))),
0068       topFolderName_(iConfig.getParameter<std::string>("TopFolderName")) {}
0069 
0070 //
0071 // Begin Run
0072 //
0073 
0074 void SiPixelMonitorRecHitsSoAAlpaka::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0075   tkGeom_ = &iSetup.getData(geomToken_);
0076   tTopo_ = &iSetup.getData(topoToken_);
0077 }
0078 
0079 //
0080 // -- Analyze
0081 //
0082 
0083 void SiPixelMonitorRecHitsSoAAlpaka::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0084   const auto& rhsoaHandle = iEvent.getHandle(tokenSoAHits_);
0085   if (!rhsoaHandle.isValid()) {
0086     edm::LogWarning("SiPixelMonitorRecHitsSoAAlpaka") << "No RecHits SoA found \n returning!";
0087     return;
0088   }
0089   auto const& rhsoa = *rhsoaHandle;
0090   auto const& soa2d = rhsoa.const_view();
0091 
0092   uint32_t nHits_ = soa2d.metadata().size();
0093   hnHits->Fill(nHits_);
0094   auto detIds = tkGeom_->detUnitIds();
0095   for (uint32_t i = 0; i < nHits_; i++) {
0096     DetId id = detIds[soa2d[i].detectorIndex()];
0097     float xG = soa2d[i].xGlobal();
0098     float yG = soa2d[i].yGlobal();
0099     float zG = soa2d[i].zGlobal();
0100     float rG = soa2d[i].rGlobal();
0101     float fphi = short2phi(soa2d[i].iphi());
0102     uint32_t charge = soa2d[i].chargeAndStatus().charge;
0103     int16_t sizeX = std::ceil(float(std::abs(soa2d[i].clusterSizeX()) / 8.));
0104     int16_t sizeY = std::ceil(float(std::abs(soa2d[i].clusterSizeY()) / 8.));
0105     hBFposZP->Fill(zG, fphi);
0106     int16_t ysign = yG >= 0 ? 1 : -1;
0107     hBFposZR->Fill(zG, rG * ysign);
0108     switch (id.subdetId()) {
0109       case PixelSubdetector::PixelBarrel:
0110         hBposXY->Fill(xG, yG);
0111         hBposZP->Fill(zG, fphi);
0112         hBcharge->Fill(charge);
0113         hBsizex->Fill(sizeX);
0114         hBsizey->Fill(sizeY);
0115         hBposZPL[tTopo_->pxbLayer(id) - 1]->Fill(zG, fphi);
0116         hBchargeL[tTopo_->pxbLayer(id) - 1]->Fill(charge);
0117         hBsizexL[tTopo_->pxbLayer(id) - 1]->Fill(sizeX);
0118         hBsizeyL[tTopo_->pxbLayer(id) - 1]->Fill(sizeY);
0119         break;
0120       case PixelSubdetector::PixelEndcap:
0121         hFposXY->Fill(xG, yG);
0122         hFposZP->Fill(zG, fphi);
0123         hFcharge->Fill(charge);
0124         hFsizex->Fill(sizeX);
0125         hFsizey->Fill(sizeY);
0126         hFposXYD[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(xG, yG);
0127         hFchargeD[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(charge);
0128         hFsizexD[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeX);
0129         hFsizeyD[tTopo_->pxfSide(id) - 1][tTopo_->pxfDisk(id) - 1]->Fill(sizeY);
0130         break;
0131     }
0132   }
0133 }
0134 
0135 //
0136 // -- Book Histograms
0137 //
0138 
0139 void SiPixelMonitorRecHitsSoAAlpaka::bookHistograms(DQMStore::IBooker& iBook,
0140                                                     edm::Run const& iRun,
0141                                                     edm::EventSetup const& iSetup) {
0142   iBook.cd();
0143   iBook.setCurrentFolder(topFolderName_);
0144 
0145   // clang-format off
0146   //Global
0147   hnHits = iBook.book1D("nHits", "RecHits per event;RecHits;#events", 200, 0, 5000);
0148   hBFposZP = iBook.book2D("recHitsGlobalPosZP", "RecHits position Global;Z;#phi", 1000, -60, 60, 200,-3.2,3.2);
0149   hBFposZR = iBook.book2D("recHitsGlobalPosZR", "RecHits position Global;Z;R", 1000, -60, 60, 200,-20,20);
0150   //Barrel
0151   hBposXY = iBook.book2D("recHitsBarrelPosXY", "RecHits position Barrel;X;Y", 200, -20, 20, 200,-20,20);
0152   hBposZP = iBook.book2D("recHitsBarrelPosZP", "RecHits position Barrel;Z;#phi", 300, -30, 30, 200,-3.2,3.2);
0153   hBcharge = iBook.book1D("recHitsBarrelCharge", "RecHits Charge Barrel;Charge;#events", 250, 0, 100000);
0154   hBsizex = iBook.book1D("recHitsBarrelSizex", "RecHits SizeX Barrel;SizeX;#events", 50, 0, 50);
0155   hBsizey = iBook.book1D("recHitsBarrelSizey", "RecHits SizeY Barrel;SizeY;#events", 50, 0, 50);
0156   //Barrel Layer
0157   for(unsigned int il=0;il<tkGeom_->numberOfLayers(PixelSubdetector::PixelBarrel);il++){
0158     hBposZPL[il] = iBook.book2D(Form("recHitsBLay%dPosZP",il+1), Form("RecHits position Barrel Layer%d;Z;#phi",il+1), 300, -30, 30, 200,-3.2,3.2);
0159     hBchargeL[il] = iBook.book1D(Form("recHitsBLay%dCharge",il+1), Form("RecHits Charge Barrel Layer%d;Charge;#events",il+1), 250, 0, 100000);
0160     hBsizexL[il] = iBook.book1D(Form("recHitsBLay%dSizex",il+1), Form("RecHits SizeX Barrel Layer%d;SizeX;#events",il+1), 50, 0, 50);
0161     hBsizeyL[il] = iBook.book1D(Form("recHitsBLay%dSizey",il+1), Form("RecHits SizeY Barrel Layer%d;SizeY;#events",il+1), 50, 0, 50);
0162   }
0163   //Endcaps
0164   hFposXY = iBook.book2D("recHitsEndcapsPosXY", "RecHits position Endcaps;X;Y", 200, -20, 20, 200,-20, 20);
0165   hFposZP = iBook.book2D("recHitsEndcapsPosZP", "RecHits position Endcaps;Z;#phi", 600, -60, 60, 200,-3.2,3.2);
0166   hFcharge = iBook.book1D("recHitsEndcapsCharge", "RecHits Charge Endcaps;Charge;#events", 250, 0, 100000);
0167   hFsizex = iBook.book1D("recHitsEndcapsSizex", "RecHits SizeX Endcaps;SizeX;#events", 50, 0, 50);
0168   hFsizey = iBook.book1D("recHitsEndcapsSizey", "RecHits SizeY Endcaps;SizeY;#events", 50, 0, 50);
0169   //Endcaps Disk
0170   for(int is=0;is<2;is++){
0171     int sign=is==0? -1:1;
0172     for(unsigned int id=0;id<tkGeom_->numberOfLayers(PixelSubdetector::PixelEndcap);id++){
0173       hFposXYD[is][id] = iBook.book2D(Form("recHitsFDisk%+dPosXY",id*sign+sign), Form("RecHits position Endcaps Disk%+d;X;Y",id*sign+sign), 200, -20, 20, 200,-20,20);
0174       hFchargeD[is][id] = iBook.book1D(Form("recHitsFDisk%+dCharge",id*sign+sign), Form("RecHits Charge Endcaps Disk%+d;Charge;#events",id*sign+sign), 250, 0, 100000);
0175       hFsizexD[is][id] = iBook.book1D(Form("recHitsFDisk%+dSizex",id*sign+sign), Form("RecHits SizeX Endcaps Disk%+d;SizeX;#events",id*sign+sign), 50, 0, 50);
0176       hFsizeyD[is][id] = iBook.book1D(Form("recHitsFDisk%+dSizey",id*sign+sign), Form("RecHits SizeY Endcaps Disk%+d;SizeY;#events",id*sign+sign), 50, 0, 50);
0177     }
0178   }
0179 }
0180 
0181 void SiPixelMonitorRecHitsSoAAlpaka::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0182   // monitorpixelRecHitsSoA
0183   edm::ParameterSetDescription desc;
0184   desc.add<edm::InputTag>("pixelHitsSrc", edm::InputTag("siPixelRecHitsPreSplittingAlpaka"));
0185   desc.add<std::string>("TopFolderName", "SiPixelHeterogeneous/PixelRecHitsAlpaka");
0186   descriptions.addWithDefaultLabel(desc);
0187 }
0188 
0189 using SiPixelPhase1MonitorRecHitsSoAAlpaka = SiPixelMonitorRecHitsSoAAlpaka;
0190 using SiPixelPhase2MonitorRecHitsSoAAlpaka = SiPixelMonitorRecHitsSoAAlpaka;
0191 using SiPixelHIonPhase1MonitorRecHitsSoAAlpaka = SiPixelMonitorRecHitsSoAAlpaka;
0192 
0193 #include "FWCore/Framework/interface/MakerMacros.h"
0194 DEFINE_FWK_MODULE(SiPixelMonitorRecHitsSoAAlpaka);
0195 DEFINE_FWK_MODULE(SiPixelPhase1MonitorRecHitsSoAAlpaka);
0196 DEFINE_FWK_MODULE(SiPixelPhase2MonitorRecHitsSoAAlpaka);
0197 DEFINE_FWK_MODULE(SiPixelHIonPhase1MonitorRecHitsSoAAlpaka);
0198