pi.c 6.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010,2011,2013-2015 Université de Bordeaux
  4. * Copyright (C) 2012,2013 Inria
  5. * Copyright (C) 2010 Mehdi Juhoor
  6. * Copyright (C) 2010-2013,2015-2017 CNRS
  7. *
  8. * StarPU is free software; you can redistribute it and/or modify
  9. * it under the terms of the GNU Lesser General Public License as published by
  10. * the Free Software Foundation; either version 2.1 of the License, or (at
  11. * your option) any later version.
  12. *
  13. * StarPU is distributed in the hope that it will be useful, but
  14. * WITHOUT ANY WARRANTY; without even the implied warranty of
  15. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  16. *
  17. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  18. */
  19. /*
  20. * This computes Pi by using drawing random coordinates (thanks to the sobol
  21. * generator) and check whether they fall within one quarter of a circle. The
  22. * proportion gives an approximation of Pi. For each task, we draw a number of
  23. * coordinates, and we gather the number of successful draws.
  24. *
  25. * TODO: use curandGenerateUniform instead of the sobol generator, like pi_redux.c does
  26. */
  27. #include "SobolQRNG/sobol.h"
  28. #include "SobolQRNG/sobol_gold.h"
  29. #include "pi.h"
  30. #ifdef STARPU_USE_CUDA
  31. void cuda_kernel(void **descr, void *cl_arg);
  32. #endif
  33. #define FPRINTF(ofile, fmt, ...) do { if (!getenv("STARPU_SSILENT")) {fprintf(ofile, fmt, ## __VA_ARGS__); }} while(0)
  34. /* default value */
  35. static unsigned ntasks = 1024;
  36. static unsigned long long nshot_per_task = 16*1024*1024ULL;
  37. void cpu_kernel(void *descr[], void *cl_arg)
  38. {
  39. (void)cl_arg;
  40. unsigned *directions = (unsigned *)STARPU_VECTOR_GET_PTR(descr[0]);
  41. unsigned nx = nshot_per_task;
  42. TYPE *random_numbers = malloc(2*nx*sizeof(TYPE));
  43. sobolCPU(2*nx/n_dimensions, n_dimensions, directions, random_numbers);
  44. TYPE *random_numbers_x = &random_numbers[0];
  45. TYPE *random_numbers_y = &random_numbers[nx];
  46. unsigned current_cnt = 0;
  47. unsigned i;
  48. for (i = 0; i < nx; i++)
  49. {
  50. TYPE x = random_numbers_x[i];
  51. TYPE y = random_numbers_y[i];
  52. TYPE dist = (x*x + y*y);
  53. unsigned success = (dist <= 1.0);
  54. current_cnt += success;
  55. }
  56. unsigned *cnt = (unsigned *)STARPU_VECTOR_GET_PTR(descr[1]);
  57. *cnt = current_cnt;
  58. free(random_numbers);
  59. }
  60. /* The amount of work does not depend on the data size at all :) */
  61. static size_t size_base(struct starpu_task *task, unsigned nimpl)
  62. {
  63. (void)task;
  64. (void)nimpl;
  65. return nshot_per_task;
  66. }
  67. static void parse_args(int argc, char **argv)
  68. {
  69. int i;
  70. for (i = 1; i < argc; i++)
  71. {
  72. if (strcmp(argv[i], "-ntasks") == 0)
  73. {
  74. char *argptr;
  75. ntasks = strtol(argv[++i], &argptr, 10);
  76. }
  77. if (strcmp(argv[i], "-nshot") == 0)
  78. {
  79. char *argptr;
  80. nshot_per_task = strtol(argv[++i], &argptr, 10);
  81. }
  82. if (strcmp(argv[i], "-h") == 0 || strcmp(argv[i], "--help") == 0)
  83. {
  84. fprintf(stderr,"Usage: %s [options...]\n", argv[0]);
  85. fprintf(stderr,"\n");
  86. fprintf(stderr,"Options:\n");
  87. fprintf(stderr,"-ntasks <n> select the number of tasks\n");
  88. fprintf(stderr,"-nshot <n> select the number of shot per task\n");
  89. exit(0);
  90. }
  91. }
  92. }
  93. static struct starpu_perfmodel model =
  94. {
  95. .type = STARPU_HISTORY_BASED,
  96. .size_base = size_base,
  97. .symbol = "monte_carlo_pi"
  98. };
  99. static struct starpu_codelet pi_cl =
  100. {
  101. .cpu_funcs = {cpu_kernel},
  102. .cpu_funcs_name = {"cpu_kernel"},
  103. #ifdef STARPU_USE_CUDA
  104. .cuda_funcs = {cuda_kernel},
  105. #endif
  106. .nbuffers = 2,
  107. .modes = {STARPU_R, STARPU_W},
  108. .model = &model
  109. };
  110. int main(int argc, char **argv)
  111. {
  112. unsigned i;
  113. int ret;
  114. parse_args(argc, argv);
  115. ret = starpu_init(NULL);
  116. if (ret == -ENODEV) return 77;
  117. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  118. /* Initialize the random number generator */
  119. unsigned *sobol_qrng_directions = malloc(n_dimensions*n_directions*sizeof(unsigned));
  120. STARPU_ASSERT(sobol_qrng_directions);
  121. initSobolDirectionVectors(n_dimensions, sobol_qrng_directions);
  122. /* Any worker may use that array now */
  123. starpu_data_handle_t sobol_qrng_direction_handle;
  124. starpu_vector_data_register(&sobol_qrng_direction_handle, STARPU_MAIN_RAM,
  125. (uintptr_t)sobol_qrng_directions, n_dimensions*n_directions, sizeof(unsigned));
  126. unsigned *cnt_array = calloc(ntasks, sizeof(unsigned));
  127. STARPU_ASSERT(cnt_array);
  128. starpu_data_handle_t cnt_array_handle;
  129. starpu_vector_data_register(&cnt_array_handle, STARPU_MAIN_RAM, (uintptr_t)cnt_array, ntasks, sizeof(unsigned));
  130. /* Use a write-through policy : when the data is modified on an
  131. * accelerator, we know that it will only be modified once and be
  132. * accessed by the CPU later on */
  133. starpu_data_set_wt_mask(cnt_array_handle, (1<<0));
  134. struct starpu_data_filter f =
  135. {
  136. .filter_func = starpu_vector_filter_block,
  137. .nchildren = ntasks
  138. };
  139. starpu_data_partition(cnt_array_handle, &f);
  140. double start;
  141. double end;
  142. start = starpu_timing_now();
  143. for (i = 0; i < ntasks; i++)
  144. {
  145. struct starpu_task *task = starpu_task_create();
  146. task->cl = &pi_cl;
  147. STARPU_ASSERT(starpu_data_get_sub_data(cnt_array_handle, 1, i));
  148. task->handles[0] = sobol_qrng_direction_handle;
  149. task->handles[1] = starpu_data_get_sub_data(cnt_array_handle, 1, i);
  150. ret = starpu_task_submit(task);
  151. STARPU_ASSERT(!ret);
  152. }
  153. starpu_task_wait_for_all();
  154. /* Get the cnt_array back in main memory */
  155. starpu_data_unpartition(cnt_array_handle, STARPU_MAIN_RAM);
  156. starpu_data_unregister(cnt_array_handle);
  157. starpu_data_unregister(sobol_qrng_direction_handle);
  158. /* Count the total number of entries */
  159. unsigned long total_cnt = 0;
  160. for (i = 0; i < ntasks; i++)
  161. total_cnt += cnt_array[i];
  162. end = starpu_timing_now();
  163. double timing = end - start;
  164. unsigned long total_shot_cnt = ntasks * nshot_per_task;
  165. /* Total surface : Pi * r^ 2 = Pi*1^2, total square surface : 2^2 = 4, probability to impact the disk: pi/4 */
  166. FPRINTF(stderr, "Pi approximation : %f (%lu / %lu)\n", ((TYPE)total_cnt*4)/(total_shot_cnt), total_cnt, total_shot_cnt);
  167. FPRINTF(stderr, "Total time : %f ms\n", timing/1000.0);
  168. FPRINTF(stderr, "Speed : %f GShot/s\n", total_shot_cnt/(1e3*timing));
  169. if (!getenv("STARPU_SSILENT")) starpu_codelet_display_stats(&pi_cl);
  170. starpu_shutdown();
  171. return 0;
  172. }