axpy_partition_gpu.h 4.3 KB

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