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