pi_redux.c 7.1 KB

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