Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:12:16

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 // CUDA plugin include(s).
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 // System include(s).
0017 #include <cassert>
0018 #include <cmath>
0019 
0020 namespace {
0021 
0022 /// Type of "other spacepoint" passed to the kernel
0023 enum OtherSPType : int {
0024   BottomSP = 0,  //< The "other" spacepoint is a bottom one
0025   TopSP = 1      //< The "other" spacepoint is a top one
0026 };
0027 
0028 }  // namespace
0029 
0030 namespace Acts {
0031 namespace Cuda {
0032 namespace Kernels {
0033 
0034 /// @name Helper functions for calculating "deltaR" between spacepoints
0035 /// @{
0036 
0037 /// Prototype for the @c Acts::Cuda::Kernels::getDeltaR functions
0038 template <int SPType>
0039 __device__ float getDeltaR(float /*middleR*/, float /*otherR*/) {
0040   // This function should *never* be called.
0041   assert(false);
0042   return 0.0f;
0043 }
0044 
0045 /// Function specialisation for middle-bottom spacepoint pairs
0046 template <>
0047 __device__ float getDeltaR<BottomSP>(float middleR, float bottomR) {
0048   return middleR - bottomR;
0049 }
0050 
0051 /// Function specialisation for middle-top spacepoint pairs
0052 template <>
0053 __device__ float getDeltaR<TopSP>(float middleR, float topR) {
0054   return topR - middleR;
0055 }
0056 
0057 /// @}
0058 
0059 /// @name Helper functions for calculating "cotTheta" between spacepoints
0060 /// @{
0061 
0062 /// Prototype for the @c Acts::Cuda::Kernels::getCotTheta functions
0063 template <int SPType>
0064 __device__ float getCotTheta(float /*middleZ*/, float /*otherZ*/,
0065                              float /*deltaR*/) {
0066   // This function should *never* be called.
0067   assert(false);
0068   return 0.0f;
0069 }
0070 
0071 /// Function specialisation for middle-bottom spacepoint pairs
0072 template <>
0073 __device__ float getCotTheta<BottomSP>(float middleZ, float bottomZ,
0074                                        float deltaR) {
0075   return (middleZ - bottomZ) / deltaR;
0076 }
0077 
0078 /// Function specialisation for middle-top spacepoint pairs
0079 template <>
0080 __device__ float getCotTheta<TopSP>(float middleZ, float topZ, float deltaR) {
0081   return (topZ - middleZ) / deltaR;
0082 }
0083 
0084 /// @}
0085 
0086 /// Kernel performing the spacepoint dublet finding
0087 ///
0088 /// @tparam SPType The "other" spacepoint type used in the call
0089 ///         (@c ::BottomSP or @c TopSP)
0090 /// @param[in] nMiddleSPs The number of middle spacepoints in @c middleSPs
0091 /// @param[in] middleSPs Properties of all of the middle spacepoints
0092 /// @param[in] nOtherSPs The number of other (bottom or top) spacepoints in
0093 ///            @c otherSPs
0094 /// @param[in] otherSPs Properties of all of the other (bottom or top)
0095 ///            spacepoints
0096 /// @param[in] deltaRMin Configuration parameter from @c Acts::SeedFinderConfig
0097 /// @param[in] deltaRMax Configuration parameter from @c Acts::SeedFinderConfig
0098 /// @param[in] cotThetaMax Configuration parameter from
0099 ///            @c Acts::SeedFinderConfig
0100 /// @param[in] collisionRegionMin Configuration parameter from
0101 ///            @c Acts::SeedFinderConfig
0102 /// @param[in] collisionRegionMax Configuration parameter from
0103 ///            @c Acts::SeedFinderConfig
0104 /// @param[out] dubletCounts 1-D array of the middle-other dublets found
0105 ///             for each middle spacepoint
0106 /// @param[out] dublets 2-D matrix of size @c nMiddleSPs x @c nOtherSPs, holding
0107 ///             the "other" spacepoint indices for the identified dublets
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   // Figure out which dublet the kernel operates on.
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   // If we're outside of bounds, stop here.
0122   if ((middleIndex >= nMiddleSPs) || (otherIndex >= nOtherSPs)) {
0123     return;
0124   }
0125 
0126   // Calculate variables used in the compatibility check.
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   // Perform the compatibility check.
0135   const bool isCompatible =
0136       ((deltaR >= deltaRMin) && (deltaR <= deltaRMax) &&
0137        (fabs(cotTheta) <= cotThetaMax) && (zOrigin >= collisionRegionMin) &&
0138        (zOrigin <= collisionRegionMax));
0139 
0140   // If they are compatible, save their indices into the output matrix.
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 }  // namespace Kernels
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   // Calculate the parallelisation for the middle<->bottom spacepoint
0165   // compatibility flagging.
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   // Launch the middle-bottom dublet finding.
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   // Calculate the parallelisation for the middle<->top spacepoint
0178   // compatibility flagging.
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   // Launch the middle-bottom dublet finding.
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 }  // namespace Details
0194 }  // namespace Cuda
0195 }  // namespace Acts