Back to home page

Project CMSSW displayed by LXR

 
 

    


Warning, /RecoLocalTracker/SiPixelRecHits/plugins/PixelRecHitGPUKernel.cu is written in an unsupported language. File is not indexed.

0001 // C++ headers
0002 #include <algorithm>
0003 #include <numeric>
0004 
0005 // CUDA runtime
0006 #include <cuda_runtime.h>
0007 
0008 // CMSSW headers
0009 #include "CUDADataFormats/SiPixelCluster/interface/gpuClusteringConstants.h"
0010 #include "HeterogeneousCore/CUDAUtilities/interface/cudaCheck.h"
0011 #include "HeterogeneousCore/CUDAUtilities/interface/device_unique_ptr.h"
0012 
0013 #include "PixelRecHitGPUKernel.h"
0014 #include "gpuPixelRecHits.h"
0015 
0016 //#define GPU_DEBUG
0017 
0018 namespace {
0019   template <typename TrackerTraits>
0020   __global__ void setHitsLayerStart(uint32_t const* __restrict__ hitsModuleStart,
0021                                     pixelCPEforGPU::ParamsOnGPUT<TrackerTraits> const* cpeParams,
0022                                     uint32_t* hitsLayerStart) {
0023     auto i = blockIdx.x * blockDim.x + threadIdx.x;
0024     constexpr auto m = TrackerTraits::numberOfLayers;
0025 
0026     assert(0 == hitsModuleStart[0]);
0027 
0028     if (i <= m) {
0029       hitsLayerStart[i] = hitsModuleStart[cpeParams->layerGeometry().layerStart[i]];
0030 #ifdef GPU_DEBUG
0031       int old = i == 0 ? 0 : hitsModuleStart[cpeParams->layerGeometry().layerStart[i - 1]];
0032       printf("LayerStart %d/%d at module %d: %d - %d\n",
0033              i,
0034              m,
0035              cpeParams->layerGeometry().layerStart[i],
0036              hitsLayerStart[i],
0037              hitsLayerStart[i] - old);
0038 #endif
0039     }
0040   }
0041 }  // namespace
0042 
0043 namespace pixelgpudetails {
0044 
0045   template <typename TrackerTraits>
0046   TrackingRecHitSoADevice<TrackerTraits> PixelRecHitGPUKernel<TrackerTraits>::makeHitsAsync(
0047       SiPixelDigisCUDA const& digis_d,
0048       SiPixelClustersCUDA const& clusters_d,
0049       BeamSpotCUDA const& bs_d,
0050       pixelCPEforGPU::ParamsOnGPUT<TrackerTraits> const* cpeParams,
0051       cudaStream_t stream) const {
0052     using namespace gpuPixelRecHits;
0053     auto nHits = clusters_d.nClusters();
0054 
0055     TrackingRecHitSoADevice<TrackerTraits> hits_d(
0056         nHits, clusters_d.offsetBPIX2(), cpeParams, clusters_d->clusModuleStart(), stream);
0057 
0058     int activeModulesWithDigis = digis_d.nModules();
0059     // protect from empty events
0060     if (activeModulesWithDigis) {
0061       int threadsPerBlock = 128;
0062       int blocks = activeModulesWithDigis;
0063 
0064 #ifdef GPU_DEBUG
0065       std::cout << "launching getHits kernel for " << blocks << " blocks" << std::endl;
0066 #endif
0067       getHits<TrackerTraits><<<blocks, threadsPerBlock, 0, stream>>>(
0068           cpeParams, bs_d.data(), digis_d.view(), digis_d.nDigis(), clusters_d.const_view(), hits_d.view());
0069       cudaCheck(cudaGetLastError());
0070 #ifdef GPU_DEBUG
0071       cudaCheck(cudaDeviceSynchronize());
0072 #endif
0073 
0074       // assuming full warp of threads is better than a smaller number...
0075       if (nHits) {
0076         setHitsLayerStart<TrackerTraits>
0077             <<<1, 32, 0, stream>>>(clusters_d->clusModuleStart(), cpeParams, hits_d.view().hitsLayerStart().data());
0078         cudaCheck(cudaGetLastError());
0079         constexpr auto nLayers = TrackerTraits::numberOfLayers;
0080         cms::cuda::fillManyFromVector(&(hits_d.view().phiBinner()),
0081                                       nLayers,
0082                                       hits_d.view().iphi(),
0083                                       hits_d.view().hitsLayerStart().data(),
0084                                       nHits,
0085                                       256,
0086                                       hits_d.view().phiBinnerStorage(),
0087                                       stream);
0088         cudaCheck(cudaGetLastError());
0089 
0090 #ifdef GPU_DEBUG
0091         cudaCheck(cudaDeviceSynchronize());
0092 #endif
0093       }
0094     }
0095 
0096 #ifdef GPU_DEBUG
0097     cudaCheck(cudaDeviceSynchronize());
0098     std::cout << "PixelRecHitGPUKernel -> DONE!" << std::endl;
0099 #endif
0100 
0101     return hits_d;
0102   }
0103 
0104   template class PixelRecHitGPUKernel<pixelTopology::Phase1>;
0105   template class PixelRecHitGPUKernel<pixelTopology::Phase2>;
0106   template class PixelRecHitGPUKernel<pixelTopology::HIonPhase1>;
0107 }  // namespace pixelgpudetails