CMS 3D CMS Logo

SiPixelPhase2DigiToCluster.cc
Go to the documentation of this file.
1 // C++ includes
2 #include <memory>
3 #include <stdexcept>
4 #include <string>
5 #include <utility>
6 #include <vector>
7 
8 #include <alpaka/alpaka.hpp>
9 
27 
29 
31 
33  public:
34  explicit SiPixelPhase2DigiToCluster(const edm::ParameterSet& iConfig);
35  ~SiPixelPhase2DigiToCluster() override = default;
36 
37  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
39 
40  private:
41  void acquire(device::Event const& iEvent, device::EventSetup const& iSetup) override;
42  void produce(device::Event& iEvent, device::EventSetup const& iSetup) override;
43 
46 
49 
51 
53  uint32_t nDigis_ = 0;
54 
56  };
57 
59  : geomToken_(esConsumes()),
60  pixelDigiToken_(consumes<edm::DetSetVector<PixelDigi>>(iConfig.getParameter<edm::InputTag>("InputDigis"))),
61  digiPutToken_(produces()),
62  clusterPutToken_(produces()),
63  clusterThresholds_{iConfig.getParameter<int32_t>("clusterThreshold_layer1"),
64  iConfig.getParameter<int32_t>("clusterThreshold_otherLayers"),
65  static_cast<float>(iConfig.getParameter<double>("ElectronPerADCGain")),
66  static_cast<int8_t>(iConfig.getParameter<int>("Phase2ReadoutMode")),
67  static_cast<uint16_t>(iConfig.getParameter<uint32_t>("Phase2DigiBaseline")),
68  static_cast<uint8_t>(iConfig.getParameter<uint32_t>("Phase2KinkADC"))} {}
69 
72 
73  desc.add<bool>("IncludeErrors", true);
74  desc.add<int32_t>("clusterThreshold_layer1",
76  desc.add<int32_t>("clusterThreshold_otherLayers", pixelClustering::clusterThresholdPhase2OtherLayers);
77  desc.add<double>("ElectronPerADCGain", 1500.);
78  desc.add<int32_t>("Phase2ReadoutMode", 3);
79  desc.add<uint32_t>("Phase2DigiBaseline", 1000);
80  desc.add<uint32_t>("Phase2KinkADC", 8);
81  desc.add<edm::InputTag>("InputDigis", edm::InputTag("simSiPixelDigis:Pixel"));
82  descriptions.addWithDefaultLabel(desc);
83  }
84 
86  auto const& input = iEvent.get(pixelDigiToken_);
87 
88  const TrackerGeometry* geom_ = &iSetup.getData(geomToken_);
89 
90  uint32_t nDigis = 0;
91 
92  for (const auto& det : input) {
93  nDigis += det.size();
94  }
95 
96  if (nDigis == 0)
97  return;
98 
99  nDigis_ = nDigis;
100  SiPixelDigisHost digis_h(nDigis_, iEvent.queue());
101 
102  nDigis = 0;
103  for (const auto& det : input) {
104  unsigned int detid = det.detId();
105  DetId detIdObject(detid);
106  const GeomDetUnit* genericDet = geom_->idToDetUnit(detIdObject);
107  auto const gind = genericDet->index();
108  for (auto const& px : det) {
109  digis_h.view()[nDigis].moduleId() = uint16_t(gind);
110 
111  digis_h.view()[nDigis].xx() = uint16_t(px.row());
112  digis_h.view()[nDigis].yy() = uint16_t(px.column());
113  digis_h.view()[nDigis].adc() = uint16_t(px.adc());
114 
115  digis_h.view()[nDigis].clus() = 0;
116 
117  digis_h.view()[nDigis].pdigi() = uint32_t(px.packedData());
118 
119  digis_h.view()[nDigis].rawIdArr() = uint32_t(detid);
120 
121  nDigis++;
122  }
123  }
124 
125  digis_d = SiPixelDigisSoACollection(nDigis, iEvent.queue());
126  alpaka::memcpy(iEvent.queue(), digis_d.buffer(), digis_h.buffer());
127 
129  }
130 
132  if (nDigis_ == 0) {
134  SiPixelDigisSoACollection digis_d_zero{nDigis_, iEvent.queue()};
135  iEvent.emplace(digiPutToken_, std::move(digis_d_zero));
136  iEvent.emplace(clusterPutToken_, std::move(clusters_d));
137  return;
138  }
139 
140  digis_d.setNModules(Algo_.nModules());
143  }
144 
145 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
146 
147 // define as framework plugin
149 DEFINE_FWK_ALPAKA_MODULE(SiPixelPhase2DigiToCluster);
device::EDPutToken< SiPixelDigisSoACollection > digiPutToken_
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
int index() const
Definition: GeomDet.h:83
constexpr uint16_t clusterThresholdPhase2LayerOne
static constexpr uint16_t numberOfModules
void makePhase2ClustersAsync(Queue &queue, const SiPixelClusterThresholds clusterThresholds, SiPixelDigisSoAView &digis_view, const uint32_t numDigis)
static std::string const input
Definition: EdmProvDump.cc:50
int iEvent
Definition: GenABIO.cc:224
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, SiPixelDigisHost, SiPixelDigisDevice< Device > > SiPixelDigisSoACollection
void acquire(device::Event const &iEvent, device::EventSetup const &iSetup) override
device::EDPutToken< SiPixelClustersSoACollection > clusterPutToken_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
T const & getData(edm::ESGetToken< T, R > const &iToken) const
Definition: EventSetup.h:32
Definition: DetId.h:17
void produce(device::Event &iEvent, device::EventSetup const &iSetup) override
HLT enums.
constexpr uint16_t clusterThresholdPhase2OtherLayers
#define DEFINE_FWK_ALPAKA_MODULE(name)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > pixelDigiToken_
def move(src, dest)
Definition: eostools.py:511
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, SiPixelClustersHost, SiPixelClustersDevice< Device > > SiPixelClustersSoACollection
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)