axpy_partition_gpu.h 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2016-2020 Université de Bordeaux, CNRS (LaBRI UMR 5800), Inria
  4. * Copyright (C) 2016 Uppsala University
  5. *
  6. * StarPU is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU Lesser General Public License as published by
  8. * the Free Software Foundation; either version 2.1 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * StarPU is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  14. *
  15. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  16. */
  17. /*
  18. * This creates two dumb vectors, splits them into chunks, and for each pair of
  19. * chunk, run axpy on them.
  20. */
  21. #pragma once
  22. __device__ static uint get_smid(void)
  23. {
  24. #if defined(__CUDACC__)
  25. uint ret;
  26. asm("mov.u32 %0, %smid;" : "=r"(ret) );
  27. return ret;
  28. #else
  29. return 0;
  30. #endif
  31. }
  32. #define __P_HKARGS dimGrid, active_blocks ,occupancy, block_assignment_d, mapping_start
  33. #define __P_KARGS dim3 blocks, int active_blocks, int occupancy, unsigned int* block_assignment, int mapping_start
  34. #define __P_DARGS blocks,blockid
  35. #define __P_BEGIN \
  36. __shared__ unsigned int block_start; \
  37. int smid = get_smid(); \
  38. if(threadIdx.x == 0 && threadIdx.y == 0 && threadIdx.z == 0) \
  39. { \
  40. block_start = atomicDec(&block_assignment[smid],0xDEADBEEF); \
  41. } \
  42. __syncthreads(); \
  43. \
  44. if(block_start > active_blocks) \
  45. { \
  46. return; \
  47. }
  48. #define __P_LOOPXY \
  49. dim3 blockid; \
  50. blockid.z = 0; \
  51. \
  52. int gridDim_sum = blocks.x*blocks.y; \
  53. int startBlock = block_start + (smid - mapping_start) * occupancy; \
  54. \
  55. for(int blockid_sum = startBlock; blockid_sum < gridDim_sum; blockid_sum +=active_blocks) \
  56. { \
  57. blockid.x = blockid_sum % blocks.x; \
  58. blockid.y = blockid_sum / blocks.x;
  59. #define __P_LOOPEND }
  60. // Needed if shared memory is used
  61. #define __P_LOOPEND_SAFE __syncthreads(); }
  62. #define __P_LOOPX \
  63. dim3 blockid; \
  64. blockid.z = 0; \
  65. blockid.y = 0; \
  66. int gridDim_sum = blocks.x; \
  67. int startBlock = (smid-mapping_start) + block_start*(active_blocks/occupancy); \
  68. \
  69. for(int blockid_sum = startBlock; blockid_sum < gridDim_sum; blockid_sum +=active_blocks) \
  70. { \
  71. blockid.x = blockid_sum;
  72. // int startBlock = block_start + (smid - mapping_start) * occupancy; \
  73. //////////// HOST side functions
  74. template <typename F>
  75. static void buildPartitionedBlockMapping(F cudaFun, int threads, int shmem, int mapping_start, int allocation,
  76. int &width, int &active_blocks, unsigned int *block_assignment_d,cudaStream_t current_stream =
  77. #ifdef cudaStreamPerThread
  78. cudaStreamPerThread
  79. #else
  80. NULL
  81. #endif
  82. )
  83. {
  84. int occupancy;
  85. int nb_SM = 13; //TODO: replace with call
  86. int mapping_end = mapping_start + allocation - 1; // exclusive
  87. unsigned int block_assignment[15];
  88. #if CUDART_VERSION >= 6050
  89. cudaOccupancyMaxActiveBlocksPerMultiprocessor(&occupancy,cudaFun,threads,shmem);
  90. #else
  91. occupancy = 4;
  92. #endif
  93. width = occupancy * nb_SM; // Physical wrapper grid size. Fits GPU exactly
  94. active_blocks = occupancy*allocation; // The total number of blocks doing work
  95. for(int i = 0; i < mapping_start; i++)
  96. block_assignment[i] = (unsigned) -1;
  97. for(int i = mapping_start; i <= mapping_end; i++)
  98. {
  99. block_assignment[i] = occupancy - 1;
  100. }
  101. for(int i = mapping_end+1; i < nb_SM; i++)
  102. block_assignment[i] = (unsigned) -1;
  103. cudaMemcpyAsync((void*)block_assignment_d,block_assignment,sizeof(block_assignment),cudaMemcpyHostToDevice, current_stream);
  104. //cudaMemcpy((void*)block_assignment_d,block_assignment,sizeof(block_assignment),cudaMemcpyHostToDevice);
  105. //cudaDeviceSynchronize();
  106. }
  107. #define __P_HOSTSETUP(KERNEL,GRIDDIM,BLOCKSIZE,SHMEMSIZE,MAPPING_START,MAPPING_END,STREAM) \
  108. unsigned int* block_assignment_d; cudaMalloc((void**) &block_assignment_d,15*sizeof(unsigned int)); \
  109. int width = 0; \
  110. int active_blocks = 0; \
  111. buildPartitionedBlockMapping(KERNEL,BLOCKSIZE,SHMEMSIZE,(MAPPING_START),(MAPPING_END)-(MAPPING_START), \
  112. width, active_blocks, block_assignment_d,STREAM); \
  113. int occupancy = active_blocks/((MAPPING_END)-(MAPPING_START)); \
  114. dim3 dimGrid = (GRIDDIM);\
  115. int mapping_start = (MAPPING_START);