File indexing completed on 2025-07-03 00:42:37
0001 #include <memory>
0002 #include <vector>
0003
0004 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0007 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0008 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0009 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/global/EDProducer.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0017 #include "FWCore/Utilities/interface/InputTag.h"
0018 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0019 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0020 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0021 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0022 #include "RecoLocalTracker/SiPixelRecHits/interface/pixelCPEforDevice.h"
0023
0024 class SiPixelRecHitFromSoAAlpaka : public edm::global::EDProducer<> {
0025 using hindex_type = uint32_t;
0026 using HMSstorage = typename std::vector<hindex_type>;
0027
0028 public:
0029 explicit SiPixelRecHitFromSoAAlpaka(const edm::ParameterSet& iConfig);
0030 ~SiPixelRecHitFromSoAAlpaka() override = default;
0031
0032 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0033
0034
0035 using HitsOnHost = ::reco::TrackingRecHitHost;
0036
0037 private:
0038 void produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0039
0040 const uint32_t maxHitsInModules_;
0041 const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0042 const edm::EDGetTokenT<HitsOnHost> hitsToken_;
0043 const edm::EDGetTokenT<SiPixelClusterCollectionNew> clusterToken_;
0044 const edm::EDPutTokenT<SiPixelRecHitCollection> rechitsPutToken_;
0045 const edm::EDPutTokenT<HMSstorage> hostPutToken_;
0046 };
0047
0048 SiPixelRecHitFromSoAAlpaka::SiPixelRecHitFromSoAAlpaka(const edm::ParameterSet& iConfig)
0049 : maxHitsInModules_(iConfig.getParameter<uint32_t>("maxHitsInModules")),
0050 geomToken_(esConsumes()),
0051 hitsToken_(consumes(iConfig.getParameter<edm::InputTag>("pixelRecHitSrc"))),
0052 clusterToken_(consumes(iConfig.getParameter<edm::InputTag>("src"))),
0053 rechitsPutToken_(produces()),
0054 hostPutToken_(produces()) {}
0055
0056 void SiPixelRecHitFromSoAAlpaka::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0057 edm::ParameterSetDescription desc;
0058 desc.add<uint32_t>("maxHitsInModules", phase1PixelTopology::maxNumClustersPerModules)
0059 ->setComment("Max number of hits in a single module");
0060 desc.add<edm::InputTag>("pixelRecHitSrc", edm::InputTag("siPixelRecHitsPreSplittingAlpaka"));
0061 desc.add<edm::InputTag>("src", edm::InputTag("siPixelClustersPreSplitting"));
0062 descriptions.addWithDefaultLabel(desc);
0063 }
0064
0065 void SiPixelRecHitFromSoAAlpaka::produce(edm::StreamID streamID,
0066 edm::Event& iEvent,
0067 const edm::EventSetup& iSetup) const {
0068 auto const& hits = iEvent.get(hitsToken_);
0069 auto hitsView = hits.view();
0070 auto modulesView = hits.view<::reco::HitModuleSoA>();
0071 auto nHits = hitsView.metadata().size();
0072 auto nModules = modulesView.metadata().size();
0073 LogDebug("SiPixelRecHitFromSoAAlpaka") << "converting " << nHits << " hits in max " << nModules << " modules";
0074
0075
0076 SiPixelRecHitCollection output;
0077 output.reserve(nModules, nHits);
0078
0079 HMSstorage hmsp(nModules + 1);
0080
0081 if (0 == nHits) {
0082 hmsp.clear();
0083 iEvent.emplace(rechitsPutToken_, std::move(output));
0084 iEvent.emplace(hostPutToken_, std::move(hmsp));
0085 return;
0086 }
0087
0088
0089 for (unsigned int idx = 0; idx < hmsp.size(); ++idx) {
0090 hmsp[idx] = modulesView.moduleStart()[idx];
0091 }
0092 iEvent.emplace(hostPutToken_, std::move(hmsp));
0093
0094 auto xl = hitsView.xLocal();
0095 auto yl = hitsView.yLocal();
0096 auto xe = hitsView.xerrLocal();
0097 auto ye = hitsView.yerrLocal();
0098
0099 TrackerGeometry const& geom = iSetup.getData(geomToken_);
0100
0101 auto const hclusters = iEvent.getHandle(clusterToken_);
0102
0103 int numberOfDetUnits = 0;
0104 int numberOfClusters = 0;
0105 for (auto const& dsv : *hclusters) {
0106 numberOfDetUnits++;
0107 unsigned int detid = dsv.detId();
0108 DetId detIdObject(detid);
0109 const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
0110 auto gind = genericDet->index();
0111 const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
0112 assert(pixDet);
0113 SiPixelRecHitCollection::FastFiller recHitsOnDetUnit(output, detid);
0114 auto fc = modulesView.moduleStart()[gind];
0115 auto lc = modulesView.moduleStart()[gind + 1];
0116 auto nhits = lc - fc;
0117
0118 assert(lc > fc);
0119 LogDebug("SiPixelRecHitFromSoAAlpaka") << "in det " << gind << ": conv " << nhits << " hits from " << dsv.size()
0120 << " legacy clusters" << ' ' << fc << ',' << lc << "\n";
0121 if (nhits > maxHitsInModules_)
0122 edm::LogWarning("SiPixelRecHitFromSoAAlpaka")
0123 .format("Too many clusters {} in module {}. Only the first {} hits will be converted",
0124 nhits,
0125 gind,
0126 maxHitsInModules_);
0127
0128 nhits = std::min(nhits, maxHitsInModules_);
0129
0130 LogDebug("SiPixelRecHitFromSoAAlpaka") << "in det " << gind << "conv " << nhits << " hits from " << dsv.size()
0131 << " legacy clusters" << ' ' << lc << ',' << fc;
0132
0133 if (0 == nhits)
0134 continue;
0135 auto jnd = [&](int k) { return fc + k; };
0136 assert(nhits <= dsv.size());
0137 if (nhits != dsv.size()) {
0138 edm::LogWarning("GPUHits2CPU") << "nhits!= nclus " << nhits << ' ' << dsv.size();
0139 }
0140 for (auto const& clust : dsv) {
0141 assert(clust.originalId() >= 0);
0142 assert(clust.originalId() < dsv.size());
0143 if (clust.originalId() >= nhits)
0144 continue;
0145 auto ij = jnd(clust.originalId());
0146 LocalPoint lp(xl[ij], yl[ij]);
0147 LocalError le(xe[ij], 0, ye[ij]);
0148 SiPixelRecHitQuality::QualWordType rqw = 0;
0149
0150 numberOfClusters++;
0151
0152
0153
0154
0155
0156
0157
0158
0159
0160 edm::Ref<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster> cluster = edmNew::makeRefTo(hclusters, &clust);
0161
0162 recHitsOnDetUnit.emplace_back(lp, le, rqw, *genericDet, cluster);
0163
0164
0165 LogDebug("SiPixelRecHitFromSoAAlpaka") << "cluster " << numberOfClusters << " at " << lp << ' ' << le;
0166
0167 }
0168
0169
0170 LogDebug("SiPixelRecHitFromSoAAlpaka") << "found " << recHitsOnDetUnit.size() << " RecHits on " << detid;
0171
0172 }
0173
0174 LogDebug("SiPixelRecHitFromSoAAlpaka") << "found " << numberOfDetUnits << " dets, " << numberOfClusters
0175 << " clusters";
0176
0177 iEvent.emplace(rechitsPutToken_, std::move(output));
0178 }
0179
0180 using SiPixelRecHitFromSoAAlpakaPhase1 = SiPixelRecHitFromSoAAlpaka;
0181 using SiPixelRecHitFromSoAAlpakaPhase2 = SiPixelRecHitFromSoAAlpaka;
0182 using SiPixelRecHitFromSoAAlpakaHIonPhase1 = SiPixelRecHitFromSoAAlpaka;
0183
0184 #include "FWCore/Framework/interface/MakerMacros.h"
0185 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpaka);
0186
0187 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpakaPhase1);
0188 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpakaPhase2);
0189 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpakaHIonPhase1);