Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:21:08

0001 // C++ includes
0002 #include <memory>
0003 #include <stdexcept>
0004 #include <string>
0005 #include <utility>
0006 #include <vector>
0007 
0008 #include <alpaka/alpaka.hpp>
0009 
0010 #include "DataFormats/Common/interface/DetSetVector.h"
0011 #include "DataFormats/SiPixelClusterSoA/interface/alpaka/SiPixelClustersSoACollection.h"
0012 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0013 #include "DataFormats/SiPixelDigiSoA/interface/alpaka/SiPixelDigiErrorsSoACollection.h"
0014 #include "DataFormats/SiPixelDigiSoA/interface/alpaka/SiPixelDigisSoACollection.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0018 #include "FWCore/Utilities/interface/ESGetToken.h"
0019 #include "FWCore/Utilities/interface/InputTag.h"
0020 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0021 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0022 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/EDPutToken.h"
0023 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/Event.h"
0024 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/stream/SynchronizingEDProducer.h"
0025 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0026 #include "RecoLocalTracker/SiPixelClusterizer/interface/SiPixelClusterThresholds.h"
0027 
0028 #include "SiPixelRawToClusterKernel.h"
0029 
0030 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0031 
0032   class SiPixelPhase2DigiToCluster : public stream::SynchronizingEDProducer<> {
0033   public:
0034     explicit SiPixelPhase2DigiToCluster(const edm::ParameterSet& iConfig);
0035     ~SiPixelPhase2DigiToCluster() override = default;
0036 
0037     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0038     using Algo = pixelDetails::SiPixelRawToClusterKernel<pixelTopology::Phase2>;
0039 
0040   private:
0041     void acquire(device::Event const& iEvent, device::EventSetup const& iSetup) override;
0042     void produce(device::Event& iEvent, device::EventSetup const& iSetup) override;
0043 
0044     const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0045     const edm::EDGetTokenT<edm::DetSetVector<PixelDigi>> pixelDigiToken_;
0046 
0047     device::EDPutToken<SiPixelDigisSoACollection> digiPutToken_;
0048     device::EDPutToken<SiPixelClustersSoACollection> clusterPutToken_;
0049 
0050     Algo Algo_;
0051 
0052     const SiPixelClusterThresholds clusterThresholds_;
0053     uint32_t nDigis_ = 0;
0054 
0055     SiPixelDigisSoACollection digis_d;
0056   };
0057 
0058   SiPixelPhase2DigiToCluster::SiPixelPhase2DigiToCluster(const edm::ParameterSet& iConfig)
0059       : geomToken_(esConsumes()),
0060         pixelDigiToken_(consumes<edm::DetSetVector<PixelDigi>>(iConfig.getParameter<edm::InputTag>("InputDigis"))),
0061         digiPutToken_(produces()),
0062         clusterPutToken_(produces()),
0063         clusterThresholds_{iConfig.getParameter<int32_t>("clusterThreshold_layer1"),
0064                            iConfig.getParameter<int32_t>("clusterThreshold_otherLayers"),
0065                            static_cast<float>(iConfig.getParameter<double>("ElectronPerADCGain")),
0066                            static_cast<int8_t>(iConfig.getParameter<int>("Phase2ReadoutMode")),
0067                            static_cast<uint16_t>(iConfig.getParameter<uint32_t>("Phase2DigiBaseline")),
0068                            static_cast<uint8_t>(iConfig.getParameter<uint32_t>("Phase2KinkADC"))} {}
0069 
0070   void SiPixelPhase2DigiToCluster::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0071     edm::ParameterSetDescription desc;
0072 
0073     desc.add<bool>("IncludeErrors", true);
0074     desc.add<int32_t>("clusterThreshold_layer1",
0075                       pixelClustering::clusterThresholdPhase2LayerOne);  //FIXME (fix the CUDA)
0076     desc.add<int32_t>("clusterThreshold_otherLayers", pixelClustering::clusterThresholdPhase2OtherLayers);
0077     desc.add<double>("ElectronPerADCGain", 1500.);
0078     desc.add<int32_t>("Phase2ReadoutMode", 3);
0079     desc.add<uint32_t>("Phase2DigiBaseline", 1000);
0080     desc.add<uint32_t>("Phase2KinkADC", 8);
0081     desc.add<edm::InputTag>("InputDigis", edm::InputTag("simSiPixelDigis:Pixel"));
0082     descriptions.addWithDefaultLabel(desc);
0083   }
0084 
0085   void SiPixelPhase2DigiToCluster::acquire(device::Event const& iEvent, device::EventSetup const& iSetup) {
0086     auto const& input = iEvent.get(pixelDigiToken_);
0087 
0088     const TrackerGeometry* geom_ = &iSetup.getData(geomToken_);
0089 
0090     uint32_t nDigis = 0;
0091 
0092     for (const auto& det : input) {
0093       nDigis += det.size();
0094     }
0095 
0096     if (nDigis == 0)
0097       return;
0098 
0099     nDigis_ = nDigis;
0100     SiPixelDigisHost digis_h(nDigis_, iEvent.queue());
0101 
0102     nDigis = 0;
0103     for (const auto& det : input) {
0104       unsigned int detid = det.detId();
0105       DetId detIdObject(detid);
0106       const GeomDetUnit* genericDet = geom_->idToDetUnit(detIdObject);
0107       auto const gind = genericDet->index();
0108       for (auto const& px : det) {
0109         digis_h.view()[nDigis].moduleId() = uint16_t(gind);
0110 
0111         digis_h.view()[nDigis].xx() = uint16_t(px.row());
0112         digis_h.view()[nDigis].yy() = uint16_t(px.column());
0113         digis_h.view()[nDigis].adc() = uint16_t(px.adc());
0114 
0115         digis_h.view()[nDigis].clus() = 0;
0116 
0117         digis_h.view()[nDigis].pdigi() = uint32_t(px.packedData());
0118 
0119         digis_h.view()[nDigis].rawIdArr() = uint32_t(detid);
0120 
0121         nDigis++;
0122       }
0123     }
0124 
0125     digis_d = SiPixelDigisSoACollection(nDigis, iEvent.queue());
0126     alpaka::memcpy(iEvent.queue(), digis_d.buffer(), digis_h.buffer());
0127 
0128     Algo_.makePhase2ClustersAsync(iEvent.queue(), clusterThresholds_, digis_d.view(), nDigis_);
0129   }
0130 
0131   void SiPixelPhase2DigiToCluster::produce(device::Event& iEvent, device::EventSetup const& iSetup) {
0132     if (nDigis_ == 0) {
0133       SiPixelClustersSoACollection clusters_d{pixelTopology::Phase2::numberOfModules, iEvent.queue()};
0134       SiPixelDigisSoACollection digis_d_zero{nDigis_, iEvent.queue()};
0135       iEvent.emplace(digiPutToken_, std::move(digis_d_zero));
0136       iEvent.emplace(clusterPutToken_, std::move(clusters_d));
0137       return;
0138     }
0139 
0140     digis_d.setNModules(Algo_.nModules());
0141     iEvent.emplace(digiPutToken_, std::move(digis_d));
0142     iEvent.emplace(clusterPutToken_, Algo_.getClusters());
0143   }
0144 
0145 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0146 
0147 // define as framework plugin
0148 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0149 DEFINE_FWK_ALPAKA_MODULE(SiPixelPhase2DigiToCluster);