CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelRecHitSoAFromLegacy.cc
Go to the documentation of this file.
1 #include <cuda_runtime.h>
2 
26 
27 #include "gpuPixelRecHits.h"
28 
30 public:
32  ~SiPixelRecHitSoAFromLegacy() override = default;
33 
34  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
35 
36  using HitModuleStart = std::array<uint32_t, gpuClustering::maxNumModules + 1>;
38 
39 private:
40  void produce(edm::StreamID streamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
41 
48  const bool convert2Legacy_;
49 };
50 
52  : geomToken_(esConsumes()),
53  cpeToken_(esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("CPE")))),
54  bsGetToken_{consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))},
55  clusterToken_{consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("src"))},
56  tokenHit_{produces<TrackingRecHit2DCPU>()},
57  tokenModuleStart_{produces<HMSstorage>()},
58  convert2Legacy_(iConfig.getParameter<bool>("convertToLegacy")) {
59  if (convert2Legacy_)
60  produces<SiPixelRecHitCollectionNew>();
61 }
62 
65 
66  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
67  desc.add<edm::InputTag>("src", edm::InputTag("siPixelClustersPreSplitting"));
68  desc.add<std::string>("CPE", "PixelCPEFast");
69  desc.add<bool>("convertToLegacy", false);
70  descriptions.addWithDefaultLabel(desc);
71 }
72 
74  const TrackerGeometry* geom_ = &es.getData(geomToken_);
75  PixelCPEFast const* fcpe = dynamic_cast<const PixelCPEFast*>(&es.getData(cpeToken_));
76  if (not fcpe) {
77  throw cms::Exception("Configuration") << "SiPixelRecHitSoAFromLegacy can only use a CPE of type PixelCPEFast";
78  }
79  auto const& cpeView = fcpe->getCPUProduct();
80 
81  const reco::BeamSpot& bs = iEvent.get(bsGetToken_);
82 
83  BeamSpotPOD bsHost;
84  bsHost.x = bs.x0();
85  bsHost.y = bs.y0();
86  bsHost.z = bs.z0();
87 
89  iEvent.getByToken(clusterToken_, hclusters);
90  auto const& input = *hclusters;
91 
92  // allocate a buffer for the indices of the clusters
93  auto hmsp = std::make_unique<uint32_t[]>(gpuClustering::maxNumModules + 1);
94  // hitsModuleStart is a non-owning pointer to the buffer
95  auto hitsModuleStart = hmsp.get();
96  // wrap the buffer in a HostProduct
97  auto hms = std::make_unique<HMSstorage>(std::move(hmsp));
98  // move the HostProduct to the Event, without reallocating the buffer or affecting hitsModuleStart
99  iEvent.put(tokenModuleStart_, std::move(hms));
100 
101  // legacy output
102  auto legacyOutput = std::make_unique<SiPixelRecHitCollectionNew>();
103 
104  // storage
105  std::vector<uint16_t> xx;
106  std::vector<uint16_t> yy;
107  std::vector<uint16_t> adc;
108  std::vector<uint16_t> moduleInd;
109  std::vector<int32_t> clus;
110 
111  std::vector<edm::Ref<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster>> clusterRef;
112 
113  constexpr uint32_t maxHitsInModule = gpuClustering::maxHitsInModule();
114 
115  HitModuleStart moduleStart_; // index of the first pixel of each module
116  HitModuleStart clusInModule_;
117  memset(&clusInModule_, 0, sizeof(HitModuleStart)); // needed??
118  assert(gpuClustering::maxNumModules + 1 == clusInModule_.size());
119  assert(0 == clusInModule_[gpuClustering::maxNumModules]);
120  uint32_t moduleId_;
121  moduleStart_[1] = 0; // we run sequentially....
122 
124  moduleStart_.data(), clusInModule_.data(), &moduleId_, hitsModuleStart};
125 
126  // fill cluster arrays
127  int numberOfClusters = 0;
128  for (auto const& dsv : input) {
129  unsigned int detid = dsv.detId();
130  DetId detIdObject(detid);
131  const GeomDetUnit* genericDet = geom_->idToDetUnit(detIdObject);
132  auto gind = genericDet->index();
133  assert(gind < gpuClustering::maxNumModules);
134  auto const nclus = dsv.size();
135  clusInModule_[gind] = nclus;
136  numberOfClusters += nclus;
137  }
138  hitsModuleStart[0] = 0;
139  for (int i = 1, n = clusInModule_.size(); i < n; ++i)
140  hitsModuleStart[i] = hitsModuleStart[i - 1] + clusInModule_[i - 1];
141  assert(numberOfClusters == int(hitsModuleStart[gpuClustering::maxNumModules]));
142 
143  // output SoA
144  // element 96 is the start of BPIX2 (i.e. the number of clusters in BPIX1)
145  auto output =
146  std::make_unique<TrackingRecHit2DCPU>(numberOfClusters, hitsModuleStart[96], &cpeView, hitsModuleStart, nullptr);
147 
148  if (0 == numberOfClusters) {
149  iEvent.put(std::move(output));
150  if (convert2Legacy_)
151  iEvent.put(std::move(legacyOutput));
152  return;
153  }
154 
155  if (convert2Legacy_)
156  legacyOutput->reserve(gpuClustering::maxNumModules, numberOfClusters);
157 
158  int numberOfDetUnits = 0;
159  int numberOfHits = 0;
160  for (auto const& dsv : input) {
161  numberOfDetUnits++;
162  unsigned int detid = dsv.detId();
163  DetId detIdObject(detid);
164  const GeomDetUnit* genericDet = geom_->idToDetUnit(detIdObject);
165  auto const gind = genericDet->index();
166  assert(gind < gpuClustering::maxNumModules);
167  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
168  assert(pixDet);
169  auto const nclus = dsv.size();
170  assert(clusInModule_[gind] == nclus);
171  if (0 == nclus)
172  continue; // is this really possible?
173 
174  auto const fc = hitsModuleStart[gind];
175  auto const lc = hitsModuleStart[gind + 1];
176  assert(lc > fc);
177  LogDebug("SiPixelRecHitSoAFromLegacy") << "in det " << gind << ": conv " << nclus << " hits from " << dsv.size()
178  << " legacy clusters" << ' ' << fc << ',' << lc;
179  assert((lc - fc) == nclus);
180  if (nclus > maxHitsInModule)
181  printf(
182  "WARNING: too many clusters %d in Module %d. Only first %d Hits converted\n", nclus, gind, maxHitsInModule);
183 
184  // fill digis
185  xx.clear();
186  yy.clear();
187  adc.clear();
188  moduleInd.clear();
189  clus.clear();
190  clusterRef.clear();
191  moduleId_ = gind;
192  uint32_t ic = 0;
193  uint32_t ndigi = 0;
194  for (auto const& clust : dsv) {
195  assert(clust.size() > 0);
196  for (int i = 0, nd = clust.size(); i < nd; ++i) {
197  auto px = clust.pixel(i);
198  xx.push_back(px.x);
199  yy.push_back(px.y);
200  adc.push_back(px.adc);
201  moduleInd.push_back(gind);
202  clus.push_back(ic);
203  ++ndigi;
204  }
205  assert(clust.originalId() == ic); // make sure hits and clus are in sync
206  if (convert2Legacy_)
207  clusterRef.emplace_back(edmNew::makeRefTo(hclusters, &clust));
208  ic++;
209  }
210  assert(nclus == ic);
211  assert(clus.size() == ndigi);
212  numberOfHits += nclus;
213  // filled creates view
214  SiPixelDigisCUDA::DeviceConstView digiView{xx.data(), yy.data(), adc.data(), moduleInd.data(), clus.data()};
215  assert(digiView.adc(0) != 0);
216  // we run on blockId.x==0
217  gpuPixelRecHits::getHits(&cpeView, &bsHost, &digiView, ndigi, &clusterView, output->view());
218  for (auto h = fc; h < lc; ++h)
219  if (h - fc < maxHitsInModule)
220  assert(gind == output->view()->detectorIndex(h));
221  else
222  assert(gpuClustering::invalidModuleId == output->view()->detectorIndex(h));
223  if (convert2Legacy_) {
224  SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(*legacyOutput, detid);
225  for (auto h = fc; h < lc; ++h) {
226  auto ih = h - fc;
227  if (ih >= maxHitsInModule)
228  break;
229  assert(ih < clusterRef.size());
230  LocalPoint lp(output->view()->xLocal(h), output->view()->yLocal(h));
231  LocalError le(output->view()->xerrLocal(h), 0, output->view()->yerrLocal(h));
233  SiPixelRecHit hit(lp, le, rqw, *genericDet, clusterRef[ih]);
234  recHitsOnDetUnit.push_back(hit);
235  }
236  }
237  }
238  assert(numberOfHits == numberOfClusters);
239 
240  // fill data structure to support CA
241  for (auto i = 0U; i < phase1PixelTopology::numberOfLayers + 1; ++i) {
242  output->hitsLayerStart()[i] = hitsModuleStart[cpeView.layerGeometry().layerStart[i]];
243  }
244  cms::cuda::fillManyFromVector(output->phiBinner(),
246  output->iphi(),
247  output->hitsLayerStart(),
248  numberOfHits,
249  256,
250  output->phiBinnerStorage());
251 
252  LogDebug("SiPixelRecHitSoAFromLegacy") << "created HitSoa for " << numberOfClusters << " clusters in "
253  << numberOfDetUnits << " Dets";
254  iEvent.put(std::move(output));
255  if (convert2Legacy_)
256  iEvent.put(std::move(legacyOutput));
257 }
258 
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)
double z0() const
z coordinate
Definition: BeamSpot.h:65
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void push_back(data_type const &d)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
constexpr uint32_t numberOfLayers
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr uint32_t maxHitsInModule()
const edm::EDGetTokenT< reco::BeamSpot > bsGetToken_
assert(be >=bs)
~SiPixelRecHitSoAFromLegacy() override=default
static std::string const input
Definition: EdmProvDump.cc:47
bool getData(T &iHolder) const
Definition: EventSetup.h:128
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
printf("params %d %f %f %f\n", minT, eps, errmax, chi2max)
def move
Definition: eostools.py:511
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
constexpr uint16_t maxNumModules
int index() const
Definition: GeomDet.h:83
ParameterDescriptionBase * add(U const &iLabel, T const &value)
constexpr uint16_t invalidModuleId
Definition: DetId.h:17
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::array< uint32_t, gpuClustering::maxNumModules+1 > HitModuleStart
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > cpeToken_
Pixel cluster – collection of neighboring pixels above threshold.
const edm::EDPutTokenT< HMSstorage > tokenModuleStart_
double y0() const
y coordinate
Definition: BeamSpot.h:63
const edm::EDPutTokenT< TrackingRecHit2DCPU > tokenHit_
const edm::EDGetTokenT< SiPixelClusterCollectionNew > clusterToken_
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
SiPixelRecHitSoAFromLegacy(const edm::ParameterSet &iConfig)
uint16_t *__restrict__ uint16_t const *__restrict__ adc
Our base class.
Definition: SiPixelRecHit.h:23
#define LogDebug(id)
double x0() const
x coordinate
Definition: BeamSpot.h:61