File indexing completed on 2025-01-18 09:12:16
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "Acts/Plugins/Cuda/Seeding2/Details/FindDublets.hpp"
0011 #include "Acts/Plugins/Cuda/Seeding2/Details/Types.hpp"
0012
0013 #include "../Utilities/ErrorCheck.cuh"
0014 #include "../Utilities/MatrixMacros.hpp"
0015
0016
0017 #include <cassert>
0018 #include <cmath>
0019
0020 namespace {
0021
0022
0023 enum OtherSPType : int {
0024 BottomSP = 0,
0025 TopSP = 1
0026 };
0027
0028 }
0029
0030 namespace Acts {
0031 namespace Cuda {
0032 namespace Kernels {
0033
0034
0035
0036
0037
0038 template <int SPType>
0039 __device__ float getDeltaR(float , float ) {
0040
0041 assert(false);
0042 return 0.0f;
0043 }
0044
0045
0046 template <>
0047 __device__ float getDeltaR<BottomSP>(float middleR, float bottomR) {
0048 return middleR - bottomR;
0049 }
0050
0051
0052 template <>
0053 __device__ float getDeltaR<TopSP>(float middleR, float topR) {
0054 return topR - middleR;
0055 }
0056
0057
0058
0059
0060
0061
0062
0063 template <int SPType>
0064 __device__ float getCotTheta(float , float ,
0065 float ) {
0066
0067 assert(false);
0068 return 0.0f;
0069 }
0070
0071
0072 template <>
0073 __device__ float getCotTheta<BottomSP>(float middleZ, float bottomZ,
0074 float deltaR) {
0075 return (middleZ - bottomZ) / deltaR;
0076 }
0077
0078
0079 template <>
0080 __device__ float getCotTheta<TopSP>(float middleZ, float topZ, float deltaR) {
0081 return (topZ - middleZ) / deltaR;
0082 }
0083
0084
0085
0086
0087
0088
0089
0090
0091
0092
0093
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105
0106
0107
0108
0109 template <int SPType>
0110 __global__ void findDublets(std::size_t nMiddleSPs,
0111 const Details::SpacePoint* middleSPs,
0112 std::size_t nOtherSPs,
0113 const Details::SpacePoint* otherSPs,
0114 float deltaRMin, float deltaRMax, float cotThetaMax,
0115 float collisionRegionMin, float collisionRegionMax,
0116 unsigned int* dubletCounts, std::size_t* dublets) {
0117
0118 const std::size_t middleIndex = blockIdx.x * blockDim.x + threadIdx.x;
0119 const std::size_t otherIndex = blockIdx.y * blockDim.y + threadIdx.y;
0120
0121
0122 if ((middleIndex >= nMiddleSPs) || (otherIndex >= nOtherSPs)) {
0123 return;
0124 }
0125
0126
0127 const float deltaR = getDeltaR<SPType>(middleSPs[middleIndex].radius,
0128 otherSPs[otherIndex].radius);
0129 const float cotTheta = getCotTheta<SPType>(middleSPs[middleIndex].z,
0130 otherSPs[otherIndex].z, deltaR);
0131 const float zOrigin =
0132 middleSPs[middleIndex].z - middleSPs[middleIndex].radius * cotTheta;
0133
0134
0135 const bool isCompatible =
0136 ((deltaR >= deltaRMin) && (deltaR <= deltaRMax) &&
0137 (fabs(cotTheta) <= cotThetaMax) && (zOrigin >= collisionRegionMin) &&
0138 (zOrigin <= collisionRegionMax));
0139
0140
0141 if (isCompatible) {
0142 const unsigned int dubletRow = atomicAdd(dubletCounts + middleIndex, 1);
0143 ACTS_CUDA_MATRIX2D_ELEMENT(dublets, nMiddleSPs, nOtherSPs, middleIndex,
0144 dubletRow) = otherIndex;
0145 }
0146 return;
0147 }
0148
0149 }
0150
0151 namespace Details {
0152
0153 void findDublets(std::size_t maxBlockSize, std::size_t nBottomSPs,
0154 const device_array<SpacePoint>& bottomSPs,
0155 std::size_t nMiddleSPs,
0156 const device_array<SpacePoint>& middleSPs, std::size_t nTopSPs,
0157 const device_array<SpacePoint>& topSPs, float deltaRMin,
0158 float deltaRMax, float cotThetaMax, float collisionRegionMin,
0159 float collisionRegionMax,
0160 device_array<unsigned int>& middleBottomCounts,
0161 device_array<std::size_t>& middleBottomDublets,
0162 device_array<unsigned int>& middleTopCounts,
0163 device_array<std::size_t>& middleTopDublets) {
0164
0165
0166 const dim3 blockSizeMB(1, maxBlockSize);
0167 const dim3 numBlocksMB((nMiddleSPs + blockSizeMB.x - 1) / blockSizeMB.x,
0168 (nBottomSPs + blockSizeMB.y - 1) / blockSizeMB.y);
0169
0170
0171 Kernels::findDublets<BottomSP><<<numBlocksMB, blockSizeMB>>>(
0172 nMiddleSPs, middleSPs.get(), nBottomSPs, bottomSPs.get(), deltaRMin,
0173 deltaRMax, cotThetaMax, collisionRegionMin, collisionRegionMax,
0174 middleBottomCounts.get(), middleBottomDublets.get());
0175 ACTS_CUDA_ERROR_CHECK(cudaGetLastError());
0176
0177
0178
0179 const dim3 blockSizeMT(1, maxBlockSize);
0180 const dim3 numBlocksMT((nMiddleSPs + blockSizeMT.x - 1) / blockSizeMT.x,
0181 (nTopSPs + blockSizeMT.y - 1) / blockSizeMT.y);
0182
0183
0184 Kernels::findDublets<TopSP><<<numBlocksMT, blockSizeMT>>>(
0185 nMiddleSPs, middleSPs.get(), nTopSPs, topSPs.get(), deltaRMin, deltaRMax,
0186 cotThetaMax, collisionRegionMin, collisionRegionMax,
0187 middleTopCounts.get(), middleTopDublets.get());
0188 ACTS_CUDA_ERROR_CHECK(cudaGetLastError());
0189 ACTS_CUDA_ERROR_CHECK(cudaDeviceSynchronize());
0190 return;
0191 }
0192
0193 }
0194 }
0195 }