pi_redux.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010-2020 Université de Bordeaux, CNRS (LaBRI UMR 5800), Inria
  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)
  56. {
  57. (void)arg;
  58. #ifdef STARPU_HAVE_CURAND
  59. curandStatus_t res;
  60. #endif
  61. int workerid = starpu_worker_get_id_check();
  62. switch (starpu_worker_get_type(workerid))
  63. {
  64. case STARPU_CPU_WORKER:
  65. case STARPU_MIC_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. (void)task;
  93. (void)nimpl;
  94. return nshot_per_task;
  95. }
  96. static void parse_args(int argc, char **argv)
  97. {
  98. int i;
  99. for (i = 1; i < argc; i++)
  100. {
  101. if (strcmp(argv[i], "-ntasks") == 0)
  102. {
  103. char *argptr;
  104. ntasks = strtol(argv[++i], &argptr, 10);
  105. }
  106. if (strcmp(argv[i], "-nshot") == 0)
  107. {
  108. char *argptr;
  109. nshot_per_task = strtol(argv[++i], &argptr, 10);
  110. }
  111. if (strcmp(argv[i], "-noredux") == 0)
  112. {
  113. use_redux = 0;
  114. }
  115. if (strcmp(argv[i], "-warmup") == 0)
  116. {
  117. do_warmup = 1;
  118. ntasks_warmup = 8; /* arbitrary number of warmup tasks */
  119. }
  120. if (strcmp(argv[i], "-h") == 0 || strcmp(argv[i], "--help") == 0)
  121. {
  122. fprintf(stderr, "Usage: %s [-ntasks n] [-noredux] [-warmup] [-h]\n", argv[0]);
  123. exit(-1);
  124. }
  125. }
  126. }
  127. /*
  128. * Monte-carlo kernel
  129. */
  130. void pi_func_cpu(void *descr[], void *cl_arg)
  131. {
  132. (void)cl_arg;
  133. int workerid = starpu_worker_get_id_check();
  134. unsigned short *worker_xsub;
  135. worker_xsub = &xsubi[PADDING*workerid];
  136. starpu_drand48_data *buffer;
  137. buffer = &randbuffer[PADDING*workerid];
  138. unsigned long local_cnt = 0;
  139. /* Fill the scratchpad with random numbers */
  140. unsigned i;
  141. for (i = 0; i < nshot_per_task; i++)
  142. {
  143. double randx, randy;
  144. starpu_erand48_r(worker_xsub, buffer, &randx);
  145. starpu_erand48_r(worker_xsub, buffer, &randy);
  146. double x = (2.0*randx - 1.0);
  147. double y = (2.0*randy - 1.0);
  148. double dist = x*x + y*y;
  149. if (dist < 1.0)
  150. local_cnt++;
  151. }
  152. /* Put the contribution of that task into the counter */
  153. unsigned long *cnt = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  154. *cnt = *cnt + local_cnt;
  155. }
  156. extern void pi_redux_cuda_kernel(float *x, float *y, unsigned n, unsigned long *shot_cnt);
  157. #ifdef STARPU_HAVE_CURAND
  158. static void pi_func_cuda(void *descr[], void *cl_arg)
  159. {
  160. (void)cl_arg;
  161. curandStatus_t res;
  162. int workerid = starpu_worker_get_id_check();
  163. /* CURAND is a bit silly: it assumes that any error is fatal. Calling
  164. * cudaGetLastError resets the last error value. */
  165. (void) cudaGetLastError();
  166. /* Fill the scratchpad with random numbers. Note that both x and y
  167. * arrays are in stored the same vector. */
  168. float *scratchpad_xy = (float *)STARPU_VECTOR_GET_PTR(descr[0]);
  169. res = curandGenerateUniform(curandgens[workerid], scratchpad_xy, 2*nshot_per_task);
  170. STARPU_ASSERT(res == CURAND_STATUS_SUCCESS);
  171. float *x = &scratchpad_xy[0];
  172. float *y = &scratchpad_xy[nshot_per_task];
  173. unsigned long *shot_cnt = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  174. pi_redux_cuda_kernel(x, y, nshot_per_task, shot_cnt);
  175. }
  176. #endif
  177. static struct starpu_perfmodel pi_model =
  178. {
  179. .type = STARPU_HISTORY_BASED,
  180. .size_base = size_base,
  181. .symbol = "monte_carlo_pi_scratch"
  182. };
  183. static struct starpu_codelet pi_cl =
  184. {
  185. .cpu_funcs = {pi_func_cpu},
  186. .cpu_funcs_name = {"pi_func_cpu"},
  187. #ifdef STARPU_HAVE_CURAND
  188. .cuda_funcs = {pi_func_cuda},
  189. #endif
  190. .nbuffers = 2,
  191. .modes = {STARPU_SCRATCH, STARPU_RW},
  192. .model = &pi_model
  193. };
  194. static struct starpu_perfmodel pi_model_redux =
  195. {
  196. .type = STARPU_HISTORY_BASED,
  197. .size_base = size_base,
  198. .symbol = "monte_carlo_pi_scratch_redux"
  199. };
  200. static struct starpu_codelet pi_cl_redux =
  201. {
  202. .cpu_funcs = {pi_func_cpu},
  203. .cpu_funcs_name = {"pi_func_cpu"},
  204. #ifdef STARPU_HAVE_CURAND
  205. .cuda_funcs = {pi_func_cuda},
  206. #endif
  207. .nbuffers = 2,
  208. .modes = {STARPU_SCRATCH, STARPU_REDUX},
  209. .model = &pi_model_redux
  210. };
  211. /*
  212. * Codelets to implement reduction
  213. */
  214. void init_cpu_func(void *descr[], void *cl_arg)
  215. {
  216. (void)cl_arg;
  217. unsigned long *val = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  218. *val = 0;
  219. }
  220. #ifdef STARPU_HAVE_CURAND
  221. static void init_cuda_func(void *descr[], void *cl_arg)
  222. {
  223. (void)cl_arg;
  224. unsigned long *val = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  225. cudaMemsetAsync(val, 0, sizeof(unsigned long), starpu_cuda_get_local_stream());
  226. }
  227. #endif
  228. static struct starpu_codelet init_codelet =
  229. {
  230. .cpu_funcs = {init_cpu_func},
  231. .cpu_funcs_name = {"init_cpu_func"},
  232. #ifdef STARPU_HAVE_CURAND
  233. .cuda_funcs = {init_cuda_func},
  234. .cuda_flags = {STARPU_CUDA_ASYNC},
  235. #endif
  236. .modes = {STARPU_W},
  237. .nbuffers = 1
  238. };
  239. #ifdef STARPU_HAVE_CURAND
  240. /* Dummy implementation of the addition of two unsigned longs in CUDA */
  241. static void redux_cuda_func(void *descr[], void *cl_arg)
  242. {
  243. (void)cl_arg;
  244. unsigned long *d_a = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  245. unsigned long *d_b = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  246. unsigned long h_a, h_b;
  247. cudaMemcpyAsync(&h_a, d_a, sizeof(h_a), cudaMemcpyDeviceToHost, starpu_cuda_get_local_stream());
  248. cudaMemcpyAsync(&h_b, d_b, sizeof(h_b), cudaMemcpyDeviceToHost, starpu_cuda_get_local_stream());
  249. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  250. h_a += h_b;
  251. cudaMemcpyAsync(d_a, &h_a, sizeof(h_a), cudaMemcpyHostToDevice, starpu_cuda_get_local_stream());
  252. }
  253. #endif
  254. void redux_cpu_func(void *descr[], void *cl_arg)
  255. {
  256. (void)cl_arg;
  257. unsigned long *a = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[0]);
  258. unsigned long *b = (unsigned long *)STARPU_VARIABLE_GET_PTR(descr[1]);
  259. *a = *a + *b;
  260. }
  261. static struct starpu_codelet redux_codelet =
  262. {
  263. .cpu_funcs = {redux_cpu_func},
  264. .cpu_funcs_name = {"redux_cpu_func"},
  265. #ifdef STARPU_HAVE_CURAND
  266. .cuda_funcs = {redux_cuda_func},
  267. .cuda_flags = {STARPU_CUDA_ASYNC},
  268. #endif
  269. .modes = {STARPU_RW, STARPU_R},
  270. .nbuffers = 2
  271. };
  272. /*
  273. * Main program
  274. */
  275. int main(int argc, char **argv)
  276. {
  277. unsigned i;
  278. int ret;
  279. /* Not supported yet */
  280. if (starpu_get_env_number_default("STARPU_GLOBAL_ARBITER", 0) > 0)
  281. return 77;
  282. parse_args(argc, argv);
  283. ret = starpu_init(NULL);
  284. if (ret == -ENODEV) return 77;
  285. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  286. /* Launch a Random Number Generator (RNG) on each worker */
  287. starpu_execute_on_each_worker(init_rng, NULL, STARPU_CPU|STARPU_CUDA);
  288. /* Create a scratchpad data */
  289. starpu_data_handle_t xy_scratchpad_handle;
  290. starpu_vector_data_register(&xy_scratchpad_handle, -1, (uintptr_t)NULL,
  291. 2*nshot_per_task, sizeof(float));
  292. /* Create a variable that will be used to count the number of shots
  293. * that actually hit the unit circle when shooting randomly in
  294. * [-1,1]^2. */
  295. unsigned long shot_cnt = 0;
  296. starpu_data_handle_t shot_cnt_handle;
  297. starpu_variable_data_register(&shot_cnt_handle, STARPU_MAIN_RAM,
  298. (uintptr_t)&shot_cnt, sizeof(shot_cnt));
  299. starpu_data_set_reduction_methods(shot_cnt_handle,
  300. &redux_codelet, &init_codelet);
  301. double start;
  302. double end;
  303. for (i = 0; i < ntasks_warmup; i++)
  304. {
  305. struct starpu_task *task = starpu_task_create();
  306. task->cl = use_redux?&pi_cl_redux:&pi_cl;
  307. task->handles[0] = xy_scratchpad_handle;
  308. task->handles[1] = shot_cnt_handle;
  309. ret = starpu_task_submit(task);
  310. STARPU_ASSERT(!ret);
  311. }
  312. start = starpu_timing_now();
  313. for (i = 0; i < ntasks; i++)
  314. {
  315. struct starpu_task *task = starpu_task_create();
  316. task->cl = use_redux?&pi_cl_redux:&pi_cl;
  317. task->handles[0] = xy_scratchpad_handle;
  318. task->handles[1] = shot_cnt_handle;
  319. ret = starpu_task_submit(task);
  320. STARPU_ASSERT(!ret);
  321. }
  322. starpu_data_unregister(shot_cnt_handle);
  323. starpu_data_unregister(xy_scratchpad_handle);
  324. end = starpu_timing_now();
  325. double timing = end - start;
  326. /* Total surface : Pi * r^ 2 = Pi*1^2, total square surface : 2^2 = 4,
  327. * probability to impact the disk: pi/4 */
  328. unsigned long total = (ntasks + ntasks_warmup)*nshot_per_task;
  329. double pi_approx = ((double)shot_cnt*4.0)/total;
  330. FPRINTF(stderr, "Reductions? %s\n", use_redux?"yes":"no");
  331. FPRINTF(stderr, "Pi approximation : %f (%lu / %lu)\n", pi_approx, shot_cnt, total);
  332. FPRINTF(stderr, "Error %e \n", pi_approx - PI);
  333. FPRINTF(stderr, "Total time : %f ms\n", timing/1000.0);
  334. FPRINTF(stderr, "Speed : %f GShot/s\n", total/(1e3*timing));
  335. starpu_shutdown();
  336. if (fabs(pi_approx - PI) > 1.0)
  337. return 1;
  338. return 0;
  339. }