Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:23

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 template <typename TrackerTraits>
0025 class SiPixelRecHitFromSoAAlpaka : public edm::global::EDProducer<> {
0026   using HitModuleStartArray = typename TrackingRecHitSoA<TrackerTraits>::HitModuleStartArray;
0027   using hindex_type = typename TrackerTraits::hindex_type;
0028   using HMSstorage = typename std::vector<hindex_type>;
0029 
0030 public:
0031   explicit SiPixelRecHitFromSoAAlpaka(const edm::ParameterSet& iConfig);
0032   ~SiPixelRecHitFromSoAAlpaka() override = default;
0033 
0034   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035 
0036   // Data has been implicitly copied from Device to Host by the framework
0037   using HitsOnHost = TrackingRecHitHost<TrackerTraits>;
0038 
0039 private:
0040   void produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
0041 
0042   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0043   const edm::EDGetTokenT<HitsOnHost> hitsToken_;                      // Alpaka hits
0044   const edm::EDGetTokenT<SiPixelClusterCollectionNew> clusterToken_;  // legacy clusters
0045   const edm::EDPutTokenT<SiPixelRecHitCollection> rechitsPutToken_;   // legacy rechits
0046   const edm::EDPutTokenT<HMSstorage> hostPutToken_;
0047 };
0048 
0049 template <typename TrackerTraits>
0050 SiPixelRecHitFromSoAAlpaka<TrackerTraits>::SiPixelRecHitFromSoAAlpaka(const edm::ParameterSet& iConfig)
0051     : geomToken_(esConsumes()),
0052       hitsToken_(consumes(iConfig.getParameter<edm::InputTag>("pixelRecHitSrc"))),
0053       clusterToken_(consumes(iConfig.getParameter<edm::InputTag>("src"))),
0054       rechitsPutToken_(produces()),
0055       hostPutToken_(produces()) {}
0056 
0057 template <typename TrackerTraits>
0058 void SiPixelRecHitFromSoAAlpaka<TrackerTraits>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0059   edm::ParameterSetDescription desc;
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 template <typename TrackerTraits>
0066 void SiPixelRecHitFromSoAAlpaka<TrackerTraits>::produce(edm::StreamID streamID,
0067                                                         edm::Event& iEvent,
0068                                                         const edm::EventSetup& iSetup) const {
0069   auto const& hits = iEvent.get(hitsToken_);
0070   auto nHits = hits.view().metadata().size();
0071   LogDebug("SiPixelRecHitFromSoAAlpaka") << "converting " << nHits << " Hits";
0072 
0073   // allocate a buffer for the indices of the clusters
0074   constexpr auto nMaxModules = TrackerTraits::numberOfModules;
0075 
0076   SiPixelRecHitCollection output;
0077   output.reserve(nMaxModules, nHits);
0078 
0079   HMSstorage hmsp(nMaxModules + 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   // fill content of HMSstorage product, and put it into the Event
0089   for (unsigned int idx = 0; idx < hmsp.size(); ++idx) {
0090     hmsp[idx] = hits.view().hitsModuleStart()[idx];
0091   }
0092   iEvent.emplace(hostPutToken_, std::move(hmsp));
0093 
0094   auto xl = hits.view().xLocal();
0095   auto yl = hits.view().yLocal();
0096   auto xe = hits.view().xerrLocal();
0097   auto ye = hits.view().yerrLocal();
0098 
0099   TrackerGeometry const& geom = iSetup.getData(geomToken_);
0100 
0101   auto const hclusters = iEvent.getHandle(clusterToken_);
0102 
0103   constexpr uint32_t maxHitsInModule = pixelClustering::maxHitsInModule();
0104 
0105   int numberOfDetUnits = 0;
0106   int numberOfClusters = 0;
0107   for (auto const& dsv : *hclusters) {
0108     numberOfDetUnits++;
0109     unsigned int detid = dsv.detId();
0110     DetId detIdObject(detid);
0111     const GeomDetUnit* genericDet = geom.idToDetUnit(detIdObject);
0112     auto gind = genericDet->index();
0113     const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
0114     assert(pixDet);
0115     SiPixelRecHitCollection::FastFiller recHitsOnDetUnit(output, detid);
0116     auto fc = hits.view().hitsModuleStart()[gind];
0117     auto lc = hits.view().hitsModuleStart()[gind + 1];
0118     auto nhits = lc - fc;
0119 
0120     assert(lc > fc);
0121     LogDebug("SiPixelRecHitFromSoAAlpaka") << "in det " << gind << ": conv " << nhits << " hits from " << dsv.size()
0122                                            << " legacy clusters" << ' ' << fc << ',' << lc << "\n";
0123     if (nhits > maxHitsInModule)
0124       edm::LogWarning("SiPixelRecHitFromSoAAlpaka")
0125           .format("Too many clusters {} in module {}. Only the first {} hits will be converted",
0126                   nhits,
0127                   gind,
0128                   maxHitsInModule);
0129 
0130     nhits = std::min(nhits, maxHitsInModule);
0131 
0132     LogDebug("SiPixelRecHitFromSoAAlpaka") << "in det " << gind << "conv " << nhits << " hits from " << dsv.size()
0133                                            << " legacy clusters" << ' ' << lc << ',' << fc;
0134 
0135     if (0 == nhits)
0136       continue;
0137     auto jnd = [&](int k) { return fc + k; };
0138     assert(nhits <= dsv.size());
0139     if (nhits != dsv.size()) {
0140       edm::LogWarning("GPUHits2CPU") << "nhits!= nclus " << nhits << ' ' << dsv.size();
0141     }
0142     for (auto const& clust : dsv) {
0143       assert(clust.originalId() >= 0);
0144       assert(clust.originalId() < dsv.size());
0145       if (clust.originalId() >= nhits)
0146         continue;
0147       auto ij = jnd(clust.originalId());
0148       LocalPoint lp(xl[ij], yl[ij]);
0149       LocalError le(xe[ij], 0, ye[ij]);
0150       SiPixelRecHitQuality::QualWordType rqw = 0;
0151 
0152       numberOfClusters++;
0153 
0154       /* cpu version....  (for reference)
0155     std::tuple<LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType> tuple = cpe_->getParameters( clust, *genericDet );
0156     LocalPoint lp( std::get<0>(tuple) );
0157     LocalError le( std::get<1>(tuple) );
0158     SiPixelRecHitQuality::QualWordType rqw( std::get<2>(tuple) );
0159     */
0160 
0161       // Create a persistent edm::Ref to the cluster
0162       edm::Ref<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster> cluster = edmNew::makeRefTo(hclusters, &clust);
0163       // Make a RecHit and add it to the DetSet
0164       recHitsOnDetUnit.emplace_back(lp, le, rqw, *genericDet, cluster);
0165       // =============================
0166 
0167       LogDebug("SiPixelRecHitFromSoAAlpaka") << "cluster " << numberOfClusters << " at " << lp << ' ' << le;
0168 
0169     }  //  <-- End loop on Clusters
0170 
0171     //  LogDebug("SiPixelRecHitGPU")
0172     LogDebug("SiPixelRecHitFromSoAAlpaka") << "found " << recHitsOnDetUnit.size() << " RecHits on " << detid;
0173 
0174   }  //    <-- End loop on DetUnits
0175 
0176   LogDebug("SiPixelRecHitFromSoAAlpaka") << "found " << numberOfDetUnits << " dets, " << numberOfClusters
0177                                          << " clusters";
0178 
0179   iEvent.emplace(rechitsPutToken_, std::move(output));
0180 }
0181 
0182 using SiPixelRecHitFromSoAAlpakaPhase1 = SiPixelRecHitFromSoAAlpaka<pixelTopology::Phase1>;
0183 using SiPixelRecHitFromSoAAlpakaPhase2 = SiPixelRecHitFromSoAAlpaka<pixelTopology::Phase2>;
0184 using SiPixelRecHitFromSoAAlpakaHIonPhase1 = SiPixelRecHitFromSoAAlpaka<pixelTopology::HIonPhase1>;
0185 
0186 #include "FWCore/Framework/interface/MakerMacros.h"
0187 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpakaPhase1);
0188 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpakaPhase2);
0189 DEFINE_FWK_MODULE(SiPixelRecHitFromSoAAlpakaHIonPhase1);