pi.c 6.0 KB

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