pi_redux.c 7.3 KB

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