pi_redux.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010-2015 Université de Bordeaux
  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 computes Pi by using drawing random coordinates (thanks to the sobol
  18. * generator) and check whether they fall within one quarter of a circle. The
  19. * proportion gives an approximation of Pi. For each task, we draw a number of
  20. * coordinates, and we gather the number of successful draws.
  21. *
  22. * This version uses reduction to optimize gathering the number of successful
  23. * draws.
  24. */
  25. #include <starpu.h>
  26. #include <stdlib.h>
  27. #include <math.h>
  28. #define FPRINTF(ofile, fmt, ...) do { if (!getenv("STARPU_SSILENT")) {fprintf(ofile, fmt, ## __VA_ARGS__); }} while(0)
  29. #define PI 3.14159265358979323846
  30. #if defined(STARPU_USE_CUDA) && !defined(STARPU_HAVE_CURAND)
  31. #warning CURAND is required to run that example on CUDA devices
  32. #endif
  33. #ifdef STARPU_HAVE_CURAND
  34. #include <cuda.h>
  35. #include <curand.h>
  36. #endif
  37. static unsigned long long nshot_per_task = 16*1024*1024ULL;
  38. /* default value */
  39. static unsigned long ntasks = 1024;
  40. static unsigned long ntasks_warmup = 0;
  41. static unsigned use_redux = 1;
  42. static unsigned do_warmup = 0;
  43. /*
  44. * Initialization of the Random Number Generators (RNG)
  45. */
  46. #ifdef STARPU_HAVE_CURAND
  47. /* RNG for the CURAND library */
  48. static curandGenerator_t curandgens[STARPU_NMAXWORKERS];
  49. #endif
  50. /* state for the erand48 function : note the huge padding to avoid false-sharing */
  51. #define PADDING 1024
  52. static unsigned short xsubi[STARPU_NMAXWORKERS*PADDING];
  53. static starpu_drand48_data randbuffer[STARPU_NMAXWORKERS*PADDING];
  54. /* Function to initialize the random number generator in the current worker */
  55. static void init_rng(void *arg STARPU_ATTRIBUTE_UNUSED)
  56. {
  57. #ifdef STARPU_HAVE_CURAND
  58. curandStatus_t res;
  59. #endif
  60. int workerid = starpu_worker_get_id();
  61. switch (starpu_worker_get_type(workerid))
  62. {
  63. case STARPU_CPU_WORKER:
  64. case STARPU_MIC_WORKER:
  65. case STARPU_SCC_WORKER:
  66. /* create a seed */
  67. starpu_srand48_r((long int)workerid, &randbuffer[PADDING*workerid]);
  68. xsubi[0 + PADDING*workerid] = (unsigned short)workerid;
  69. xsubi[1 + PADDING*workerid] = (unsigned short)workerid;
  70. xsubi[2 + PADDING*workerid] = (unsigned short)workerid;
  71. break;
  72. #ifdef STARPU_HAVE_CURAND
  73. case STARPU_CUDA_WORKER:
  74. /* Create a RNG */
  75. res = curandCreateGenerator(&curandgens[workerid],
  76. CURAND_RNG_PSEUDO_DEFAULT);
  77. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  78. /* Seed it with worker's id */
  79. res = curandSetPseudoRandomGeneratorSeed(curandgens[workerid],
  80. (unsigned long long)workerid);
  81. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  82. break;
  83. #endif
  84. default:
  85. STARPU_ABORT();
  86. break;
  87. }
  88. }
  89. /* The amount of work does not depend on the data size at all :) */
  90. static size_t size_base(struct starpu_task *task, unsigned nimpl)
  91. {
  92. return nshot_per_task;
  93. }
  94. static void parse_args(int argc, char **argv)
  95. {
  96. int i;
  97. for (i = 1; i < argc; i++)
  98. {
  99. if (strcmp(argv[i], "-ntasks") == 0)
  100. {
  101. char *argptr;
  102. ntasks = strtol(argv[++i], &argptr, 10);
  103. }
  104. if (strcmp(argv[i], "-nshot") == 0)
  105. {
  106. char *argptr;
  107. nshot_per_task = strtol(argv[++i], &argptr, 10);
  108. }
  109. if (strcmp(argv[i], "-noredux") == 0)
  110. {
  111. use_redux = 0;
  112. }
  113. if (strcmp(argv[i], "-warmup") == 0)
  114. {
  115. do_warmup = 1;
  116. ntasks_warmup = 8; /* arbitrary number of warmup tasks */
  117. }
  118. if (strcmp(argv[i], "-h") == 0 || strcmp(argv[i], "--help") == 0)
  119. {
  120. fprintf(stderr, "Usage: %s [-ntasks n] [-noredux] [-warmup] [-h]\n", argv[0]);
  121. exit(-1);
  122. }
  123. }
  124. }
  125. /*
  126. * Monte-carlo kernel
  127. */
  128. void pi_func_cpu(void *descr[], void *cl_arg STARPU_ATTRIBUTE_UNUSED)
  129. {
  130. int workerid = starpu_worker_get_id();
  131. unsigned short *worker_xsub;
  132. worker_xsub = &xsubi[PADDING*workerid];
  133. struct drand48_data *buffer;
  134. buffer = &randbuffer[PADDING*workerid];
  135. unsigned long local_cnt = 0;
  136. /* Fill the scratchpad with random numbers */
  137. unsigned i;
  138. for (i = 0; i < nshot_per_task; i++)
  139. {
  140. double randx, randy;
  141. starpu_erand48_r(worker_xsub, buffer, &randx);
  142. starpu_erand48_r(worker_xsub, buffer, &randy);
  143. double x = (2.0*randx - 1.0);
  144. double y = (2.0*randy - 1.0);
  145. double dist = x*x + y*y;
  146. if (dist < 1.0)
  147. local_cnt++;
  148. }
  149. /* Put the contribution of that task into the counter */
  150. unsigned long *cnt = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  151. *cnt = *cnt + local_cnt;
  152. }
  153. extern void pi_redux_cuda_kernel(float *x, float *y, unsigned n, unsigned long *shot_cnt);
  154. #ifdef STARPU_HAVE_CURAND
  155. static void pi_func_cuda(void *descr[], void *cl_arg STARPU_ATTRIBUTE_UNUSED)
  156. {
  157. curandStatus_t res;
  158. int workerid = starpu_worker_get_id();
  159. /* CURAND is a bit silly: it assumes that any error is fatal. Calling
  160. * cudaGetLastError resets the last error value. */
  161. (void) cudaGetLastError();
  162. /* Fill the scratchpad with random numbers. Note that both x and y
  163. * arrays are in stored the same vector. */
  164. float *scratchpad_xy = (float *)STARPU_VECTOR_GET_PTR(descr[0]);
  165. res = curandGenerateUniform(curandgens[workerid], scratchpad_xy, 2*nshot_per_task);
  166. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  167. float *x = &scratchpad_xy[0];
  168. float *y = &scratchpad_xy[nshot_per_task];
  169. unsigned long *shot_cnt = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  170. pi_redux_cuda_kernel(x, y, nshot_per_task, shot_cnt);
  171. }
  172. #endif
  173. static struct starpu_perfmodel pi_model =
  174. {
  175. .type = STARPU_HISTORY_BASED,
  176. .size_base = size_base,
  177. .symbol = "monte_carlo_pi_scratch"
  178. };
  179. static struct starpu_codelet pi_cl =
  180. {
  181. .cpu_funcs = {pi_func_cpu},
  182. .cpu_funcs_name = {"pi_func_cpu"},
  183. #ifdef STARPU_HAVE_CURAND
  184. .cuda_funcs = {pi_func_cuda},
  185. #endif
  186. .nbuffers = 2,
  187. .modes = {STARPU_SCRATCH, STARPU_RW},
  188. .model = &pi_model
  189. };
  190. static struct starpu_perfmodel pi_model_redux =
  191. {
  192. .type = STARPU_HISTORY_BASED,
  193. .size_base = size_base,
  194. .symbol = "monte_carlo_pi_scratch_redux"
  195. };
  196. static struct starpu_codelet pi_cl_redux =
  197. {
  198. .cpu_funcs = {pi_func_cpu},
  199. .cpu_funcs_name = {"pi_func_cpu"},
  200. #ifdef STARPU_HAVE_CURAND
  201. .cuda_funcs = {pi_func_cuda},
  202. #endif
  203. .nbuffers = 2,
  204. .modes = {STARPU_SCRATCH, STARPU_REDUX},
  205. .model = &pi_model_redux
  206. };
  207. /*
  208. * Codelets to implement reduction
  209. */
  210. void init_cpu_func(void *descr[], void *cl_arg)
  211. {
  212. unsigned long *val = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  213. *val = 0;
  214. }
  215. #ifdef STARPU_HAVE_CURAND
  216. static void init_cuda_func(void *descr[], void *cl_arg)
  217. {
  218. unsigned long *val = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  219. cudaMemsetAsync(val, 0, sizeof(unsigned long), starpu_cuda_get_local_stream());
  220. }
  221. #endif
  222. static struct starpu_codelet init_codelet =
  223. {
  224. .cpu_funcs = {init_cpu_func},
  225. .cpu_funcs_name = {"init_cpu_func"},
  226. #ifdef STARPU_HAVE_CURAND
  227. .cuda_funcs = {init_cuda_func},
  228. .cuda_flags = {STARPU_CUDA_ASYNC},
  229. #endif
  230. .modes = {STARPU_W},
  231. .nbuffers = 1
  232. };
  233. #ifdef STARPU_HAVE_CURAND
  234. /* Dummy implementation of the addition of two unsigned longs in CUDA */
  235. static void redux_cuda_func(void *descr[], void *cl_arg)
  236. {
  237. unsigned long *d_a = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  238. unsigned long *d_b = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  239. unsigned long h_a, h_b;
  240. cudaMemcpyAsync(&h_a, d_a, sizeof(h_a), cudaMemcpyDeviceToHost, starpu_cuda_get_local_stream());
  241. cudaMemcpyAsync(&h_b, d_b, sizeof(h_b), cudaMemcpyDeviceToHost, starpu_cuda_get_local_stream());
  242. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  243. h_a += h_b;
  244. cudaMemcpyAsync(d_a, &h_a, sizeof(h_a), cudaMemcpyHostToDevice, starpu_cuda_get_local_stream());
  245. }
  246. #endif
  247. void redux_cpu_func(void *descr[], void *cl_arg)
  248. {
  249. unsigned long *a = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  250. unsigned long *b = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  251. *a = *a + *b;
  252. }
  253. static struct starpu_codelet redux_codelet =
  254. {
  255. .cpu_funcs = {redux_cpu_func},
  256. .cpu_funcs_name = {"redux_cpu_func"},
  257. #ifdef STARPU_HAVE_CURAND
  258. .cuda_funcs = {redux_cuda_func},
  259. .cuda_flags = {STARPU_CUDA_ASYNC},
  260. #endif
  261. .modes = {STARPU_RW, STARPU_R},
  262. .nbuffers = 2
  263. };
  264. /*
  265. * Main program
  266. */
  267. int main(int argc, char **argv)
  268. {
  269. unsigned i;
  270. int ret;
  271. /* Not supported yet */
  272. if (starpu_get_env_number_default("STARPU_GLOBAL_ARBITER", 0) > 0)
  273. return 77;
  274. parse_args(argc, argv);
  275. ret = starpu_init(NULL);
  276. if (ret == -ENODEV) return 77;
  277. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  278. /* Launch a Random Number Generator (RNG) on each worker */
  279. starpu_execute_on_each_worker(init_rng, NULL, STARPU_CPU|STARPU_CUDA);
  280. /* Create a scratchpad data */
  281. starpu_data_handle_t xy_scratchpad_handle;
  282. starpu_vector_data_register(&xy_scratchpad_handle, -1, (uintptr_t)NULL,
  283. 2*nshot_per_task, sizeof(float));
  284. /* Create a variable that will be used to count the number of shots
  285. * that actually hit the unit circle when shooting randomly in
  286. * [-1,1]^2. */
  287. unsigned long shot_cnt = 0;
  288. starpu_data_handle_t shot_cnt_handle;
  289. starpu_variable_data_register(&shot_cnt_handle, STARPU_MAIN_RAM,
  290. (uintptr_t)&shot_cnt, sizeof(shot_cnt));
  291. starpu_data_set_reduction_methods(shot_cnt_handle,
  292. &redux_codelet, &init_codelet);
  293. double start;
  294. double end;
  295. for (i = 0; i < ntasks_warmup; i++)
  296. {
  297. struct starpu_task *task = starpu_task_create();
  298. task->cl = use_redux?&pi_cl_redux:&pi_cl;
  299. task->handles[0] = xy_scratchpad_handle;
  300. task->handles[1] = shot_cnt_handle;
  301. ret = starpu_task_submit(task);
  302. STARPU_ASSERT(!ret);
  303. }
  304. start = starpu_timing_now();
  305. for (i = 0; i < ntasks; i++)
  306. {
  307. struct starpu_task *task = starpu_task_create();
  308. task->cl = use_redux?&pi_cl_redux:&pi_cl;
  309. task->handles[0] = xy_scratchpad_handle;
  310. task->handles[1] = shot_cnt_handle;
  311. ret = starpu_task_submit(task);
  312. STARPU_ASSERT(!ret);
  313. }
  314. starpu_data_unregister(shot_cnt_handle);
  315. starpu_data_unregister(xy_scratchpad_handle);
  316. end = starpu_timing_now();
  317. double timing = end - start;
  318. /* Total surface : Pi * r^ 2 = Pi*1^2, total square surface : 2^2 = 4,
  319. * probability to impact the disk: pi/4 */
  320. unsigned long total = (ntasks + ntasks_warmup)*nshot_per_task;
  321. double pi_approx = ((double)shot_cnt*4.0)/total;
  322. FPRINTF(stderr, "Reductions? %s\n", use_redux?"yes":"no");
  323. FPRINTF(stderr, "Pi approximation : %f (%ld / %ld)\n", pi_approx, shot_cnt, total);
  324. FPRINTF(stderr, "Error %e \n", pi_approx - PI);
  325. FPRINTF(stderr, "Total time : %f ms\n", timing/1000.0);
  326. FPRINTF(stderr, "Speed : %f GShot/s\n", total/(1e3*timing));
  327. starpu_shutdown();
  328. if (fabs(pi_approx - PI) > 1.0)
  329. return 1;
  330. return 0;
  331. }