gpu_nbody.cu 3.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2019 Mael Keryell
  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. #include <stdio.h>
  17. #include <stdint.h>
  18. #include <math.h>
  19. #include <starpu.h>
  20. struct Params
  21. {
  22. unsigned taskx;
  23. unsigned epsilon;
  24. };
  25. __global__ void gpuNbodyKernel(double *P, double *subA, double *M,
  26. uint32_t nxP, uint32_t nxA, uint32_t nxM,
  27. uint32_t ldP, uint32_t ldA,
  28. struct Params params)
  29. {
  30. uint32_t id, i, j, k;
  31. double dx, dy, modul;
  32. id = blockIdx.x * blockDim.x + threadIdx.x;
  33. i = id % nxA;
  34. j = id / nxA;
  35. if (j >= 1){
  36. return;
  37. }
  38. double sumaccx;
  39. double sumaccy;
  40. for (k = 0; k < nxP; k++){
  41. if (k != id + nxA*params.taskx){
  42. dx = P[k] - P[id + nxA*params.taskx];
  43. dy = P[k + ldP] - P[id + nxA*params.taskx + ldP];
  44. modul = dx * dx + dy * dy;
  45. sumaccx = 6.674e-11 * M[k] * dx / pow(modul + params.epsilon, 3);
  46. sumaccy = 6.674e-11 * M[k] * dy / pow(modul + params.epsilon, 3);
  47. }
  48. }
  49. subA[i] = sumaccx;
  50. subA[i + ldA] = sumaccy;
  51. // P[id + nxA * params.taskx] = subA[i];
  52. // subA[i] = 0;
  53. // subA[i + ldA] = 1;
  54. }
  55. #define THREADS_PER_BLOCK 64
  56. extern "C" void gpu_nbody(void * descr[], void * args)
  57. {
  58. double *d_P, *d_subA, *d_M;
  59. uint32_t nxP, nxA, nxM;
  60. uint32_t ldA, ldP;
  61. uint32_t nblocks;
  62. struct Params *params = (struct Params *) args;
  63. d_P = (double *) STARPU_MATRIX_GET_PTR(descr[0]);
  64. d_subA = (double *) STARPU_MATRIX_GET_PTR(descr[1]);
  65. d_M = (double *) STARPU_MATRIX_GET_PTR(descr[2]);
  66. nxP = STARPU_MATRIX_GET_NX(descr[0]);
  67. nxA = STARPU_MATRIX_GET_NX(descr[1]);
  68. nxM = STARPU_MATRIX_GET_NX(descr[2]);
  69. ldP = STARPU_MATRIX_GET_LD(descr[0]);
  70. ldA = STARPU_MATRIX_GET_LD(descr[1]);
  71. nblocks = (nxA + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
  72. gpuNbodyKernel
  73. <<< nblocks, THREADS_PER_BLOCK, 0, starpu_cuda_get_local_stream()
  74. >>> (d_P, d_subA, d_M, nxP, nxA, nxM, ldP, ldA, *params);
  75. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  76. }
  77. __global__ void gpuNbody2Kernel(double *d_subP, double *d_subV, double *d_subA,
  78. uint32_t nxP, uint32_t nxV, uint32_t nxA,
  79. uint32_t ldP, uint32_t ldV, uint32_t ldA,
  80. struct Params params)
  81. {
  82. uint32_t id, i, j;
  83. id = blockIdx.x * blockDim.x + threadIdx.x;
  84. i = id % nxP;
  85. j = id / nxP;
  86. if (j >= 1){
  87. return;
  88. }
  89. d_subV[i] = d_subV[i] + 3600*d_subA[i];
  90. d_subV[i + ldV] = d_subV[i + ldV] + 3600*d_subA[i + ldA];
  91. d_subP[i] = d_subP[i] + 3600*d_subV[i];
  92. d_subP[i + ldP] = d_subP[i + ldP] + 3600*d_subV[i + ldV];
  93. }
  94. extern "C" void gpu_nbody2(void * descr[], void *args)
  95. {
  96. double *d_subP, *d_subV, *d_subA;
  97. uint32_t nxP, nxV, nxA;
  98. uint32_t ldP, ldV, ldA;
  99. uint32_t nblocks;
  100. struct Params *params = (struct Params *) args;
  101. d_subP = (double *) STARPU_MATRIX_GET_PTR(descr[0]);
  102. d_subV = (double *) STARPU_MATRIX_GET_PTR(descr[1]);
  103. d_subA = (double *) STARPU_MATRIX_GET_PTR(descr[2]);
  104. nxP = STARPU_MATRIX_GET_NX(descr[0]);
  105. nxV = STARPU_MATRIX_GET_NX(descr[1]);
  106. nxA = STARPU_MATRIX_GET_NX(descr[2]);
  107. ldP = STARPU_MATRIX_GET_LD(descr[0]);
  108. ldV = STARPU_MATRIX_GET_LD(descr[1]);
  109. ldA = STARPU_MATRIX_GET_LD(descr[2]);
  110. nblocks = (nxA + THREADS_PER_BLOCK - 1) / THREADS_PER_BLOCK;
  111. gpuNbody2Kernel
  112. <<< nblocks, THREADS_PER_BLOCK, 0, starpu_cuda_get_local_stream()
  113. >>> (d_subP, d_subV, d_subA, nxP, nxV, nxA, ldP, ldV, ldA, *params);
  114. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  115. }