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