CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelRecHitFromCUDA.cc
Go to the documentation of this file.
1 #include <cuda_runtime.h>
2 
3 #include <fmt/printf.h>
4 
26 
27 class SiPixelRecHitFromCUDA : public edm::stream::EDProducer<edm::ExternalWork> {
28 public:
30  ~SiPixelRecHitFromCUDA() override = default;
31 
32  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
33 
35 
36 private:
37  void acquire(edm::Event const& iEvent,
38  edm::EventSetup const& iSetup,
39  edm::WaitingTaskWithArenaHolder waitingTaskHolder) override;
40  void produce(edm::Event& iEvent, edm::EventSetup const& iSetup) override;
41 
47 
48  uint32_t nHits_;
51 };
52 
54  : geomToken_(esConsumes()),
55  hitsToken_(
56  consumes<cms::cuda::Product<TrackingRecHit2DGPU>>(iConfig.getParameter<edm::InputTag>("pixelRecHitSrc"))),
57  clusterToken_(consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("src"))),
58  rechitsPutToken_(produces<SiPixelRecHitCollection>()),
59  hostPutToken_(produces<HMSstorage>()) {}
60 
63  desc.add<edm::InputTag>("pixelRecHitSrc", edm::InputTag("siPixelRecHitsPreSplittingCUDA"));
64  desc.add<edm::InputTag>("src", edm::InputTag("siPixelClustersPreSplitting"));
65  descriptions.addWithDefaultLabel(desc);
66 }
67 
69  edm::EventSetup const& iSetup,
70  edm::WaitingTaskWithArenaHolder waitingTaskHolder) {
71  cms::cuda::Product<TrackingRecHit2DGPU> const& inputDataWrapped = iEvent.get(hitsToken_);
72  cms::cuda::ScopedContextAcquire ctx{inputDataWrapped, std::move(waitingTaskHolder)};
73  auto const& inputData = ctx.get(inputDataWrapped);
74 
75  nHits_ = inputData.nHits();
76 
77  LogDebug("SiPixelRecHitFromCUDA") << "converting " << nHits_ << " Hits";
78 
79  if (0 == nHits_)
80  return;
81  store32_ = inputData.localCoordToHostAsync(ctx.stream());
82  hitsModuleStart_ = inputData.hitsModuleStartToHostAsync(ctx.stream());
83 }
84 
86  // allocate a buffer for the indices of the clusters
87  auto hmsp = std::make_unique<uint32_t[]>(gpuClustering::maxNumModules + 1);
88 
90  output.reserve(gpuClustering::maxNumModules, nHits_);
91 
92  if (0 == nHits_) {
93  iEvent.emplace(rechitsPutToken_, std::move(output));
94  iEvent.emplace(hostPutToken_, std::move(hmsp));
95  return;
96  }
97  output.reserve(gpuClustering::maxNumModules, nHits_);
98 
100  // wrap the buffer in a HostProduct, and move it to the Event, without reallocating the buffer or affecting hitsModuleStart
101  iEvent.emplace(hostPutToken_, std::move(hmsp));
102 
103  auto xl = store32_.get();
104  auto yl = xl + nHits_;
105  auto xe = yl + nHits_;
106  auto ye = xe + nHits_;
107 
108  const TrackerGeometry* geom = &es.getData(geomToken_);
109 
111  auto const& input = *hclusters;
112 
113  constexpr uint32_t maxHitsInModule = gpuClustering::maxHitsInModule();
114 
115  int numberOfDetUnits = 0;
116  int numberOfClusters = 0;
117  for (auto const& dsv : input) {
118  numberOfDetUnits++;
119  unsigned int detid = dsv.detId();
120  DetId detIdObject(detid);
121  const GeomDetUnit* genericDet = geom->idToDetUnit(detIdObject);
122  auto gind = genericDet->index();
123  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
124  assert(pixDet);
125  SiPixelRecHitCollection::FastFiller recHitsOnDetUnit(output, detid);
126  auto fc = hitsModuleStart_[gind];
127  auto lc = hitsModuleStart_[gind + 1];
128  auto nhits = lc - fc;
129 
130  assert(lc > fc);
131  LogDebug("SiPixelRecHitFromCUDA") << "in det " << gind << ": conv " << nhits << " hits from " << dsv.size()
132  << " legacy clusters" << ' ' << fc << ',' << lc;
133  if (nhits > maxHitsInModule)
134  edm::LogWarning("SiPixelRecHitFromCUDA") << fmt::sprintf(
135  "Too many clusters %d in module %d. Only the first %d hits will be converted", nhits, gind, maxHitsInModule);
136  nhits = std::min(nhits, maxHitsInModule);
137 
138  LogDebug("SiPixelRecHitFromCUDA") << "in det " << gind << "conv " << nhits << " hits from " << dsv.size()
139  << " legacy clusters" << ' ' << lc << ',' << fc;
140 
141  if (0 == nhits)
142  continue;
143  auto jnd = [&](int k) { return fc + k; };
144  assert(nhits <= dsv.size());
145  if (nhits != dsv.size()) {
146  edm::LogWarning("GPUHits2CPU") << "nhits!= nclus " << nhits << ' ' << dsv.size();
147  }
148  for (auto const& clust : dsv) {
149  assert(clust.originalId() >= 0);
150  assert(clust.originalId() < dsv.size());
151  if (clust.originalId() >= nhits)
152  continue;
153  auto ij = jnd(clust.originalId());
154  LocalPoint lp(xl[ij], yl[ij]);
155  LocalError le(xe[ij], 0, ye[ij]);
157 
158  numberOfClusters++;
159 
160  /* cpu version.... (for reference)
161  std::tuple<LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType> tuple = cpe_->getParameters( clust, *genericDet );
162  LocalPoint lp( std::get<0>(tuple) );
163  LocalError le( std::get<1>(tuple) );
164  SiPixelRecHitQuality::QualWordType rqw( std::get<2>(tuple) );
165  */
166 
167  // Create a persistent edm::Ref to the cluster
169  // Make a RecHit and add it to the DetSet
170  recHitsOnDetUnit.emplace_back(lp, le, rqw, *genericDet, cluster);
171  // =============================
172 
173  LogDebug("SiPixelRecHitFromCUDA") << "cluster " << numberOfClusters << " at " << lp << ' ' << le;
174 
175  } // <-- End loop on Clusters
176 
177  // LogDebug("SiPixelRecHitGPU")
178  LogDebug("SiPixelRecHitFromCUDA") << "found " << recHitsOnDetUnit.size() << " RecHits on " << detid;
179 
180  } // <-- End loop on DetUnits
181 
182  LogDebug("SiPixelRecHitFromCUDA") << "found " << numberOfDetUnits << " dets, " << numberOfClusters << " clusters";
183 
184  iEvent.emplace(rechitsPutToken_, std::move(output));
185 }
186 
const edm::EDGetTokenT< SiPixelClusterCollectionNew > clusterToken_
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~SiPixelRecHitFromCUDA() override=default
constexpr uint32_t maxHitsInModule()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
assert(be >=bs)
static std::string const input
Definition: EdmProvDump.cc:47
bool getData(T &iHolder) const
Definition: EventSetup.h:128
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
int iEvent
Definition: GenABIO.cc:224
cms::cuda::host::unique_ptr< float[]> store32_
def move
Definition: eostools.py:511
SiPixelRecHitFromCUDA(const edm::ParameterSet &iConfig)
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
constexpr uint16_t maxNumModules
void acquire(edm::Event const &iEvent, edm::EventSetup const &iSetup, edm::WaitingTaskWithArenaHolder waitingTaskHolder) override
T min(T a, T b)
Definition: MathUtil.h:58
int index() const
Definition: GeomDet.h:83
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void produce(edm::Event &iEvent, edm::EventSetup const &iSetup) override
Definition: DetId.h:17
const edm::EDGetTokenT< cms::cuda::Product< TrackingRecHit2DGPU > > hitsToken_
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:433
cms::cuda::host::unique_ptr< uint32_t[]> hitsModuleStart_
Pixel cluster – collection of neighboring pixels above threshold.
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
void emplace_back(Args &&...args)
std::unique_ptr< T, impl::HostDeleter > unique_ptr
const edm::EDPutTokenT< HMSstorage > hostPutToken_
Log< level::Warning, false > LogWarning
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::EDPutTokenT< SiPixelRecHitCollection > rechitsPutToken_
#define LogDebug(id)