pi_redux.c 7.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284
  1. /*
  2. * StarPU
  3. * Copyright (C) Université Bordeaux 1, CNRS 2008-2010 (see AUTHORS file)
  4. *
  5. * This program 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. * This program 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 <starpu.h>
  17. #include <starpu_cuda.h>
  18. #include <stdlib.h>
  19. #define PI 3.14159265358979323846
  20. #if defined(STARPU_USE_CUDA) && defined(STARPU_HAVE_CURAND)
  21. #error CURAND is required to run that example on CUDA devices
  22. #endif
  23. #ifdef STARPU_USE_CUDA
  24. #include <cuda.h>
  25. #include <curand.h>
  26. #endif
  27. #define NSHOT_PER_TASK (1024*1024)
  28. /* default value */
  29. static unsigned ntasks = 1024;
  30. /*
  31. * Initialization of the Random Number Generators (RNG)
  32. */
  33. #ifdef STARPU_USE_CUDA
  34. /* RNG for the CURAND library */
  35. static curandGenerator_t curandgens[STARPU_NMAXWORKERS];
  36. #endif
  37. /* state for the erand48 function */
  38. static unsigned short xsubi[3*STARPU_NMAXWORKERS];
  39. /* Function to initialize the random number generator in the current worker */
  40. static void init_rng(void *arg __attribute__((unused)))
  41. {
  42. #ifdef STARPU_USE_CUDA
  43. curandStatus_t res;
  44. #endif
  45. int workerid = starpu_worker_get_id();
  46. switch (starpu_worker_get_type(workerid)) {
  47. case STARPU_CPU_WORKER:
  48. /* create a seed */
  49. xsubi[0 + 3*workerid] = (unsigned short)workerid;
  50. xsubi[1 + 3*workerid] = (unsigned short)workerid;
  51. xsubi[2 + 3*workerid] = (unsigned short)workerid;
  52. break;
  53. #ifdef STARPU_USE_CUDA
  54. case STARPU_CUDA_WORKER:
  55. /* Create a RNG */
  56. res = curandCreateGenerator(&curandgens[workerid],
  57. CURAND_RNG_PSEUDO_DEFAULT);
  58. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  59. /* Seed it with worker's id */
  60. res = curandSetPseudoRandomGeneratorSeed(curandgens[workerid],
  61. (unsigned long long)workerid);
  62. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  63. break;
  64. #endif
  65. default:
  66. STARPU_ABORT();
  67. break;
  68. }
  69. }
  70. static void parse_args(int argc, char **argv)
  71. {
  72. int i;
  73. for (i = 1; i < argc; i++) {
  74. if (strcmp(argv[i], "-ntasks") == 0) {
  75. char *argptr;
  76. ntasks = strtol(argv[++i], &argptr, 10);
  77. }
  78. }
  79. }
  80. /*
  81. * Monte-carlo kernel
  82. */
  83. static void pi_func_cpu(void *descr[], void *cl_arg __attribute__ ((unused)))
  84. {
  85. int workerid = starpu_worker_get_id();
  86. unsigned short *worker_xsub;
  87. worker_xsub = &xsubi[3*workerid];
  88. unsigned long local_cnt = 0;
  89. /* Fill the scratchpad with random numbers */
  90. int i;
  91. for (i = 0; i < NSHOT_PER_TASK; i++)
  92. {
  93. float x = (float)(2.0*erand48(worker_xsub) - 1.0);
  94. float y = (float)(2.0*erand48(worker_xsub) - 1.0);
  95. float dist = x*x + y*y;
  96. if (dist < 1.0)
  97. local_cnt++;
  98. }
  99. /* Put the contribution of that task into the counter */
  100. unsigned long *cnt = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  101. *cnt = *cnt + local_cnt;
  102. }
  103. extern void pi_redux_cuda_kernel(float *x, float *y, unsigned n, unsigned long *shot_cnt);
  104. #ifdef STARPU_USE_CUDA
  105. static void pi_func_cuda(void *descr[], void *cl_arg __attribute__ ((unused)))
  106. {
  107. cudaError_t cures;
  108. curandStatus_t res;
  109. int workerid = starpu_worker_get_id();
  110. /* CURAND is a bit silly: it assumes that any error is fatal. Calling
  111. * cudaGetLastError resets the last error value. */
  112. cures = cudaGetLastError();
  113. // if (cures)
  114. // STARPU_CUDA_REPORT_ERROR(cures);
  115. /* Fill the scratchpad with random numbers. Note that both x and y
  116. * arrays are in stored the same vector. */
  117. float *scratchpad_xy = (float *)STARPU_VECTOR_GET_PTR(descr[0]);
  118. res = curandGenerateUniform(curandgens[workerid], scratchpad_xy, 2*NSHOT_PER_TASK);
  119. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  120. float *x = &scratchpad_xy[0];
  121. float *y = &scratchpad_xy[NSHOT_PER_TASK];
  122. unsigned long *shot_cnt = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  123. pi_redux_cuda_kernel(x, y, NSHOT_PER_TASK, shot_cnt);
  124. }
  125. #endif
  126. static struct starpu_codelet_t pi_cl = {
  127. .where = STARPU_CPU|STARPU_CUDA,
  128. .cpu_func = pi_func_cpu,
  129. #ifdef STARPU_USE_CUDA
  130. .cuda_func = pi_func_cuda,
  131. #endif
  132. .nbuffers = 2,
  133. .model = NULL
  134. };
  135. /*
  136. * Codelets to implement reduction
  137. */
  138. static void init_cpu_func(void *descr[], void *cl_arg)
  139. {
  140. unsigned long *val = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  141. *val = 0;
  142. }
  143. #ifdef STARPU_USE_CUDA
  144. static void init_cuda_func(void *descr[], void *cl_arg)
  145. {
  146. unsigned long *val = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  147. cudaMemset(val, 0, sizeof(unsigned long));
  148. cudaThreadSynchronize();
  149. }
  150. #endif
  151. static struct starpu_codelet_t init_codelet = {
  152. .where = STARPU_CPU|STARPU_CUDA,
  153. .cpu_func = init_cpu_func,
  154. #ifdef STARPU_USE_CUDA
  155. .cuda_func = init_cuda_func,
  156. #endif
  157. .nbuffers = 1
  158. };
  159. void redux_cpu_func(void *descr[], void *cl_arg)
  160. {
  161. unsigned long *a = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  162. unsigned long *b = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  163. *a = *a + *b;
  164. };
  165. static struct starpu_codelet_t redux_codelet = {
  166. .where = STARPU_CPU,
  167. .cpu_func = redux_cpu_func,
  168. .nbuffers = 2
  169. };
  170. /*
  171. * Main program
  172. */
  173. int main(int argc, char **argv)
  174. {
  175. unsigned i;
  176. parse_args(argc, argv);
  177. starpu_init(NULL);
  178. /* Launch a Random Number Generator (RNG) on each worker */
  179. starpu_execute_on_each_worker(init_rng, NULL, STARPU_CPU|STARPU_CUDA);
  180. /* Create a scratchpad data */
  181. starpu_data_handle xy_scratchpad_handle;
  182. starpu_vector_data_register(&xy_scratchpad_handle, -1, (uintptr_t)NULL,
  183. 2*NSHOT_PER_TASK, sizeof(float));
  184. /* Create a variable that will be used to count the number of shots
  185. * that actually hit the unit circle when shooting randomly in
  186. * [-1,1]^2. */
  187. unsigned long shot_cnt = 0;
  188. starpu_data_handle shot_cnt_handle;
  189. starpu_variable_data_register(&shot_cnt_handle, 0,
  190. (uintptr_t)&shot_cnt, sizeof(shot_cnt));
  191. starpu_data_set_reduction_methods(shot_cnt_handle,
  192. &redux_codelet, &init_codelet);
  193. struct timeval start;
  194. struct timeval end;
  195. gettimeofday(&start, NULL);
  196. for (i = 0; i < ntasks; i++)
  197. {
  198. struct starpu_task *task = starpu_task_create();
  199. task->cl = &pi_cl;
  200. task->buffers[0].handle = xy_scratchpad_handle;
  201. task->buffers[0].mode = STARPU_SCRATCH;
  202. task->buffers[1].handle = shot_cnt_handle;
  203. task->buffers[1].mode = STARPU_REDUX;
  204. int ret = starpu_task_submit(task);
  205. STARPU_ASSERT(!ret);
  206. }
  207. starpu_data_acquire(shot_cnt_handle, STARPU_R);
  208. gettimeofday(&end, NULL);
  209. double timing = (double)((end.tv_sec - start.tv_sec)*1000000 + (end.tv_usec - start.tv_usec));
  210. /* Total surface : Pi * r^ 2 = Pi*1^2, total square surface : 2^2 = 4,
  211. * probability to impact the disk: pi/4 */
  212. unsigned long total = ntasks*NSHOT_PER_TASK;
  213. double pi_approx = ((double)shot_cnt*4.0)/total;
  214. starpu_data_release(shot_cnt_handle);
  215. fprintf(stderr, "Pi approximation : %lf (%ld / %ld)\n", pi_approx, shot_cnt, total);
  216. fprintf(stderr, "Error %le \n", pi_approx - PI);
  217. fprintf(stderr, "Total time : %f ms\n", timing/1000.0);
  218. fprintf(stderr, "Speed : %f GShot/s\n", total/(1e3*timing));
  219. starpu_data_unregister(shot_cnt_handle);
  220. starpu_shutdown();
  221. if (abs(pi_approx - PI) > 1.0)
  222. return 1;
  223. return 0;
  224. }