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 
25 
27 
29 
31  public:
32  explicit SiPixelPhase2DigiToCluster(const edm::ParameterSet& iConfig);
33  ~SiPixelPhase2DigiToCluster() override = default;
34 
35  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
37 
38  private:
39  void acquire(device::Event const& iEvent, device::EventSetup const& iSetup) override;
40  void produce(device::Event& iEvent, device::EventSetup const& iSetup) override;
41 
44 
48 
50 
51  const bool includeErrors_;
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  includeErrors_(iConfig.getParameter<bool>("IncludeErrors")),
64  clusterThresholds_{iConfig.getParameter<int32_t>("clusterThreshold_layer1"),
65  iConfig.getParameter<int32_t>("clusterThreshold_otherLayers"),
66  static_cast<float>(iConfig.getParameter<double>("ElectronPerADCGain")),
67  static_cast<int8_t>(iConfig.getParameter<int>("Phase2ReadoutMode")),
68  static_cast<uint16_t>(iConfig.getParameter<uint32_t>("Phase2DigiBaseline")),
69  static_cast<uint8_t>(iConfig.getParameter<uint32_t>("Phase2KinkADC"))} {
70  if (includeErrors_) {
71  digiErrorPutToken_ = produces();
72  }
73  }
74 
77 
78  desc.add<bool>("IncludeErrors", true);
79  desc.add<int32_t>("clusterThreshold_layer1",
81  desc.add<int32_t>("clusterThreshold_otherLayers", pixelClustering::clusterThresholdPhase2OtherLayers);
82  desc.add<double>("ElectronPerADCGain", 1500.);
83  desc.add<int32_t>("Phase2ReadoutMode", 3);
84  desc.add<uint32_t>("Phase2DigiBaseline", 1000);
85  desc.add<uint32_t>("Phase2KinkADC", 8);
86  desc.add<edm::InputTag>("InputDigis", edm::InputTag("simSiPixelDigis:Pixel"));
87  descriptions.addWithDefaultLabel(desc);
88  }
89 
91  auto const& input = iEvent.get(pixelDigiToken_);
92 
93  const TrackerGeometry* geom_ = &iSetup.getData(geomToken_);
94 
95  uint32_t nDigis = 0;
96 
97  for (const auto& det : input) {
98  nDigis += det.size();
99  }
100 
101  if (nDigis_ == 0)
102  return;
103 
104  SiPixelDigisHost digis_h(nDigis, iEvent.queue());
105  nDigis_ = nDigis;
106 
107  nDigis = 0;
108  for (const auto& det : input) {
109  unsigned int detid = det.detId();
110  DetId detIdObject(detid);
111  const GeomDetUnit* genericDet = geom_->idToDetUnit(detIdObject);
112  auto const gind = genericDet->index();
113  for (auto const& px : det) {
114  digis_h.view()[nDigis].moduleId() = uint16_t(gind);
115 
116  digis_h.view()[nDigis].xx() = uint16_t(px.row());
117  digis_h.view()[nDigis].yy() = uint16_t(px.column());
118  digis_h.view()[nDigis].adc() = uint16_t(px.adc());
119 
120  digis_h.view()[nDigis].pdigi() = uint32_t(px.packedData());
121 
122  digis_h.view()[nDigis].rawIdArr() = uint32_t(detid);
123 
124  nDigis++;
125  }
126  }
127 
128  digis_d = SiPixelDigisSoACollection(nDigis, iEvent.queue());
129  alpaka::memcpy(iEvent.queue(), digis_d.buffer(), digis_h.buffer());
130 
132  }
133 
135  if (nDigis_ == 0) {
138  iEvent.emplace(clusterPutToken_, std::move(clusters_d));
139  if (includeErrors_) {
141  }
142  return;
143  }
144 
145  digis_d.setNModulesDigis(Algo_.nModules(), nDigis_);
146 
149  if (includeErrors_) {
151  }
152  }
153 
154 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
155 
156 // define as framework plugin
158 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.
device::EDPutToken< SiPixelDigiErrorsSoACollection > digiErrorPutToken_
int index() const
Definition: GeomDet.h:83
constexpr uint16_t clusterThresholdPhase2LayerOne
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
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, SiPixelDigiErrorsHost, SiPixelDigiErrorsDevice< Device > > SiPixelDigiErrorsSoACollection
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 constexpr uint16_t numberOfModules
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)