CMS 3D CMS Logo

SiPixelDigisClustersFromSoAAlpaka.cc
Go to the documentation of this file.
1 #include <memory>
2 
20 
21 // local include(s)
22 #include "PixelClusterizerBase.h"
23 
24 //#define EDM_ML_DEBUG
25 //#define GPU_DEBUG
26 
27 template <typename TrackerTraits>
29 public:
31  ~SiPixelDigisClustersFromSoAAlpaka() override = default;
32 
33  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
34 
35 private:
36  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
37 
40  const SiPixelClusterThresholds clusterThresholds_; // Cluster threshold in electrons
41  const bool produceDigis_;
42  const bool storeDigis_;
43 
46 };
47 
48 template <typename TrackerTraits>
50  : topoToken_(esConsumes()),
51  digisHostToken_(consumes(iConfig.getParameter<edm::InputTag>("src"))),
52  clusterThresholds_(iConfig.getParameter<int>("clusterThreshold_layer1"),
53  iConfig.getParameter<int>("clusterThreshold_otherLayers")),
54  produceDigis_(iConfig.getParameter<bool>("produceDigis")),
55  storeDigis_(produceDigis_ && iConfig.getParameter<bool>("storeDigis")),
56  clustersPutToken_(produces()) {
57  if (produceDigis_)
58  digisPutToken_ = produces<edm::DetSetVector<PixelDigi>>();
59 }
60 
61 template <typename TrackerTraits>
64  desc.add<edm::InputTag>("src", edm::InputTag("siPixelDigisSoA"));
65  desc.add<int>("clusterThreshold_layer1", pixelClustering::clusterThresholdLayerOne);
66  desc.add<int>("clusterThreshold_otherLayers", pixelClustering::clusterThresholdOtherLayers);
67  desc.add<bool>("produceDigis", true);
68  desc.add<bool>("storeDigis", true);
69 
70  descriptions.addWithDefaultLabel(desc);
71 }
72 
73 template <typename TrackerTraits>
76  const edm::EventSetup& iSetup) const {
77  const auto& digisHost = iEvent.get(digisHostToken_);
78  const auto& digisView = digisHost.const_view();
79  const uint32_t nDigis = digisHost.nDigis();
80 
81  const auto& ttopo = iSetup.getData(topoToken_);
82  constexpr auto maxModules = TrackerTraits::numberOfModules;
83 
84  std::unique_ptr<edm::DetSetVector<PixelDigi>> outputDigis;
85  if (produceDigis_)
86  outputDigis = std::make_unique<edm::DetSetVector<PixelDigi>>();
87  if (storeDigis_)
88  outputDigis->reserve(maxModules);
89  auto outputClusters = std::make_unique<SiPixelClusterCollectionNew>();
90  outputClusters->reserve(maxModules, nDigis / 2);
91 
92  edm::DetSet<PixelDigi>* detDigis = nullptr;
93  uint32_t detId = 0;
94 
95  for (uint32_t i = 0; i < nDigis; i++) {
96  // check for uninitialized digis
97  // this is set in RawToDigi_kernel in SiPixelRawToClusterGPUKernel.cu
98  if (digisView[i].rawIdArr() == 0)
99  continue;
100 
101  // check for noisy/dead pixels (electrons set to 0)
102  if (digisView[i].adc() == 0)
103  continue;
104 
105  detId = digisView[i].rawIdArr();
106  if (storeDigis_) {
107  detDigis = &outputDigis->find_or_insert(detId);
108 
109  if ((*detDigis).empty())
110  (*detDigis).data.reserve(64); // avoid the first relocations
111  }
112 
113  break;
114  }
115 
116  int32_t nclus = -1;
118 #ifdef EDM_ML_DEBUG
119  auto totClustersFilled = 0;
120 #endif
121 
122  auto fillClusters = [&](uint32_t detId) {
123  if (nclus < 0)
124  return; // this in reality should never happen
126  auto layer = (DetId(detId).subdetId() == 1) ? ttopo.pxbLayer(detId) : 0;
127  auto clusterThreshold = clusterThresholds_.getThresholdForLayerOnCondition(layer == 1);
128  for (int32_t ic = 0; ic < nclus + 1; ++ic) {
129  auto const& acluster = aclusters[ic];
130  // in any case we cannot go out of sync with gpu...
131  if (acluster.charge < clusterThreshold)
132  edm::LogWarning("SiPixelDigisClustersFromSoAAlpaka")
133  << "cluster below charge Threshold "
134  << "Layer/DetId/clusId " << layer << '/' << detId << '/' << ic << " size/charge " << acluster.isize << '/'
135  << acluster.charge << "\n";
136  // sort by row (x)
137  spc.emplace_back(acluster.isize, acluster.adc, acluster.x, acluster.y, acluster.xmin, acluster.ymin, ic);
138  aclusters[ic].clear();
139 #ifdef EDM_ML_DEBUG
140  ++totClustersFilled;
141  const auto& cluster{spc.back()};
142  // LogDebug("SiPixelDigisClustersFromSoAAlpaka")
143  std::cout << "putting in this cluster " << ic << " " << cluster.charge() << " " << cluster.pixelADC().size()
144  << "\n";
145 #endif
146  std::push_heap(spc.begin(), spc.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
147  return cl1.minPixelRow() < cl2.minPixelRow();
148  });
149  }
150  nclus = -1;
151  // sort by row (x)
152  std::sort_heap(spc.begin(), spc.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) {
153  return cl1.minPixelRow() < cl2.minPixelRow();
154  });
155  if (spc.empty())
156  spc.abort();
157  };
158 
159 #ifdef GPU_DEBUG
160  std::cout << "Dumping all digis. nDigis = " << nDigis << std::endl;
161 #endif
162  for (uint32_t i = 0; i < nDigis; i++) {
163 #ifdef GPU_DEBUG
164  PixelDigi dig2{digisView[i].pdigi()};
165  std::cout << i << ";" << digisView[i].rawIdArr() << ";" << digisView[i].clus() << ";" << digisView[i].pdigi() << ";"
166  << digisView[i].adc() << ";" << dig2.row() << ";" << dig2.column() << std::endl;
167 #endif
168 
169  // check for uninitialized digis
170  if (digisView[i].rawIdArr() == 0)
171  continue;
172  // check for noisy/dead pixels (electrons set to 0)
173  if (digisView[i].adc() == 0)
174  continue;
175  if (digisView[i].clus() >= -pixelClustering::invalidClusterId)
176  continue; // not in cluster; TODO add an assert for the size
177  if (digisView[i].clus() == pixelClustering::invalidModuleId)
178  continue; // from clusters killed by charge cut
179 #ifdef EDM_ML_DEBUG
180  assert(digisView[i].rawIdArr() > 109999);
181 #endif
182  if (detId != digisView[i].rawIdArr()) {
183 #ifdef GPU_DEBUG
184  std::cout << ">> Closed module --" << detId << "; nclus = " << nclus << std::endl;
185 #endif
186  // new module
187  fillClusters(detId);
188 #ifdef EDM_ML_DEBUG
189  assert(nclus == -1);
190 #endif
191  detId = digisView[i].rawIdArr();
192  if (storeDigis_) {
193  detDigis = &outputDigis->find_or_insert(detId);
194  if ((*detDigis).empty())
195  (*detDigis).data.reserve(64); // avoid the first relocations
196  else {
197  edm::LogWarning("SiPixelDigisClustersFromSoAAlpaka")
198  << "Problem det present twice in input! " << (*detDigis).detId();
199  }
200  }
201  }
202  PixelDigi dig{digisView[i].pdigi()};
203 
204  if (storeDigis_)
205  (*detDigis).data.emplace_back(dig);
206  // fill clusters
207 #ifdef EDM_ML_DEBUG
208  assert(digisView[i].clus() >= 0);
209  assert(digisView[i].clus() < static_cast<int>(TrackerTraits::maxNumClustersPerModules));
210 #endif
211  nclus = std::max(digisView[i].clus(), nclus);
212  auto row = dig.row();
213  auto col = dig.column();
214  SiPixelCluster::PixelPos pix(row, col);
215  aclusters[digisView[i].clus()].add(pix, digisView[i].adc());
216  }
217 
218  // fill final clusters
219  if (detId > 0)
220  fillClusters(detId);
221 
222 #ifdef EDM_ML_DEBUG
223  LogDebug("SiPixelDigisClustersFromSoAAlpaka") << "filled " << totClustersFilled << " clusters";
224 #endif
225 
226  if (produceDigis_)
227  iEvent.put(digisPutToken_, std::move(outputDigis));
228 
229  iEvent.put(clustersPutToken_, std::move(outputClusters));
230 }
231 
233 
236 
239 
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr uint16_t clusterThresholdOtherLayers
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
assert(be >=bs)
constexpr uint16_t numberOfModules
int minPixelRow() const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > const topoToken_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< SiPixelDigisHost > const digisHostToken_
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
constexpr uint16_t clusterThresholdLayerOne
edm::EDPutTokenT< edm::DetSetVector< PixelDigi > > digisPutToken_
Definition: DetId.h:17
constexpr int invalidClusterId
SiPixelDigisClustersFromSoAAlpaka(const edm::ParameterSet &iConfig)
collection_type data
Definition: DetSet.h:80
Pixel cluster – collection of neighboring pixels above threshold.
constexpr int32_t maxNumClustersPerModules
HLT enums.
void emplace_back(Args &&... args)
col
Definition: cuy.py:1009
bool add(SiPixelCluster::PixelPos const &p, uint16_t const iadc)
Log< level::Warning, false > LogWarning
~SiPixelDigisClustersFromSoAAlpaka() override=default
def move(src, dest)
Definition: eostools.py:511
edm::EDPutTokenT< SiPixelClusterCollectionNew > clustersPutToken_
constexpr uint16_t invalidModuleId
uint16_t *__restrict__ uint16_t const *__restrict__ adc
#define LogDebug(id)