Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-15 02:28:15

0001 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0002 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0003 #include "HeterogeneousCore/AlpakaInterface/interface/CopyToDevice.h"
0004 
0005 #include "LSTEvent.h"
0006 
0007 #include "Hit.h"
0008 #include "Kernels.h"
0009 #include "MiniDoublet.h"
0010 #include "PixelQuintuplet.h"
0011 #include "PixelTriplet.h"
0012 #include "Quintuplet.h"
0013 #include "Segment.h"
0014 #include "TrackCandidate.h"
0015 #include "Triplet.h"
0016 
0017 using Device = ALPAKA_ACCELERATOR_NAMESPACE::Device;
0018 using Queue = ALPAKA_ACCELERATOR_NAMESPACE::Queue;
0019 using Acc1D = ALPAKA_ACCELERATOR_NAMESPACE::Acc1D;
0020 using Acc3D = ALPAKA_ACCELERATOR_NAMESPACE::Acc3D;
0021 
0022 using namespace ALPAKA_ACCELERATOR_NAMESPACE::lst;
0023 
0024 void LSTEvent::initSync() {
0025   alpaka::wait(queue_);  // other calls can be asynchronous
0026 
0027   //reset the arrays
0028   for (int i = 0; i < 6; i++) {
0029     n_minidoublets_by_layer_barrel_[i] = 0;
0030     n_segments_by_layer_barrel_[i] = 0;
0031     n_triplets_by_layer_barrel_[i] = 0;
0032     n_quintuplets_by_layer_barrel_[i] = 0;
0033     if (i < 5) {
0034       n_minidoublets_by_layer_endcap_[i] = 0;
0035       n_segments_by_layer_endcap_[i] = 0;
0036       n_triplets_by_layer_endcap_[i] = 0;
0037       n_quintuplets_by_layer_endcap_[i] = 0;
0038     }
0039   }
0040 }
0041 
0042 void LSTEvent::resetEventSync() {
0043   alpaka::wait(queue_);  // synchronize to reset consistently
0044   //reset the arrays
0045   for (int i = 0; i < 6; i++) {
0046     n_minidoublets_by_layer_barrel_[i] = 0;
0047     n_segments_by_layer_barrel_[i] = 0;
0048     n_triplets_by_layer_barrel_[i] = 0;
0049     n_quintuplets_by_layer_barrel_[i] = 0;
0050     if (i < 5) {
0051       n_minidoublets_by_layer_endcap_[i] = 0;
0052       n_segments_by_layer_endcap_[i] = 0;
0053       n_triplets_by_layer_endcap_[i] = 0;
0054       n_quintuplets_by_layer_endcap_[i] = 0;
0055     }
0056   }
0057   lstInputDC_ = nullptr;
0058   hitsDC_.reset();
0059   rangesDC_.reset();
0060   miniDoubletsDC_.reset();
0061   segmentsDC_.reset();
0062   pixelSegmentsDC_.reset();
0063   tripletsDC_.reset();
0064   quintupletsDC_.reset();
0065   trackCandidatesDC_.reset();
0066   pixelTripletsDC_.reset();
0067   pixelQuintupletsDC_.reset();
0068 
0069   lstInputHC_.reset();
0070   hitsHC_.reset();
0071   rangesHC_.reset();
0072   miniDoubletsHC_.reset();
0073   segmentsHC_.reset();
0074   pixelSegmentsHC_.reset();
0075   tripletsHC_.reset();
0076   quintupletsHC_.reset();
0077   pixelTripletsHC_.reset();
0078   pixelQuintupletsHC_.reset();
0079   trackCandidatesHC_.reset();
0080   modulesHC_.reset();
0081 }
0082 
0083 void LSTEvent::addInputToEvent(LSTInputDeviceCollection const* lstInputDC) {
0084   lstInputDC_ = lstInputDC;
0085 
0086   pixelSize_ = lstInputDC_->sizes()[1];
0087   pixelModuleIndex_ = pixelMapping_.pixelModuleIndex;
0088 }
0089 
0090 void LSTEvent::addHitToEvent() {
0091   if (!hitsDC_) {
0092     int nHits = lstInputDC_->sizes()[0];
0093     std::array<int, 2> const hits_sizes{{nHits, static_cast<int>(nModules_)}};
0094     hitsDC_.emplace(hits_sizes, queue_);
0095     auto buf = hitsDC_->buffer();
0096     alpaka::memset(queue_, buf, 0xff);
0097   }
0098 
0099   if (!rangesDC_) {
0100     rangesDC_.emplace(nLowerModules_ + 1, queue_);
0101     auto buf = rangesDC_->buffer();
0102     alpaka::memset(queue_, buf, 0xff);
0103   }
0104 
0105   auto const hit_loop_workdiv = cms::alpakatools::make_workdiv<Acc1D>(max_blocks, 256);
0106 
0107   alpaka::exec<Acc1D>(queue_,
0108                       hit_loop_workdiv,
0109                       HitLoopKernel{},
0110                       Endcap,
0111                       TwoS,
0112                       nModules_,
0113                       nEndCapMap_,
0114                       endcapGeometry_.const_view(),
0115                       modules_.const_view<ModulesSoA>(),
0116                       lstInputDC_->const_view<HitsBaseSoA>(),
0117                       hitsDC_->view<HitsExtendedSoA>(),
0118                       hitsDC_->view<HitsRangesSoA>());
0119 
0120   auto const module_ranges_workdiv = cms::alpakatools::make_workdiv<Acc1D>(max_blocks, 256);
0121 
0122   alpaka::exec<Acc1D>(queue_,
0123                       module_ranges_workdiv,
0124                       ModuleRangesKernel{},
0125                       modules_.const_view<ModulesSoA>(),
0126                       hitsDC_->view<HitsRangesSoA>(),
0127                       nLowerModules_);
0128 }
0129 
0130 void LSTEvent::addPixelSegmentToEventStart() {
0131   if (pixelSize_ == n_max_pixel_segments_per_module) {
0132     lstWarning(
0133         "\
0134           *********************************************************\n\
0135           * Warning: Pixel line segments may be truncated.        *\n\
0136           * You need to increase n_max_pixel_segments_per_module. *\n\
0137           *********************************************************");
0138   }
0139 
0140   if (!pixelSegmentsDC_) {
0141     pixelSegmentsDC_.emplace(pixelSize_, queue_);
0142   }
0143 }
0144 
0145 void LSTEvent::addPixelSegmentToEventFinalize() {
0146   auto const addPixelSegmentToEvent_workdiv = cms::alpakatools::make_workdiv<Acc1D>(max_blocks, 256);
0147 
0148   alpaka::exec<Acc1D>(queue_,
0149                       addPixelSegmentToEvent_workdiv,
0150                       AddPixelSegmentToEventKernel{},
0151                       modules_.const_view<ModulesSoA>(),
0152                       rangesDC_->const_view(),
0153                       lstInputDC_->const_view<HitsBaseSoA>(),
0154                       hitsDC_->view<HitsExtendedSoA>(),
0155                       lstInputDC_->const_view<PixelSeedsSoA>(),
0156                       miniDoubletsDC_->view<MiniDoubletsSoA>(),
0157                       segmentsDC_->view<SegmentsSoA>(),
0158                       pixelSegmentsDC_->view(),
0159                       pixelModuleIndex_,
0160                       pixelSize_);
0161 }
0162 
0163 void LSTEvent::createMiniDoublets() {
0164   if (!miniDoubletsDC_) {
0165     // Create a view for the element nLowerModules_ inside rangesOccupancy->miniDoubletModuleOccupancy
0166     auto rangesOccupancy = rangesDC_->view();
0167     auto dst_view_miniDoubletModuleOccupancy =
0168         cms::alpakatools::make_device_view(queue_, rangesOccupancy.miniDoubletModuleOccupancy()[nLowerModules_]);
0169 
0170     // Create a host buffer for a value to be passed to the device
0171     auto pixelMaxMDs_buf_h = cms::alpakatools::make_host_buffer<int>(queue_);
0172     *pixelMaxMDs_buf_h.data() = n_max_pixel_md_per_modules;
0173 
0174     alpaka::memcpy(queue_, dst_view_miniDoubletModuleOccupancy, pixelMaxMDs_buf_h);
0175 
0176     auto dst_view_miniDoubletModuleOccupancyPix =
0177         cms::alpakatools::make_device_view(queue_, rangesOccupancy.miniDoubletModuleOccupancy()[pixelModuleIndex_]);
0178 
0179     alpaka::memcpy(queue_, dst_view_miniDoubletModuleOccupancyPix, pixelMaxMDs_buf_h);
0180 
0181     auto const createMDArrayRangesGPU_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0182 
0183     alpaka::exec<Acc1D>(queue_,
0184                         createMDArrayRangesGPU_workDiv,
0185                         CreateMDArrayRangesGPU{},
0186                         modules_.const_view<ModulesSoA>(),
0187                         hitsDC_->const_view<HitsRangesSoA>(),
0188                         rangesDC_->view(),
0189                         ptCut_);
0190 
0191     auto nTotalMDs_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0192     auto nTotalMDs_buf_d = cms::alpakatools::make_device_view(queue_, rangesOccupancy.nTotalMDs());
0193     alpaka::memcpy(queue_, nTotalMDs_buf_h, nTotalMDs_buf_d);
0194     alpaka::wait(queue_);  // wait to get the data before manipulation
0195 
0196     *nTotalMDs_buf_h.data() += n_max_pixel_md_per_modules;
0197     unsigned int nTotalMDs = *nTotalMDs_buf_h.data();
0198 
0199     std::array<int, 2> const mds_sizes{{static_cast<int>(nTotalMDs), static_cast<int>(nLowerModules_ + 1)}};
0200     miniDoubletsDC_.emplace(mds_sizes, queue_);
0201 
0202     auto mdsOccupancy = miniDoubletsDC_->view<MiniDoubletsOccupancySoA>();
0203     auto nMDs_view = cms::alpakatools::make_device_view(queue_, mdsOccupancy.nMDs(), mdsOccupancy.metadata().size());
0204     auto totOccupancyMDs_view =
0205         cms::alpakatools::make_device_view(queue_, mdsOccupancy.totOccupancyMDs(), mdsOccupancy.metadata().size());
0206     alpaka::memset(queue_, nMDs_view, 0u);
0207     alpaka::memset(queue_, totOccupancyMDs_view, 0u);
0208   }
0209 
0210   unsigned int mdSize = pixelSize_ * 2;
0211   auto src_view_mdSize = cms::alpakatools::make_host_view(mdSize);
0212 
0213   auto mdsOccupancy = miniDoubletsDC_->view<MiniDoubletsOccupancySoA>();
0214   auto dst_view_nMDs = cms::alpakatools::make_device_view(queue_, mdsOccupancy.nMDs()[pixelModuleIndex_]);
0215   alpaka::memcpy(queue_, dst_view_nMDs, src_view_mdSize);
0216 
0217   auto dst_view_totOccupancyMDs =
0218       cms::alpakatools::make_device_view(queue_, mdsOccupancy.totOccupancyMDs()[pixelModuleIndex_]);
0219   alpaka::memcpy(queue_, dst_view_totOccupancyMDs, src_view_mdSize);
0220 
0221   alpaka::wait(queue_);  // FIXME: remove synch after inputs refactored to be in pinned memory
0222 
0223   constexpr int threadsPerBlockY = 16;
0224   auto const createMiniDoublets_workDiv =
0225       cms::alpakatools::make_workdiv<Acc2D>({nLowerModules_ / threadsPerBlockY, 1}, {threadsPerBlockY, 32});
0226 
0227   alpaka::exec<Acc2D>(queue_,
0228                       createMiniDoublets_workDiv,
0229                       CreateMiniDoublets{},
0230                       modules_.const_view<ModulesSoA>(),
0231                       lstInputDC_->const_view<HitsBaseSoA>(),
0232                       hitsDC_->const_view<HitsExtendedSoA>(),
0233                       hitsDC_->const_view<HitsRangesSoA>(),
0234                       miniDoubletsDC_->view<MiniDoubletsSoA>(),
0235                       miniDoubletsDC_->view<MiniDoubletsOccupancySoA>(),
0236                       rangesDC_->const_view(),
0237                       ptCut_);
0238 
0239   auto const addMiniDoubletRangesToEventExplicit_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0240 
0241   alpaka::exec<Acc1D>(queue_,
0242                       addMiniDoubletRangesToEventExplicit_workDiv,
0243                       AddMiniDoubletRangesToEventExplicit{},
0244                       modules_.const_view<ModulesSoA>(),
0245                       miniDoubletsDC_->view<MiniDoubletsOccupancySoA>(),
0246                       rangesDC_->view(),
0247                       hitsDC_->const_view<HitsRangesSoA>());
0248 
0249   if (addObjects_) {
0250     addMiniDoubletsToEventExplicit();
0251   }
0252 }
0253 
0254 void LSTEvent::createSegmentsWithModuleMap() {
0255   if (!segmentsDC_) {
0256     auto const createSegmentArrayRanges_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0257 
0258     alpaka::exec<Acc1D>(queue_,
0259                         createSegmentArrayRanges_workDiv,
0260                         CreateSegmentArrayRanges{},
0261                         modules_.const_view<ModulesSoA>(),
0262                         rangesDC_->view(),
0263                         miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0264                         miniDoubletsDC_->const_view<MiniDoubletsOccupancySoA>(),
0265                         ptCut_);
0266 
0267     auto rangesOccupancy = rangesDC_->view();
0268     auto nTotalSegments_view_h = cms::alpakatools::make_host_view(nTotalSegments_);
0269     auto nTotalSegments_view_d = cms::alpakatools::make_device_view(queue_, rangesOccupancy.nTotalSegs());
0270     alpaka::memcpy(queue_, nTotalSegments_view_h, nTotalSegments_view_d);
0271     alpaka::wait(queue_);  // wait to get the value before manipulation
0272 
0273     nTotalSegments_ += n_max_pixel_segments_per_module;
0274 
0275     std::array<int, 2> const segments_sizes{{static_cast<int>(nTotalSegments_), static_cast<int>(nLowerModules_ + 1)}};
0276     segmentsDC_.emplace(segments_sizes, queue_);
0277 
0278     auto segmentsOccupancy = segmentsDC_->view<SegmentsOccupancySoA>();
0279     auto nSegments_view =
0280         cms::alpakatools::make_device_view(queue_, segmentsOccupancy.nSegments(), segmentsOccupancy.metadata().size());
0281     auto totOccupancySegments_view = cms::alpakatools::make_device_view(
0282         queue_, segmentsOccupancy.totOccupancySegments(), segmentsOccupancy.metadata().size());
0283     alpaka::memset(queue_, nSegments_view, 0u);
0284     alpaka::memset(queue_, totOccupancySegments_view, 0u);
0285 
0286     auto src_view_size = cms::alpakatools::make_host_view(pixelSize_);
0287 
0288     auto dst_view_segments =
0289         cms::alpakatools::make_device_view(queue_, segmentsOccupancy.nSegments()[pixelModuleIndex_]);
0290     alpaka::memcpy(queue_, dst_view_segments, src_view_size);
0291 
0292     auto dst_view_totOccupancySegments =
0293         cms::alpakatools::make_device_view(queue_, segmentsOccupancy.totOccupancySegments()[pixelModuleIndex_]);
0294     alpaka::memcpy(queue_, dst_view_totOccupancySegments, src_view_size);
0295     alpaka::wait(queue_);
0296   }
0297 
0298   auto const createSegments_workDiv = cms::alpakatools::make_workdiv<Acc3D>({nLowerModules_, 1, 1}, {1, 1, 64});
0299 
0300   alpaka::exec<Acc3D>(queue_,
0301                       createSegments_workDiv,
0302                       CreateSegments{},
0303                       modules_.const_view<ModulesSoA>(),
0304                       miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0305                       miniDoubletsDC_->const_view<MiniDoubletsOccupancySoA>(),
0306                       segmentsDC_->view<SegmentsSoA>(),
0307                       segmentsDC_->view<SegmentsOccupancySoA>(),
0308                       rangesDC_->const_view(),
0309                       ptCut_);
0310 
0311   auto const addSegmentRangesToEventExplicit_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0312 
0313   alpaka::exec<Acc1D>(queue_,
0314                       addSegmentRangesToEventExplicit_workDiv,
0315                       AddSegmentRangesToEventExplicit{},
0316                       modules_.const_view<ModulesSoA>(),
0317                       segmentsDC_->view<SegmentsOccupancySoA>(),
0318                       rangesDC_->view());
0319 
0320   if (addObjects_) {
0321     addSegmentsToEventExplicit();
0322   }
0323 }
0324 
0325 void LSTEvent::createTriplets() {
0326   if (!tripletsDC_) {
0327     auto const createTripletArrayRanges_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0328 
0329     alpaka::exec<Acc1D>(queue_,
0330                         createTripletArrayRanges_workDiv,
0331                         CreateTripletArrayRanges{},
0332                         modules_.const_view<ModulesSoA>(),
0333                         rangesDC_->view(),
0334                         segmentsDC_->const_view<SegmentsSoA>(),
0335                         segmentsDC_->const_view<SegmentsOccupancySoA>(),
0336                         ptCut_);
0337 
0338     // TODO: Why are we pulling this back down only to put it back on the device in a new struct?
0339     auto rangesOccupancy = rangesDC_->view();
0340     auto maxTriplets_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0341     auto maxTriplets_buf_d = cms::alpakatools::make_device_view(queue_, rangesOccupancy.nTotalTrips());
0342     alpaka::memcpy(queue_, maxTriplets_buf_h, maxTriplets_buf_d);
0343     alpaka::wait(queue_);  // wait to get the value before using it
0344 
0345     std::array<int, 2> const triplets_sizes{
0346         {static_cast<int>(*maxTriplets_buf_h.data()), static_cast<int>(nLowerModules_)}};
0347     tripletsDC_.emplace(triplets_sizes, queue_);
0348 
0349     auto tripletsOccupancy = tripletsDC_->view<TripletsOccupancySoA>();
0350     auto nTriplets_view =
0351         cms::alpakatools::make_device_view(queue_, tripletsOccupancy.nTriplets(), tripletsOccupancy.metadata().size());
0352     alpaka::memset(queue_, nTriplets_view, 0u);
0353     auto totOccupancyTriplets_view = cms::alpakatools::make_device_view(
0354         queue_, tripletsOccupancy.totOccupancyTriplets(), tripletsOccupancy.metadata().size());
0355     alpaka::memset(queue_, totOccupancyTriplets_view, 0u);
0356     auto triplets = tripletsDC_->view<TripletsSoA>();
0357     auto partOfPT5_view = cms::alpakatools::make_device_view(queue_, triplets.partOfPT5(), triplets.metadata().size());
0358     alpaka::memset(queue_, partOfPT5_view, 0u);
0359     auto partOfT5_view = cms::alpakatools::make_device_view(queue_, triplets.partOfT5(), triplets.metadata().size());
0360     alpaka::memset(queue_, partOfT5_view, 0u);
0361     auto partOfPT3_view = cms::alpakatools::make_device_view(queue_, triplets.partOfPT3(), triplets.metadata().size());
0362     alpaka::memset(queue_, partOfPT3_view, 0u);
0363   }
0364 
0365   uint16_t nonZeroModules = 0;
0366   unsigned int max_InnerSeg = 0;
0367 
0368   // Allocate and copy nSegments from device to host (only nLowerModules in OT, not the +1 with pLSs)
0369   auto nSegments_buf_h = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nLowerModules_);
0370   auto nSegments_buf_d = cms::alpakatools::make_device_view(
0371       queue_, segmentsDC_->const_view<SegmentsOccupancySoA>().nSegments(), nLowerModules_);
0372   alpaka::memcpy(queue_, nSegments_buf_h, nSegments_buf_d, nLowerModules_);
0373 
0374   // ... same for module_nConnectedModules
0375   // FIXME: replace by ES host data
0376   auto modules = modules_.const_view<ModulesSoA>();
0377   auto module_nConnectedModules_buf_h = cms::alpakatools::make_host_buffer<uint16_t[]>(queue_, nLowerModules_);
0378   auto module_nConnectedModules_buf_d =
0379       cms::alpakatools::make_device_view(queue_, modules.nConnectedModules(), nLowerModules_);  // only lower modules
0380   alpaka::memcpy(queue_, module_nConnectedModules_buf_h, module_nConnectedModules_buf_d, nLowerModules_);
0381 
0382   alpaka::wait(queue_);  // wait for nSegments and module_nConnectedModules before using
0383 
0384   auto const* nSegments = nSegments_buf_h.data();
0385   auto const* module_nConnectedModules = module_nConnectedModules_buf_h.data();
0386 
0387   // Allocate host index and fill it directly
0388   auto index_buf_h = cms::alpakatools::make_host_buffer<uint16_t[]>(queue_, nLowerModules_);
0389   auto* index = index_buf_h.data();
0390 
0391   for (uint16_t innerLowerModuleIndex = 0; innerLowerModuleIndex < nLowerModules_; innerLowerModuleIndex++) {
0392     uint16_t nConnectedModules = module_nConnectedModules[innerLowerModuleIndex];
0393     unsigned int nInnerSegments = nSegments[innerLowerModuleIndex];
0394     if (nConnectedModules != 0 and nInnerSegments != 0) {
0395       index[nonZeroModules] = innerLowerModuleIndex;
0396       nonZeroModules++;
0397     }
0398     max_InnerSeg = std::max(max_InnerSeg, nInnerSegments);
0399   }
0400 
0401   // Allocate and copy to device index
0402   auto index_gpu_buf = cms::alpakatools::make_device_buffer<uint16_t[]>(queue_, nLowerModules_);
0403   alpaka::memcpy(queue_, index_gpu_buf, index_buf_h, nonZeroModules);
0404 
0405   auto const createTriplets_workDiv = cms::alpakatools::make_workdiv<Acc3D>({nonZeroModules, 1, 1}, {1, 16, 16});
0406 
0407   alpaka::exec<Acc3D>(queue_,
0408                       createTriplets_workDiv,
0409                       CreateTriplets{},
0410                       modules_.const_view<ModulesSoA>(),
0411                       miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0412                       segmentsDC_->const_view<SegmentsSoA>(),
0413                       segmentsDC_->const_view<SegmentsOccupancySoA>(),
0414                       tripletsDC_->view<TripletsSoA>(),
0415                       tripletsDC_->view<TripletsOccupancySoA>(),
0416                       rangesDC_->const_view(),
0417                       index_gpu_buf.data(),
0418                       nonZeroModules,
0419                       ptCut_);
0420 
0421   auto const addTripletRangesToEventExplicit_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0422 
0423   alpaka::exec<Acc1D>(queue_,
0424                       addTripletRangesToEventExplicit_workDiv,
0425                       AddTripletRangesToEventExplicit{},
0426                       modules_.const_view<ModulesSoA>(),
0427                       tripletsDC_->const_view<TripletsOccupancySoA>(),
0428                       rangesDC_->view());
0429 
0430   if (addObjects_) {
0431     addTripletsToEventExplicit();
0432   }
0433 }
0434 
0435 void LSTEvent::createTrackCandidates(bool no_pls_dupclean, bool tc_pls_triplets) {
0436   if (!trackCandidatesDC_) {
0437     trackCandidatesDC_.emplace(n_max_nonpixel_track_candidates + n_max_pixel_track_candidates, queue_);
0438     auto buf = trackCandidatesDC_->buffer();
0439     alpaka::memset(queue_, buf, 0u);
0440   }
0441 
0442   auto const crossCleanpT3_workDiv = cms::alpakatools::make_workdiv<Acc2D>({20, 4}, {64, 16});
0443 
0444   alpaka::exec<Acc2D>(queue_,
0445                       crossCleanpT3_workDiv,
0446                       CrossCleanpT3{},
0447                       modules_.const_view<ModulesSoA>(),
0448                       rangesDC_->const_view(),
0449                       pixelTripletsDC_->view(),
0450                       lstInputDC_->const_view<PixelSeedsSoA>(),
0451                       pixelQuintupletsDC_->const_view());
0452 
0453   auto const addpT3asTrackCandidates_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 512);
0454 
0455   alpaka::exec<Acc1D>(queue_,
0456                       addpT3asTrackCandidates_workDiv,
0457                       AddpT3asTrackCandidates{},
0458                       nLowerModules_,
0459                       pixelTripletsDC_->const_view(),
0460                       trackCandidatesDC_->view(),
0461                       lstInputDC_->const_view<PixelSeedsSoA>(),
0462                       rangesDC_->const_view());
0463 
0464   // Pull nEligibleT5Modules from the device.
0465   auto rangesOccupancy = rangesDC_->view();
0466   auto nEligibleModules_buf_h = cms::alpakatools::make_host_buffer<uint16_t>(queue_);
0467   auto nEligibleModules_buf_d = cms::alpakatools::make_device_view(queue_, rangesOccupancy.nEligibleT5Modules());
0468   alpaka::memcpy(queue_, nEligibleModules_buf_h, nEligibleModules_buf_d);
0469   alpaka::wait(queue_);  // wait to get the value before using
0470   auto const nEligibleModules = *nEligibleModules_buf_h.data();
0471 
0472   constexpr int threadsPerBlockY = 16;
0473   constexpr int threadsPerBlockX = 32;
0474   auto const removeDupQuintupletsBeforeTC_workDiv = cms::alpakatools::make_workdiv<Acc2D>(
0475       {std::max(nEligibleModules / threadsPerBlockY, 1), std::max(nEligibleModules / threadsPerBlockX, 1)}, {16, 32});
0476 
0477   alpaka::exec<Acc2D>(queue_,
0478                       removeDupQuintupletsBeforeTC_workDiv,
0479                       RemoveDupQuintupletsBeforeTC{},
0480                       quintupletsDC_->view<QuintupletsSoA>(),
0481                       quintupletsDC_->view<QuintupletsOccupancySoA>(),
0482                       rangesDC_->const_view());
0483 
0484   constexpr int threadsPerBlock = 32;
0485   auto const crossCleanT5_workDiv = cms::alpakatools::make_workdiv<Acc3D>(
0486       {(nLowerModules_ / threadsPerBlock) + 1, 1, max_blocks}, {threadsPerBlock, 1, threadsPerBlock});
0487 
0488   alpaka::exec<Acc3D>(queue_,
0489                       crossCleanT5_workDiv,
0490                       CrossCleanT5{},
0491                       modules_.const_view<ModulesSoA>(),
0492                       quintupletsDC_->view<QuintupletsSoA>(),
0493                       quintupletsDC_->const_view<QuintupletsOccupancySoA>(),
0494                       pixelQuintupletsDC_->const_view(),
0495                       pixelTripletsDC_->const_view(),
0496                       rangesDC_->const_view());
0497 
0498   auto const addT5asTrackCandidate_workDiv = cms::alpakatools::make_workdiv<Acc2D>({8, 10}, {8, 128});
0499 
0500   alpaka::exec<Acc2D>(queue_,
0501                       addT5asTrackCandidate_workDiv,
0502                       AddT5asTrackCandidate{},
0503                       nLowerModules_,
0504                       quintupletsDC_->const_view<QuintupletsSoA>(),
0505                       quintupletsDC_->const_view<QuintupletsOccupancySoA>(),
0506                       trackCandidatesDC_->view(),
0507                       rangesDC_->const_view());
0508 
0509   if (!no_pls_dupclean) {
0510     auto const checkHitspLS_workDiv = cms::alpakatools::make_workdiv<Acc2D>({max_blocks * 4, max_blocks / 4}, {16, 16});
0511 
0512     alpaka::exec<Acc2D>(queue_,
0513                         checkHitspLS_workDiv,
0514                         CheckHitspLS{},
0515                         modules_.const_view<ModulesSoA>(),
0516                         segmentsDC_->const_view<SegmentsOccupancySoA>(),
0517                         lstInputDC_->const_view<PixelSeedsSoA>(),
0518                         pixelSegmentsDC_->view(),
0519                         true);
0520   }
0521 
0522   auto const crossCleanpLS_workDiv = cms::alpakatools::make_workdiv<Acc2D>({20, 4}, {32, 16});
0523 
0524   alpaka::exec<Acc2D>(queue_,
0525                       crossCleanpLS_workDiv,
0526                       CrossCleanpLS{},
0527                       modules_.const_view<ModulesSoA>(),
0528                       rangesDC_->const_view(),
0529                       pixelTripletsDC_->const_view(),
0530                       trackCandidatesDC_->view(),
0531                       segmentsDC_->const_view<SegmentsSoA>(),
0532                       segmentsDC_->const_view<SegmentsOccupancySoA>(),
0533                       lstInputDC_->const_view<PixelSeedsSoA>(),
0534                       pixelSegmentsDC_->view(),
0535                       miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0536                       lstInputDC_->const_view<HitsBaseSoA>(),
0537                       quintupletsDC_->const_view<QuintupletsSoA>());
0538 
0539   auto const addpLSasTrackCandidate_workDiv = cms::alpakatools::make_workdiv<Acc1D>(max_blocks, 384);
0540 
0541   alpaka::exec<Acc1D>(queue_,
0542                       addpLSasTrackCandidate_workDiv,
0543                       AddpLSasTrackCandidate{},
0544                       nLowerModules_,
0545                       trackCandidatesDC_->view(),
0546                       segmentsDC_->const_view<SegmentsOccupancySoA>(),
0547                       lstInputDC_->const_view<PixelSeedsSoA>(),
0548                       pixelSegmentsDC_->const_view(),
0549                       tc_pls_triplets);
0550 
0551   // Check if either n_max_pixel_track_candidates or n_max_nonpixel_track_candidates was reached
0552   auto nTrackCanpT5Host_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0553   auto nTrackCanpT3Host_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0554   auto nTrackCanpLSHost_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0555   auto nTrackCanT5Host_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0556   alpaka::memcpy(queue_,
0557                  nTrackCanpT5Host_buf,
0558                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatespT5()));
0559   alpaka::memcpy(queue_,
0560                  nTrackCanpT3Host_buf,
0561                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatespT3()));
0562   alpaka::memcpy(queue_,
0563                  nTrackCanpLSHost_buf,
0564                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatespLS()));
0565   alpaka::memcpy(queue_,
0566                  nTrackCanT5Host_buf,
0567                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatesT5()));
0568   alpaka::wait(queue_);  // wait to get the values before using them
0569 
0570   auto nTrackCandidatespT5 = *nTrackCanpT5Host_buf.data();
0571   auto nTrackCandidatespT3 = *nTrackCanpT3Host_buf.data();
0572   auto nTrackCandidatespLS = *nTrackCanpLSHost_buf.data();
0573   auto nTrackCandidatesT5 = *nTrackCanT5Host_buf.data();
0574   if ((nTrackCandidatespT5 + nTrackCandidatespT3 + nTrackCandidatespLS == n_max_pixel_track_candidates) ||
0575       (nTrackCandidatesT5 == n_max_nonpixel_track_candidates)) {
0576     lstWarning(
0577         "\
0578         ****************************************************************************************************\n\
0579         * Track candidates were possibly truncated.                                                        *\n\
0580         * You may need to increase either n_max_pixel_track_candidates or n_max_nonpixel_track_candidates. *\n\
0581         * Run the code with the WARNINGS flag activated for more details.                                  *\n\
0582         ****************************************************************************************************");
0583   }
0584 }
0585 
0586 void LSTEvent::createPixelTriplets() {
0587   if (!pixelTripletsDC_) {
0588     pixelTripletsDC_.emplace(n_max_pixel_triplets, queue_);
0589     auto nPixelTriplets_view = cms::alpakatools::make_device_view(queue_, (*pixelTripletsDC_)->nPixelTriplets());
0590     alpaka::memset(queue_, nPixelTriplets_view, 0u);
0591     auto totOccupancyPixelTriplets_view =
0592         cms::alpakatools::make_device_view(queue_, (*pixelTripletsDC_)->totOccupancyPixelTriplets());
0593     alpaka::memset(queue_, totOccupancyPixelTriplets_view, 0u);
0594   }
0595   SegmentsOccupancy segmentsOccupancy = segmentsDC_->view<SegmentsOccupancySoA>();
0596   PixelSeedsConst pixelSeeds = lstInputDC_->const_view<PixelSeedsSoA>();
0597 
0598   auto superbins_buf = cms::alpakatools::make_host_buffer<int[]>(queue_, n_max_pixel_segments_per_module);
0599   auto pixelTypes_buf = cms::alpakatools::make_host_buffer<PixelType[]>(queue_, n_max_pixel_segments_per_module);
0600 
0601   alpaka::memcpy(queue_, superbins_buf, cms::alpakatools::make_device_view(queue_, pixelSeeds.superbin(), pixelSize_));
0602   alpaka::memcpy(
0603       queue_, pixelTypes_buf, cms::alpakatools::make_device_view(queue_, pixelSeeds.pixelType(), pixelSize_));
0604   auto const* superbins = superbins_buf.data();
0605   auto const* pixelTypes = pixelTypes_buf.data();
0606 
0607   unsigned int nInnerSegments;
0608   auto nInnerSegments_src_view = cms::alpakatools::make_host_view(nInnerSegments);
0609 
0610   // Create a sub-view for the device buffer
0611   auto dev_view_nSegments = cms::alpakatools::make_device_view(queue_, segmentsOccupancy.nSegments()[nLowerModules_]);
0612 
0613   alpaka::memcpy(queue_, nInnerSegments_src_view, dev_view_nSegments);
0614   alpaka::wait(queue_);  // wait to get nInnerSegments (also superbins and pixelTypes) before using
0615 
0616   auto connectedPixelSize_host_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nInnerSegments);
0617   auto connectedPixelIndex_host_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nInnerSegments);
0618   auto connectedPixelSize_dev_buf = cms::alpakatools::make_device_buffer<unsigned int[]>(queue_, nInnerSegments);
0619   auto connectedPixelIndex_dev_buf = cms::alpakatools::make_device_buffer<unsigned int[]>(queue_, nInnerSegments);
0620 
0621   unsigned int* connectedPixelSize_host = connectedPixelSize_host_buf.data();
0622   unsigned int* connectedPixelIndex_host = connectedPixelIndex_host_buf.data();
0623 
0624   int pixelIndexOffsetPos =
0625       pixelMapping_.connectedPixelsIndex[size_superbins - 1] + pixelMapping_.connectedPixelsSizes[size_superbins - 1];
0626   int pixelIndexOffsetNeg = pixelMapping_.connectedPixelsIndexPos[size_superbins - 1] +
0627                             pixelMapping_.connectedPixelsSizesPos[size_superbins - 1] + pixelIndexOffsetPos;
0628 
0629   // TODO: check if a map/reduction to just eligible pLSs would speed up the kernel
0630   // the current selection still leaves a significant fraction of unmatchable pLSs
0631   for (unsigned int i = 0; i < nInnerSegments; i++) {  // loop over # pLS
0632     PixelType pixelType = pixelTypes[i];               // Get pixel type for this pLS
0633     int superbin = superbins[i];                       // Get superbin for this pixel
0634     if ((superbin < 0) or (superbin >= (int)size_superbins) or
0635         ((pixelType != PixelType::kHighPt) and (pixelType != PixelType::kLowPtPosCurv) and
0636          (pixelType != PixelType::kLowPtNegCurv))) {
0637       connectedPixelSize_host[i] = 0;
0638       connectedPixelIndex_host[i] = 0;
0639       continue;
0640     }
0641 
0642     // Used pixel type to select correct size-index arrays
0643     switch (pixelType) {
0644       case PixelType::kInvalid:
0645         break;
0646       case PixelType::kHighPt:
0647         // number of connected modules to this pixel
0648         connectedPixelSize_host[i] = pixelMapping_.connectedPixelsSizes[superbin];
0649         // index to get start of connected modules for this superbin in map
0650         connectedPixelIndex_host[i] = pixelMapping_.connectedPixelsIndex[superbin];
0651         break;
0652       case PixelType::kLowPtPosCurv:
0653         // number of connected modules to this pixel
0654         connectedPixelSize_host[i] = pixelMapping_.connectedPixelsSizesPos[superbin];
0655         // index to get start of connected modules for this superbin in map
0656         connectedPixelIndex_host[i] = pixelMapping_.connectedPixelsIndexPos[superbin] + pixelIndexOffsetPos;
0657         break;
0658       case PixelType::kLowPtNegCurv:
0659         // number of connected modules to this pixel
0660         connectedPixelSize_host[i] = pixelMapping_.connectedPixelsSizesNeg[superbin];
0661         // index to get start of connected modules for this superbin in map
0662         connectedPixelIndex_host[i] = pixelMapping_.connectedPixelsIndexNeg[superbin] + pixelIndexOffsetNeg;
0663         break;
0664     }
0665   }
0666 
0667   alpaka::memcpy(queue_, connectedPixelSize_dev_buf, connectedPixelSize_host_buf, nInnerSegments);
0668   alpaka::memcpy(queue_, connectedPixelIndex_dev_buf, connectedPixelIndex_host_buf, nInnerSegments);
0669 
0670   auto const createPixelTripletsFromMap_workDiv =
0671       cms::alpakatools::make_workdiv<Acc3D>({4096, 16 /* above median of connected modules*/, 1}, {4, 1, 32});
0672 
0673   alpaka::exec<Acc3D>(queue_,
0674                       createPixelTripletsFromMap_workDiv,
0675                       CreatePixelTripletsFromMap{},
0676                       modules_.const_view<ModulesSoA>(),
0677                       modules_.const_view<ModulesPixelSoA>(),
0678                       rangesDC_->const_view(),
0679                       miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0680                       segmentsDC_->const_view<SegmentsSoA>(),
0681                       lstInputDC_->const_view<PixelSeedsSoA>(),
0682                       pixelSegmentsDC_->const_view(),
0683                       tripletsDC_->view<TripletsSoA>(),
0684                       tripletsDC_->const_view<TripletsOccupancySoA>(),
0685                       pixelTripletsDC_->view(),
0686                       connectedPixelSize_dev_buf.data(),
0687                       connectedPixelIndex_dev_buf.data(),
0688                       nInnerSegments,
0689                       ptCut_);
0690 
0691 #ifdef WARNINGS
0692   auto nPixelTriplets_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0693 
0694   alpaka::memcpy(
0695       queue_, nPixelTriplets_buf, cms::alpakatools::make_device_view(queue_, (*pixelTripletsDC_)->nPixelTriplets()));
0696   alpaka::wait(queue_);  // wait to get the value before using it
0697 
0698   std::cout << "number of pixel triplets = " << *nPixelTriplets_buf.data() << std::endl;
0699 #endif
0700 
0701   //pT3s can be cleaned here because they're not used in making pT5s!
0702   //seems like more blocks lead to conflicting writes
0703   auto const removeDupPixelTripletsFromMap_workDiv = cms::alpakatools::make_workdiv<Acc2D>({40, 1}, {16, 16});
0704 
0705   alpaka::exec<Acc2D>(
0706       queue_, removeDupPixelTripletsFromMap_workDiv, RemoveDupPixelTripletsFromMap{}, pixelTripletsDC_->view());
0707 }
0708 
0709 void LSTEvent::createQuintuplets() {
0710   auto const createEligibleModulesListForQuintuplets_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0711 
0712   alpaka::exec<Acc1D>(queue_,
0713                       createEligibleModulesListForQuintuplets_workDiv,
0714                       CreateEligibleModulesListForQuintuplets{},
0715                       modules_.const_view<ModulesSoA>(),
0716                       tripletsDC_->const_view<TripletsOccupancySoA>(),
0717                       rangesDC_->view(),
0718                       tripletsDC_->view<TripletsSoA>(),
0719                       ptCut_);
0720 
0721   auto nEligibleT5Modules_buf = cms::alpakatools::make_host_buffer<uint16_t>(queue_);
0722   auto nTotalQuintuplets_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0723   auto rangesOccupancy = rangesDC_->view();
0724   auto nEligibleT5Modules_view_d = cms::alpakatools::make_device_view(queue_, rangesOccupancy.nEligibleT5Modules());
0725   auto nTotalQuintuplets_view_d = cms::alpakatools::make_device_view(queue_, rangesOccupancy.nTotalQuints());
0726   alpaka::memcpy(queue_, nEligibleT5Modules_buf, nEligibleT5Modules_view_d);
0727   alpaka::memcpy(queue_, nTotalQuintuplets_buf, nTotalQuintuplets_view_d);
0728   alpaka::wait(queue_);  // wait for the values before using them
0729 
0730   auto nEligibleT5Modules = *nEligibleT5Modules_buf.data();
0731   auto nTotalQuintuplets = *nTotalQuintuplets_buf.data();
0732 
0733   if (!quintupletsDC_) {
0734     std::array<int, 2> const quintuplets_sizes{{static_cast<int>(nTotalQuintuplets), static_cast<int>(nLowerModules_)}};
0735     quintupletsDC_.emplace(quintuplets_sizes, queue_);
0736     auto quintupletsOccupancy = quintupletsDC_->view<QuintupletsOccupancySoA>();
0737     auto nQuintuplets_view = cms::alpakatools::make_device_view(
0738         queue_, quintupletsOccupancy.nQuintuplets(), quintupletsOccupancy.metadata().size());
0739     alpaka::memset(queue_, nQuintuplets_view, 0u);
0740     auto totOccupancyQuintuplets_view = cms::alpakatools::make_device_view(
0741         queue_, quintupletsOccupancy.totOccupancyQuintuplets(), quintupletsOccupancy.metadata().size());
0742     alpaka::memset(queue_, totOccupancyQuintuplets_view, 0u);
0743     auto quintuplets = quintupletsDC_->view<QuintupletsSoA>();
0744     auto isDup_view = cms::alpakatools::make_device_view(queue_, quintuplets.isDup(), quintuplets.metadata().size());
0745     alpaka::memset(queue_, isDup_view, 0u);
0746     auto tightCutFlag_view =
0747         cms::alpakatools::make_device_view(queue_, quintuplets.tightCutFlag(), quintuplets.metadata().size());
0748     alpaka::memset(queue_, tightCutFlag_view, 0u);
0749     auto partOfPT5_view =
0750         cms::alpakatools::make_device_view(queue_, quintuplets.partOfPT5(), quintuplets.metadata().size());
0751     alpaka::memset(queue_, partOfPT5_view, 0u);
0752   }
0753 
0754   auto const createQuintuplets_workDiv =
0755       cms::alpakatools::make_workdiv<Acc3D>({std::max((int)nEligibleT5Modules, 1), 1, 1}, {1, 8, 32});
0756 
0757   alpaka::exec<Acc3D>(queue_,
0758                       createQuintuplets_workDiv,
0759                       CreateQuintuplets{},
0760                       modules_.const_view<ModulesSoA>(),
0761                       miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0762                       segmentsDC_->const_view<SegmentsSoA>(),
0763                       tripletsDC_->view<TripletsSoA>(),
0764                       tripletsDC_->const_view<TripletsOccupancySoA>(),
0765                       quintupletsDC_->view<QuintupletsSoA>(),
0766                       quintupletsDC_->view<QuintupletsOccupancySoA>(),
0767                       rangesDC_->const_view(),
0768                       nEligibleT5Modules,
0769                       ptCut_);
0770 
0771   auto const removeDupQuintupletsAfterBuild_workDiv =
0772       cms::alpakatools::make_workdiv<Acc3D>({max_blocks, 1, 1}, {1, 16, 16});
0773 
0774   alpaka::exec<Acc3D>(queue_,
0775                       removeDupQuintupletsAfterBuild_workDiv,
0776                       RemoveDupQuintupletsAfterBuild{},
0777                       modules_.const_view<ModulesSoA>(),
0778                       quintupletsDC_->view<QuintupletsSoA>(),
0779                       quintupletsDC_->const_view<QuintupletsOccupancySoA>(),
0780                       rangesDC_->const_view());
0781 
0782   auto const addQuintupletRangesToEventExplicit_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 1024);
0783 
0784   alpaka::exec<Acc1D>(queue_,
0785                       addQuintupletRangesToEventExplicit_workDiv,
0786                       AddQuintupletRangesToEventExplicit{},
0787                       modules_.const_view<ModulesSoA>(),
0788                       quintupletsDC_->const_view<QuintupletsOccupancySoA>(),
0789                       rangesDC_->view());
0790 
0791   if (addObjects_) {
0792     addQuintupletsToEventExplicit();
0793   }
0794 }
0795 
0796 void LSTEvent::pixelLineSegmentCleaning(bool no_pls_dupclean) {
0797   if (!no_pls_dupclean) {
0798     auto const checkHitspLS_workDiv = cms::alpakatools::make_workdiv<Acc2D>({max_blocks * 4, max_blocks / 4}, {16, 16});
0799 
0800     alpaka::exec<Acc2D>(queue_,
0801                         checkHitspLS_workDiv,
0802                         CheckHitspLS{},
0803                         modules_.const_view<ModulesSoA>(),
0804                         segmentsDC_->const_view<SegmentsOccupancySoA>(),
0805                         lstInputDC_->const_view<PixelSeedsSoA>(),
0806                         pixelSegmentsDC_->view(),
0807                         false);
0808   }
0809 }
0810 
0811 void LSTEvent::createPixelQuintuplets() {
0812   if (!pixelQuintupletsDC_) {
0813     pixelQuintupletsDC_.emplace(n_max_pixel_quintuplets, queue_);
0814     auto nPixelQuintuplets_view =
0815         cms::alpakatools::make_device_view(queue_, (*pixelQuintupletsDC_)->nPixelQuintuplets());
0816     alpaka::memset(queue_, nPixelQuintuplets_view, 0u);
0817     auto totOccupancyPixelQuintuplets_view =
0818         cms::alpakatools::make_device_view(queue_, (*pixelQuintupletsDC_)->totOccupancyPixelQuintuplets());
0819     alpaka::memset(queue_, totOccupancyPixelQuintuplets_view, 0u);
0820   }
0821   if (!trackCandidatesDC_) {
0822     trackCandidatesDC_.emplace(n_max_nonpixel_track_candidates + n_max_pixel_track_candidates, queue_);
0823     auto buf = trackCandidatesDC_->buffer();
0824     alpaka::memset(queue_, buf, 0u);
0825   }
0826   SegmentsOccupancy segmentsOccupancy = segmentsDC_->view<SegmentsOccupancySoA>();
0827   PixelSeedsConst pixelSeeds = lstInputDC_->const_view<PixelSeedsSoA>();
0828 
0829   auto superbins_buf = cms::alpakatools::make_host_buffer<int[]>(queue_, n_max_pixel_segments_per_module);
0830   auto pixelTypes_buf = cms::alpakatools::make_host_buffer<PixelType[]>(queue_, n_max_pixel_segments_per_module);
0831 
0832   alpaka::memcpy(queue_, superbins_buf, cms::alpakatools::make_device_view(queue_, pixelSeeds.superbin(), pixelSize_));
0833   alpaka::memcpy(
0834       queue_, pixelTypes_buf, cms::alpakatools::make_device_view(queue_, pixelSeeds.pixelType(), pixelSize_));
0835   auto const* superbins = superbins_buf.data();
0836   auto const* pixelTypes = pixelTypes_buf.data();
0837 
0838   unsigned int nInnerSegments;
0839   auto nInnerSegments_src_view = cms::alpakatools::make_host_view(nInnerSegments);
0840 
0841   // Create a sub-view for the device buffer
0842   unsigned int totalModules = nLowerModules_ + 1;
0843   auto dev_view_nSegments_buf = cms::alpakatools::make_device_view(queue_, segmentsOccupancy.nSegments(), totalModules);
0844   auto dev_view_nSegments = cms::alpakatools::make_device_view(queue_, segmentsOccupancy.nSegments()[nLowerModules_]);
0845 
0846   alpaka::memcpy(queue_, nInnerSegments_src_view, dev_view_nSegments);
0847   alpaka::wait(queue_);  // wait to get nInnerSegments (also superbins and pixelTypes) before using
0848 
0849   auto connectedPixelSize_host_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nInnerSegments);
0850   auto connectedPixelIndex_host_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nInnerSegments);
0851   auto connectedPixelSize_dev_buf = cms::alpakatools::make_device_buffer<unsigned int[]>(queue_, nInnerSegments);
0852   auto connectedPixelIndex_dev_buf = cms::alpakatools::make_device_buffer<unsigned int[]>(queue_, nInnerSegments);
0853 
0854   auto* connectedPixelSize_host = connectedPixelSize_host_buf.data();
0855   auto* connectedPixelIndex_host = connectedPixelIndex_host_buf.data();
0856 
0857   int pixelIndexOffsetPos = pixelMapping_.connectedPixelsIndex[::size_superbins - 1] +
0858                             pixelMapping_.connectedPixelsSizes[::size_superbins - 1];
0859   int pixelIndexOffsetNeg = pixelMapping_.connectedPixelsIndexPos[::size_superbins - 1] +
0860                             pixelMapping_.connectedPixelsSizesPos[::size_superbins - 1] + pixelIndexOffsetPos;
0861 
0862   // Loop over # pLS
0863   for (unsigned int i = 0; i < nInnerSegments; i++) {
0864     PixelType pixelType = pixelTypes[i];  // Get pixel type for this pLS
0865     int superbin = superbins[i];          // Get superbin for this pixel
0866     if ((superbin < 0) or (superbin >= (int)size_superbins) or
0867         ((pixelType != PixelType::kHighPt) and (pixelType != PixelType::kLowPtPosCurv) and
0868          (pixelType != PixelType::kLowPtNegCurv))) {
0869       connectedPixelSize_host[i] = 0;
0870       connectedPixelIndex_host[i] = 0;
0871       continue;
0872     }
0873 
0874     // Used pixel type to select correct size-index arrays
0875     switch (pixelType) {
0876       case PixelType::kInvalid:
0877         break;
0878       case PixelType::kHighPt:
0879         // number of connected modules to this pixel
0880         connectedPixelSize_host[i] = pixelMapping_.connectedPixelsSizes[superbin];
0881         // index to get start of connected modules for this superbin in map
0882         connectedPixelIndex_host[i] = pixelMapping_.connectedPixelsIndex[superbin];
0883         break;
0884       case PixelType::kLowPtPosCurv:
0885         // number of connected modules to this pixel
0886         connectedPixelSize_host[i] = pixelMapping_.connectedPixelsSizesPos[superbin];
0887         // index to get start of connected modules for this superbin in map
0888         connectedPixelIndex_host[i] = pixelMapping_.connectedPixelsIndexPos[superbin] + pixelIndexOffsetPos;
0889         break;
0890       case PixelType::kLowPtNegCurv:
0891         // number of connected modules to this pixel
0892         connectedPixelSize_host[i] = pixelMapping_.connectedPixelsSizesNeg[superbin];
0893         // index to get start of connected modules for this superbin in map
0894         connectedPixelIndex_host[i] = pixelMapping_.connectedPixelsIndexNeg[superbin] + pixelIndexOffsetNeg;
0895         break;
0896     }
0897   }
0898 
0899   alpaka::memcpy(queue_, connectedPixelSize_dev_buf, connectedPixelSize_host_buf, nInnerSegments);
0900   alpaka::memcpy(queue_, connectedPixelIndex_dev_buf, connectedPixelIndex_host_buf, nInnerSegments);
0901 
0902   auto const createPixelQuintupletsFromMap_workDiv =
0903       cms::alpakatools::make_workdiv<Acc3D>({max_blocks, 16, 1}, {16, 1, 16});
0904 
0905   alpaka::exec<Acc3D>(queue_,
0906                       createPixelQuintupletsFromMap_workDiv,
0907                       CreatePixelQuintupletsFromMap{},
0908                       modules_.const_view<ModulesSoA>(),
0909                       modules_.const_view<ModulesPixelSoA>(),
0910                       miniDoubletsDC_->const_view<MiniDoubletsSoA>(),
0911                       segmentsDC_->const_view<SegmentsSoA>(),
0912                       lstInputDC_->const_view<PixelSeedsSoA>(),
0913                       pixelSegmentsDC_->view(),
0914                       tripletsDC_->view<TripletsSoA>(),
0915                       quintupletsDC_->view<QuintupletsSoA>(),
0916                       quintupletsDC_->const_view<QuintupletsOccupancySoA>(),
0917                       pixelQuintupletsDC_->view(),
0918                       connectedPixelSize_dev_buf.data(),
0919                       connectedPixelIndex_dev_buf.data(),
0920                       nInnerSegments,
0921                       rangesDC_->const_view(),
0922                       ptCut_);
0923 
0924   auto const removeDupPixelQuintupletsFromMap_workDiv =
0925       cms::alpakatools::make_workdiv<Acc2D>({max_blocks, 1}, {16, 16});
0926 
0927   alpaka::exec<Acc2D>(queue_,
0928                       removeDupPixelQuintupletsFromMap_workDiv,
0929                       RemoveDupPixelQuintupletsFromMap{},
0930                       pixelQuintupletsDC_->view());
0931 
0932   auto const addpT5asTrackCandidate_workDiv = cms::alpakatools::make_workdiv<Acc1D>(1, 256);
0933 
0934   alpaka::exec<Acc1D>(queue_,
0935                       addpT5asTrackCandidate_workDiv,
0936                       AddpT5asTrackCandidate{},
0937                       nLowerModules_,
0938                       pixelQuintupletsDC_->const_view(),
0939                       trackCandidatesDC_->view(),
0940                       lstInputDC_->const_view<PixelSeedsSoA>(),
0941                       rangesDC_->const_view());
0942 
0943 #ifdef WARNINGS
0944   auto nPixelQuintuplets_buf = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
0945 
0946   alpaka::memcpy(queue_,
0947                  nPixelQuintuplets_buf,
0948                  cms::alpakatools::make_device_view(queue_, (*pixelQuintupletsDC_)->nPixelQuintuplets()));
0949   alpaka::wait(queue_);  // wait to get the value before using it
0950 
0951   std::cout << "number of pixel quintuplets = " << *nPixelQuintuplets_buf.data() << std::endl;
0952 #endif
0953 }
0954 
0955 void LSTEvent::addMiniDoubletsToEventExplicit() {
0956   auto nMDsCPU_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nLowerModules_);
0957   auto mdsOccupancy = miniDoubletsDC_->const_view<MiniDoubletsOccupancySoA>();
0958   auto nMDs_view =
0959       cms::alpakatools::make_device_view(queue_, mdsOccupancy.nMDs(), nLowerModules_);  // exclude pixel part
0960   alpaka::memcpy(queue_, nMDsCPU_buf, nMDs_view, nLowerModules_);
0961 
0962   auto modules = modules_.const_view<ModulesSoA>();
0963 
0964   // FIXME: replace by ES host data
0965   auto module_subdets_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
0966   auto module_subdets_view =
0967       cms::alpakatools::make_device_view(queue_, modules.subdets(), nLowerModules_);  // only lower modules
0968   alpaka::memcpy(queue_, module_subdets_buf, module_subdets_view, nLowerModules_);
0969 
0970   auto module_layers_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
0971   auto module_layers_view =
0972       cms::alpakatools::make_device_view(queue_, modules.layers(), nLowerModules_);  // only lower modules
0973   alpaka::memcpy(queue_, module_layers_buf, module_layers_view, nLowerModules_);
0974 
0975   alpaka::wait(queue_);  // wait for inputs before using them
0976 
0977   auto const* nMDsCPU = nMDsCPU_buf.data();
0978   auto const* module_subdets = module_subdets_buf.data();
0979   auto const* module_layers = module_layers_buf.data();
0980 
0981   for (unsigned int i = 0; i < nLowerModules_; i++) {
0982     if (nMDsCPU[i] != 0) {
0983       if (module_subdets[i] == Barrel) {
0984         n_minidoublets_by_layer_barrel_[module_layers[i] - 1] += nMDsCPU[i];
0985       } else {
0986         n_minidoublets_by_layer_endcap_[module_layers[i] - 1] += nMDsCPU[i];
0987       }
0988     }
0989   }
0990 }
0991 
0992 void LSTEvent::addSegmentsToEventExplicit() {
0993   auto nSegmentsCPU_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nLowerModules_);
0994   auto nSegments_buf = cms::alpakatools::make_device_view(
0995       queue_, segmentsDC_->const_view<SegmentsOccupancySoA>().nSegments(), nLowerModules_);
0996   alpaka::memcpy(queue_, nSegmentsCPU_buf, nSegments_buf, nLowerModules_);
0997 
0998   auto modules = modules_.const_view<ModulesSoA>();
0999 
1000   // FIXME: replace by ES host data
1001   auto module_subdets_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
1002   auto module_subdets_view =
1003       cms::alpakatools::make_device_view(queue_, modules.subdets(), nLowerModules_);  // only lower modules
1004   alpaka::memcpy(queue_, module_subdets_buf, module_subdets_view, nLowerModules_);
1005 
1006   auto module_layers_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
1007   auto module_layers_view =
1008       cms::alpakatools::make_device_view(queue_, modules.layers(), nLowerModules_);  // only lower modules
1009   alpaka::memcpy(queue_, module_layers_buf, module_layers_view, nLowerModules_);
1010 
1011   alpaka::wait(queue_);  // wait for inputs before using them
1012 
1013   auto const* nSegmentsCPU = nSegmentsCPU_buf.data();
1014   auto const* module_subdets = module_subdets_buf.data();
1015   auto const* module_layers = module_layers_buf.data();
1016 
1017   for (unsigned int i = 0; i < nLowerModules_; i++) {
1018     if (!(nSegmentsCPU[i] == 0)) {
1019       if (module_subdets[i] == Barrel) {
1020         n_segments_by_layer_barrel_[module_layers[i] - 1] += nSegmentsCPU[i];
1021       } else {
1022         n_segments_by_layer_endcap_[module_layers[i] - 1] += nSegmentsCPU[i];
1023       }
1024     }
1025   }
1026 }
1027 
1028 void LSTEvent::addQuintupletsToEventExplicit() {
1029   auto quintupletsOccupancy = quintupletsDC_->const_view<QuintupletsOccupancySoA>();
1030   auto nQuintuplets_view =
1031       cms::alpakatools::make_device_view(queue_, quintupletsOccupancy.nQuintuplets(), nLowerModules_);
1032   auto nQuintupletsCPU_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nLowerModules_);
1033   alpaka::memcpy(queue_, nQuintupletsCPU_buf, nQuintuplets_view);
1034 
1035   auto modules = modules_.const_view<ModulesSoA>();
1036 
1037   // FIXME: replace by ES host data
1038   auto module_subdets_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
1039   auto module_subdets_view =
1040       cms::alpakatools::make_device_view(queue_, modules.subdets(), nLowerModules_);  // only lower modules
1041   alpaka::memcpy(queue_, module_subdets_buf, module_subdets_view, nLowerModules_);
1042 
1043   auto module_layers_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
1044   auto module_layers_view =
1045       cms::alpakatools::make_device_view(queue_, modules.layers(), nLowerModules_);  // only lower modules
1046   alpaka::memcpy(queue_, module_layers_buf, module_layers_view, nLowerModules_);
1047 
1048   auto module_quintupletModuleIndices_buf = cms::alpakatools::make_host_buffer<int[]>(queue_, nLowerModules_);
1049   auto rangesOccupancy = rangesDC_->view();
1050   auto quintupletModuleIndices_view_d =
1051       cms::alpakatools::make_device_view(queue_, rangesOccupancy.quintupletModuleIndices(), nLowerModules_);
1052   alpaka::memcpy(queue_, module_quintupletModuleIndices_buf, quintupletModuleIndices_view_d);
1053 
1054   alpaka::wait(queue_);  // wait for inputs before using them
1055 
1056   auto const* nQuintupletsCPU = nQuintupletsCPU_buf.data();
1057   auto const* module_subdets = module_subdets_buf.data();
1058   auto const* module_layers = module_layers_buf.data();
1059   auto const* module_quintupletModuleIndices = module_quintupletModuleIndices_buf.data();
1060 
1061   for (uint16_t i = 0; i < nLowerModules_; i++) {
1062     if (!(nQuintupletsCPU[i] == 0 or module_quintupletModuleIndices[i] == -1)) {
1063       if (module_subdets[i] == Barrel) {
1064         n_quintuplets_by_layer_barrel_[module_layers[i] - 1] += nQuintupletsCPU[i];
1065       } else {
1066         n_quintuplets_by_layer_endcap_[module_layers[i] - 1] += nQuintupletsCPU[i];
1067       }
1068     }
1069   }
1070 }
1071 
1072 void LSTEvent::addTripletsToEventExplicit() {
1073   auto tripletsOccupancy = tripletsDC_->const_view<TripletsOccupancySoA>();
1074   auto nTriplets_view = cms::alpakatools::make_device_view(queue_, tripletsOccupancy.nTriplets(), nLowerModules_);
1075   auto nTripletsCPU_buf = cms::alpakatools::make_host_buffer<unsigned int[]>(queue_, nLowerModules_);
1076   alpaka::memcpy(queue_, nTripletsCPU_buf, nTriplets_view);
1077 
1078   auto modules = modules_.const_view<ModulesSoA>();
1079 
1080   // FIXME: replace by ES host data
1081   auto module_subdets_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
1082   auto module_subdets_view =
1083       cms::alpakatools::make_device_view(queue_, modules.subdets(), nLowerModules_);  // only lower modules
1084   alpaka::memcpy(queue_, module_subdets_buf, module_subdets_view, nLowerModules_);
1085 
1086   auto module_layers_buf = cms::alpakatools::make_host_buffer<short[]>(queue_, nLowerModules_);
1087   auto module_layers_view =
1088       cms::alpakatools::make_device_view(queue_, modules.layers(), nLowerModules_);  // only lower modules
1089   alpaka::memcpy(queue_, module_layers_buf, module_layers_view, nLowerModules_);
1090 
1091   alpaka::wait(queue_);  // wait for inputs before using them
1092 
1093   auto const* nTripletsCPU = nTripletsCPU_buf.data();
1094   auto const* module_subdets = module_subdets_buf.data();
1095   auto const* module_layers = module_layers_buf.data();
1096 
1097   for (uint16_t i = 0; i < nLowerModules_; i++) {
1098     if (nTripletsCPU[i] != 0) {
1099       if (module_subdets[i] == Barrel) {
1100         n_triplets_by_layer_barrel_[module_layers[i] - 1] += nTripletsCPU[i];
1101       } else {
1102         n_triplets_by_layer_endcap_[module_layers[i] - 1] += nTripletsCPU[i];
1103       }
1104     }
1105   }
1106 }
1107 
1108 unsigned int LSTEvent::getNumberOfMiniDoublets() {
1109   unsigned int miniDoublets = 0;
1110   for (auto& it : n_minidoublets_by_layer_barrel_) {
1111     miniDoublets += it;
1112   }
1113   for (auto& it : n_minidoublets_by_layer_endcap_) {
1114     miniDoublets += it;
1115   }
1116 
1117   return miniDoublets;
1118 }
1119 
1120 unsigned int LSTEvent::getNumberOfMiniDoubletsByLayerBarrel(unsigned int layer) {
1121   return n_minidoublets_by_layer_barrel_[layer];
1122 }
1123 
1124 unsigned int LSTEvent::getNumberOfMiniDoubletsByLayerEndcap(unsigned int layer) {
1125   return n_minidoublets_by_layer_endcap_[layer];
1126 }
1127 
1128 unsigned int LSTEvent::getNumberOfSegments() {
1129   unsigned int segments = 0;
1130   for (auto& it : n_segments_by_layer_barrel_) {
1131     segments += it;
1132   }
1133   for (auto& it : n_segments_by_layer_endcap_) {
1134     segments += it;
1135   }
1136 
1137   return segments;
1138 }
1139 
1140 unsigned int LSTEvent::getNumberOfSegmentsByLayerBarrel(unsigned int layer) {
1141   return n_segments_by_layer_barrel_[layer];
1142 }
1143 
1144 unsigned int LSTEvent::getNumberOfSegmentsByLayerEndcap(unsigned int layer) {
1145   return n_segments_by_layer_endcap_[layer];
1146 }
1147 
1148 unsigned int LSTEvent::getNumberOfTriplets() {
1149   unsigned int triplets = 0;
1150   for (auto& it : n_triplets_by_layer_barrel_) {
1151     triplets += it;
1152   }
1153   for (auto& it : n_triplets_by_layer_endcap_) {
1154     triplets += it;
1155   }
1156 
1157   return triplets;
1158 }
1159 
1160 unsigned int LSTEvent::getNumberOfTripletsByLayerBarrel(unsigned int layer) {
1161   return n_triplets_by_layer_barrel_[layer];
1162 }
1163 
1164 unsigned int LSTEvent::getNumberOfTripletsByLayerEndcap(unsigned int layer) {
1165   return n_triplets_by_layer_endcap_[layer];
1166 }
1167 
1168 int LSTEvent::getNumberOfPixelTriplets() {
1169   auto nPixelTriplets_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1170 
1171   alpaka::memcpy(
1172       queue_, nPixelTriplets_buf_h, cms::alpakatools::make_device_view(queue_, (*pixelTripletsDC_)->nPixelTriplets()));
1173   alpaka::wait(queue_);
1174 
1175   return *nPixelTriplets_buf_h.data();
1176 }
1177 
1178 int LSTEvent::getNumberOfPixelQuintuplets() {
1179   auto nPixelQuintuplets_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1180 
1181   alpaka::memcpy(queue_,
1182                  nPixelQuintuplets_buf_h,
1183                  cms::alpakatools::make_device_view(queue_, (*pixelQuintupletsDC_)->nPixelQuintuplets()));
1184   alpaka::wait(queue_);
1185 
1186   return *nPixelQuintuplets_buf_h.data();
1187 }
1188 
1189 unsigned int LSTEvent::getNumberOfQuintuplets() {
1190   unsigned int quintuplets = 0;
1191   for (auto& it : n_quintuplets_by_layer_barrel_) {
1192     quintuplets += it;
1193   }
1194   for (auto& it : n_quintuplets_by_layer_endcap_) {
1195     quintuplets += it;
1196   }
1197 
1198   return quintuplets;
1199 }
1200 
1201 unsigned int LSTEvent::getNumberOfQuintupletsByLayerBarrel(unsigned int layer) {
1202   return n_quintuplets_by_layer_barrel_[layer];
1203 }
1204 
1205 unsigned int LSTEvent::getNumberOfQuintupletsByLayerEndcap(unsigned int layer) {
1206   return n_quintuplets_by_layer_endcap_[layer];
1207 }
1208 
1209 int LSTEvent::getNumberOfTrackCandidates() {
1210   auto nTrackCandidates_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1211 
1212   alpaka::memcpy(queue_,
1213                  nTrackCandidates_buf_h,
1214                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidates()));
1215   alpaka::wait(queue_);
1216 
1217   return *nTrackCandidates_buf_h.data();
1218 }
1219 
1220 int LSTEvent::getNumberOfPT5TrackCandidates() {
1221   auto nTrackCandidatesPT5_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1222 
1223   alpaka::memcpy(queue_,
1224                  nTrackCandidatesPT5_buf_h,
1225                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatespT5()));
1226   alpaka::wait(queue_);
1227 
1228   return *nTrackCandidatesPT5_buf_h.data();
1229 }
1230 
1231 int LSTEvent::getNumberOfPT3TrackCandidates() {
1232   auto nTrackCandidatesPT3_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1233 
1234   alpaka::memcpy(queue_,
1235                  nTrackCandidatesPT3_buf_h,
1236                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatespT3()));
1237   alpaka::wait(queue_);
1238 
1239   return *nTrackCandidatesPT3_buf_h.data();
1240 }
1241 
1242 int LSTEvent::getNumberOfPLSTrackCandidates() {
1243   auto nTrackCandidatesPLS_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1244 
1245   alpaka::memcpy(queue_,
1246                  nTrackCandidatesPLS_buf_h,
1247                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatespLS()));
1248   alpaka::wait(queue_);
1249 
1250   return *nTrackCandidatesPLS_buf_h.data();
1251 }
1252 
1253 int LSTEvent::getNumberOfPixelTrackCandidates() {
1254   auto nTrackCandidates_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1255   auto nTrackCandidatesT5_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1256 
1257   alpaka::memcpy(queue_,
1258                  nTrackCandidates_buf_h,
1259                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidates()));
1260   alpaka::memcpy(queue_,
1261                  nTrackCandidatesT5_buf_h,
1262                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatesT5()));
1263   alpaka::wait(queue_);
1264 
1265   return (*nTrackCandidates_buf_h.data()) - (*nTrackCandidatesT5_buf_h.data());
1266 }
1267 
1268 int LSTEvent::getNumberOfT5TrackCandidates() {
1269   auto nTrackCandidatesT5_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1270 
1271   alpaka::memcpy(queue_,
1272                  nTrackCandidatesT5_buf_h,
1273                  cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidatesT5()));
1274   alpaka::wait(queue_);
1275 
1276   return *nTrackCandidatesT5_buf_h.data();
1277 }
1278 
1279 template <typename TSoA, typename TDev>
1280 typename TSoA::ConstView LSTEvent::getInput(bool sync) {
1281   if constexpr (std::is_same_v<TDev, DevHost>) {
1282     return lstInputDC_->const_view<TSoA>();
1283   } else {
1284     // In case getTrimmedInput was called first
1285     if (!lstInputHC_ || lstInputHC_->sizes()[1] == 0) {
1286       lstInputHC_.emplace(
1287           cms::alpakatools::CopyToHost<PortableMultiCollection<TDev, HitsBaseSoA, PixelSeedsSoA>>::copyAsync(
1288               queue_, *lstInputDC_));
1289       if (sync)
1290         alpaka::wait(queue_);  // host consumers expect filled data
1291     }
1292     return lstInputHC_->const_view<TSoA>();
1293   }
1294 }
1295 template HitsBaseConst LSTEvent::getInput<HitsBaseSoA>(bool);
1296 template PixelSeedsConst LSTEvent::getInput<PixelSeedsSoA>(bool);
1297 
1298 template <typename TDev>
1299 HitsBaseConst LSTEvent::getTrimmedHitsBase(bool sync) {
1300   if constexpr (std::is_same_v<TDev, DevHost>) {
1301     return lstInputDC_->const_view<HitsBaseSoA>();
1302   } else {
1303     if (!lstInputHC_) {
1304       auto hits_d = lstInputDC_->view<HitsBaseSoA>();
1305       int nHits = hits_d.metadata().size();
1306       std::array<int, 2> const hits_sizes{{nHits, 0}};
1307       lstInputHC_.emplace(hits_sizes, queue_);
1308       auto hits_h = lstInputHC_->view<HitsBaseSoA>();
1309       auto idxs_h = cms::alpakatools::make_host_view(hits_h.idxs(), nHits);
1310       auto idxs_d = cms::alpakatools::make_device_view(queue_, hits_d.idxs(), nHits);
1311       alpaka::memcpy(queue_, idxs_h, idxs_d);
1312       if (sync)
1313         alpaka::wait(queue_);  // host consumers expect filled data
1314     }
1315     return lstInputHC_->const_view<HitsBaseSoA>();
1316   }
1317 }
1318 template HitsBaseConst LSTEvent::getTrimmedHitsBase(bool);
1319 
1320 template <typename TSoA, typename TDev>
1321 typename TSoA::ConstView LSTEvent::getHits(bool sync) {
1322   if constexpr (std::is_same_v<TDev, DevHost>) {
1323     return hitsDC_->const_view<TSoA>();
1324   } else {
1325     if (!hitsHC_) {
1326       hitsHC_.emplace(
1327           cms::alpakatools::CopyToHost<PortableMultiCollection<TDev, HitsExtendedSoA, HitsRangesSoA>>::copyAsync(
1328               queue_, *hitsDC_));
1329       if (sync)
1330         alpaka::wait(queue_);  // host consumers expect filled data
1331     }
1332     return hitsHC_->const_view<TSoA>();
1333   }
1334 }
1335 template HitsExtendedConst LSTEvent::getHits<HitsExtendedSoA>(bool);
1336 template HitsRangesConst LSTEvent::getHits<HitsRangesSoA>(bool);
1337 
1338 template <typename TDev>
1339 ObjectRangesConst LSTEvent::getRanges(bool sync) {
1340   if constexpr (std::is_same_v<TDev, DevHost>) {
1341     return rangesDC_->const_view();
1342   } else {
1343     if (!rangesHC_) {
1344       rangesHC_.emplace(
1345           cms::alpakatools::CopyToHost<PortableDeviceCollection<ObjectRangesSoA, TDev>>::copyAsync(queue_, *rangesDC_));
1346       if (sync)
1347         alpaka::wait(queue_);  // host consumers expect filled data
1348     }
1349     return rangesHC_->const_view();
1350   }
1351 }
1352 template ObjectRangesConst LSTEvent::getRanges<>(bool);
1353 
1354 template <typename TSoA, typename TDev>
1355 typename TSoA::ConstView LSTEvent::getMiniDoublets(bool sync) {
1356   if constexpr (std::is_same_v<TDev, DevHost>) {
1357     return miniDoubletsDC_->const_view<TSoA>();
1358   } else {
1359     if (!miniDoubletsHC_) {
1360       miniDoubletsHC_.emplace(
1361           cms::alpakatools::CopyToHost<
1362               PortableMultiCollection<TDev, MiniDoubletsSoA, MiniDoubletsOccupancySoA>>::copyAsync(queue_,
1363                                                                                                    *miniDoubletsDC_));
1364       if (sync)
1365         alpaka::wait(queue_);  // host consumers expect filled data
1366     }
1367     return miniDoubletsHC_->const_view<TSoA>();
1368   }
1369 }
1370 template MiniDoubletsConst LSTEvent::getMiniDoublets<MiniDoubletsSoA>(bool);
1371 template MiniDoubletsOccupancyConst LSTEvent::getMiniDoublets<MiniDoubletsOccupancySoA>(bool);
1372 
1373 template <typename TSoA, typename TDev>
1374 typename TSoA::ConstView LSTEvent::getSegments(bool sync) {
1375   if constexpr (std::is_same_v<TDev, DevHost>) {
1376     return segmentsDC_->const_view<TSoA>();
1377   } else {
1378     if (!segmentsHC_) {
1379       segmentsHC_.emplace(
1380           cms::alpakatools::CopyToHost<PortableMultiCollection<TDev, SegmentsSoA, SegmentsOccupancySoA>>::copyAsync(
1381               queue_, *segmentsDC_));
1382       if (sync)
1383         alpaka::wait(queue_);  // host consumers expect filled data
1384     }
1385     return segmentsHC_->const_view<TSoA>();
1386   }
1387 }
1388 template SegmentsConst LSTEvent::getSegments<SegmentsSoA>(bool);
1389 template SegmentsOccupancyConst LSTEvent::getSegments<SegmentsOccupancySoA>(bool);
1390 
1391 template <typename TDev>
1392 PixelSegmentsConst LSTEvent::getPixelSegments(bool sync) {
1393   if constexpr (std::is_same_v<TDev, DevHost>) {
1394     return pixelSegmentsDC_->const_view();
1395   } else {
1396     if (!pixelSegmentsHC_) {
1397       pixelSegmentsHC_.emplace(cms::alpakatools::CopyToHost<::PortableCollection<PixelSegmentsSoA, TDev>>::copyAsync(
1398           queue_, *pixelSegmentsDC_));
1399 
1400       if (sync)
1401         alpaka::wait(queue_);  // host consumers expect filled data
1402     }
1403   }
1404   return pixelSegmentsHC_->const_view();
1405 }
1406 template PixelSegmentsConst LSTEvent::getPixelSegments<>(bool);
1407 
1408 template <typename TSoA, typename TDev>
1409 typename TSoA::ConstView LSTEvent::getTriplets(bool sync) {
1410   if constexpr (std::is_same_v<TDev, DevHost>) {
1411     return tripletsDC_->const_view<TSoA>();
1412   } else {
1413     if (!tripletsHC_) {
1414       tripletsHC_.emplace(
1415           cms::alpakatools::CopyToHost<PortableMultiCollection<TDev, TripletsSoA, TripletsOccupancySoA>>::copyAsync(
1416               queue_, *tripletsDC_));
1417 
1418       if (sync)
1419         alpaka::wait(queue_);  // host consumers expect filled data
1420     }
1421   }
1422   return tripletsHC_->const_view<TSoA>();
1423 }
1424 template TripletsConst LSTEvent::getTriplets<TripletsSoA>(bool);
1425 template TripletsOccupancyConst LSTEvent::getTriplets<TripletsOccupancySoA>(bool);
1426 
1427 template <typename TSoA, typename TDev>
1428 typename TSoA::ConstView LSTEvent::getQuintuplets(bool sync) {
1429   if constexpr (std::is_same_v<TDev, DevHost>) {
1430     return quintupletsDC_->const_view<TSoA>();
1431   } else {
1432     if (!quintupletsHC_) {
1433       quintupletsHC_.emplace(
1434           cms::alpakatools::CopyToHost<PortableMultiCollection<TDev, QuintupletsSoA, QuintupletsOccupancySoA>>::copyAsync(
1435               queue_, *quintupletsDC_));
1436 
1437       if (sync)
1438         alpaka::wait(queue_);  // host consumers expect filled data
1439     }
1440   }
1441   return quintupletsHC_->const_view<TSoA>();
1442 }
1443 template QuintupletsConst LSTEvent::getQuintuplets<QuintupletsSoA>(bool);
1444 template QuintupletsOccupancyConst LSTEvent::getQuintuplets<QuintupletsOccupancySoA>(bool);
1445 
1446 template <typename TDev>
1447 PixelTripletsConst LSTEvent::getPixelTriplets(bool sync) {
1448   if constexpr (std::is_same_v<TDev, DevHost>) {
1449     return pixelTripletsDC_->const_view();
1450   } else {
1451     if (!pixelTripletsHC_) {
1452       pixelTripletsHC_.emplace(cms::alpakatools::CopyToHost<::PortableCollection<PixelTripletsSoA, TDev>>::copyAsync(
1453           queue_, *pixelTripletsDC_));
1454 
1455       if (sync)
1456         alpaka::wait(queue_);  // host consumers expect filled data
1457     }
1458   }
1459   return pixelTripletsHC_->const_view();
1460 }
1461 template PixelTripletsConst LSTEvent::getPixelTriplets<>(bool);
1462 
1463 template <typename TDev>
1464 PixelQuintupletsConst LSTEvent::getPixelQuintuplets(bool sync) {
1465   if constexpr (std::is_same_v<TDev, DevHost>) {
1466     return pixelQuintupletsDC_->const_view();
1467   } else {
1468     if (!pixelQuintupletsHC_) {
1469       pixelQuintupletsHC_.emplace(
1470           cms::alpakatools::CopyToHost<::PortableCollection<PixelQuintupletsSoA, TDev>>::copyAsync(
1471               queue_, *pixelQuintupletsDC_));
1472 
1473       if (sync)
1474         alpaka::wait(queue_);  // host consumers expect filled data
1475     }
1476   }
1477   return pixelQuintupletsHC_->const_view();
1478 }
1479 template PixelQuintupletsConst LSTEvent::getPixelQuintuplets<>(bool);
1480 
1481 const TrackCandidatesConst& LSTEvent::getTrackCandidates(bool inCMSSW, bool sync) {
1482   if (!trackCandidatesHC_) {
1483     // Get nTrackCanHost parameter to initialize host based instance
1484     auto nTrackCanHost_buf_h = cms::alpakatools::make_host_buffer<unsigned int>(queue_);
1485     alpaka::memcpy(queue_,
1486                    nTrackCanHost_buf_h,
1487                    cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->nTrackCandidates()));
1488     alpaka::wait(queue_);  // wait here before we get nTrackCanHost and trackCandidatesInCPU becomes usable
1489 
1490     auto const nTrackCanHost = *nTrackCanHost_buf_h.data();
1491     trackCandidatesHC_.emplace(nTrackCanHost, queue_);
1492 
1493     (*trackCandidatesHC_)->nTrackCandidates() = nTrackCanHost;
1494     alpaka::memcpy(queue_,
1495                    cms::alpakatools::make_host_view((*trackCandidatesHC_)->hitIndices()->data(),
1496                                                     Params_pT5::kHits * nTrackCanHost),
1497                    cms::alpakatools::make_device_view(
1498                        queue_, (*trackCandidatesDC_)->hitIndices()->data(), Params_pT5::kHits * nTrackCanHost));
1499     alpaka::memcpy(queue_,
1500                    cms::alpakatools::make_host_view((*trackCandidatesHC_)->pixelSeedIndex(), nTrackCanHost),
1501                    cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->pixelSeedIndex(), nTrackCanHost));
1502     if (not inCMSSW) {
1503       alpaka::memcpy(queue_,
1504                      cms::alpakatools::make_host_view((*trackCandidatesHC_)->logicalLayers()->data(),
1505                                                       Params_pT5::kLayers * nTrackCanHost),
1506                      cms::alpakatools::make_device_view(
1507                          queue_, (*trackCandidatesDC_)->logicalLayers()->data(), Params_pT5::kLayers * nTrackCanHost));
1508       alpaka::memcpy(
1509           queue_,
1510           cms::alpakatools::make_host_view((*trackCandidatesHC_)->directObjectIndices(), nTrackCanHost),
1511           cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->directObjectIndices(), nTrackCanHost));
1512       alpaka::memcpy(
1513           queue_,
1514           cms::alpakatools::make_host_view((*trackCandidatesHC_)->objectIndices()->data(), 2 * nTrackCanHost),
1515           cms::alpakatools::make_device_view(
1516               queue_, (*trackCandidatesDC_)->objectIndices()->data(), 2 * nTrackCanHost));
1517     }
1518     alpaka::memcpy(
1519         queue_,
1520         cms::alpakatools::make_host_view((*trackCandidatesHC_)->trackCandidateType(), nTrackCanHost),
1521         cms::alpakatools::make_device_view(queue_, (*trackCandidatesDC_)->trackCandidateType(), nTrackCanHost));
1522     if (sync)
1523       alpaka::wait(queue_);  // host consumers expect filled data
1524   }
1525   return trackCandidatesHC_.value().const_view();
1526 }
1527 
1528 template <typename TSoA, typename TDev>
1529 typename TSoA::ConstView LSTEvent::getModules(bool sync) {
1530   if constexpr (std::is_same_v<TDev, DevHost>) {
1531     return modules_.const_view<TSoA>();
1532   } else {
1533     if (!modulesHC_) {
1534       modulesHC_.emplace(
1535           cms::alpakatools::CopyToHost<PortableMultiCollection<TDev, ModulesSoA, ModulesPixelSoA>>::copyAsync(
1536               queue_, modules_));
1537       if (sync)
1538         alpaka::wait(queue_);  // host consumers expect filled data
1539     }
1540     return modulesHC_->const_view<TSoA>();
1541   }
1542 }
1543 template ModulesConst LSTEvent::getModules<ModulesSoA>(bool);
1544 template ModulesPixelConst LSTEvent::getModules<ModulesPixelSoA>(bool);