pi.c 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010-2011, 2013-2014 Université de Bordeaux
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011, 2012, 2013 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. #include "SobolQRNG/sobol.h"
  19. #include "SobolQRNG/sobol_gold.h"
  20. #include "pi.h"
  21. #ifdef STARPU_USE_CUDA
  22. void cuda_kernel(void **descr, void *cl_arg);
  23. #endif
  24. #define FPRINTF(ofile, fmt, ...) do { if (!getenv("STARPU_SSILENT")) {fprintf(ofile, fmt, ## __VA_ARGS__); }} while(0)
  25. /* default value */
  26. static unsigned ntasks = 1024;
  27. static unsigned long long nshot_per_task = 16*1024*1024ULL;
  28. void cpu_kernel(void *descr[], void *cl_arg)
  29. {
  30. unsigned *directions = (unsigned *)STARPU_VECTOR_GET_PTR(descr[0]);
  31. unsigned nx = nshot_per_task;
  32. TYPE *random_numbers = malloc(2*nx*sizeof(TYPE));
  33. sobolCPU(2*nx/n_dimensions, n_dimensions, directions, random_numbers);
  34. TYPE *random_numbers_x = &random_numbers[0];
  35. TYPE *random_numbers_y = &random_numbers[nx];
  36. unsigned current_cnt = 0;
  37. unsigned i;
  38. for (i = 0; i < nx; i++)
  39. {
  40. TYPE x = random_numbers_x[i];
  41. TYPE y = random_numbers_y[i];
  42. TYPE dist = (x*x + y*y);
  43. unsigned success = (dist <= 1.0);
  44. current_cnt += success;
  45. }
  46. unsigned *cnt = (unsigned *)STARPU_VECTOR_GET_PTR(descr[1]);
  47. *cnt = current_cnt;
  48. free(random_numbers);
  49. }
  50. /* The amount of work does not depend on the data size at all :) */
  51. static size_t size_base(struct starpu_task *task, unsigned nimpl)
  52. {
  53. return nshot_per_task;
  54. }
  55. static void parse_args(int argc, char **argv)
  56. {
  57. int i;
  58. for (i = 1; i < argc; i++)
  59. {
  60. if (strcmp(argv[i], "-ntasks") == 0)
  61. {
  62. char *argptr;
  63. ntasks = strtol(argv[++i], &argptr, 10);
  64. }
  65. if (strcmp(argv[i], "-nshot") == 0)
  66. {
  67. char *argptr;
  68. nshot_per_task = strtol(argv[++i], &argptr, 10);
  69. }
  70. if (strcmp(argv[i], "-h") == 0 || strcmp(argv[i], "--help") == 0)
  71. {
  72. fprintf(stderr,"Usage: %s [options...]\n", argv[0]);
  73. fprintf(stderr,"\n");
  74. fprintf(stderr,"Options:\n");
  75. fprintf(stderr,"-ntasks <n> select the number of tasks\n");
  76. fprintf(stderr,"-nshot <n> select the number of shot per task\n");
  77. exit(0);
  78. }
  79. }
  80. }
  81. static struct starpu_perfmodel model =
  82. {
  83. .type = STARPU_HISTORY_BASED,
  84. .size_base = size_base,
  85. .symbol = "monte_carlo_pi"
  86. };
  87. static struct starpu_codelet pi_cl =
  88. {
  89. .cpu_funcs = {cpu_kernel},
  90. .cpu_funcs_name = {"cpu_kernel"},
  91. #ifdef STARPU_USE_CUDA
  92. .cuda_funcs = {cuda_kernel},
  93. #endif
  94. .nbuffers = 2,
  95. .modes = {STARPU_R, STARPU_W},
  96. .model = &model
  97. };
  98. int main(int argc, char **argv)
  99. {
  100. unsigned i;
  101. int ret;
  102. parse_args(argc, argv);
  103. ret = starpu_init(NULL);
  104. if (ret == -ENODEV) return 77;
  105. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  106. /* Initialize the random number generator */
  107. unsigned *sobol_qrng_directions = malloc(n_dimensions*n_directions*sizeof(unsigned));
  108. STARPU_ASSERT(sobol_qrng_directions);
  109. initSobolDirectionVectors(n_dimensions, sobol_qrng_directions);
  110. /* Any worker may use that array now */
  111. starpu_data_handle_t sobol_qrng_direction_handle;
  112. starpu_vector_data_register(&sobol_qrng_direction_handle, STARPU_MAIN_RAM,
  113. (uintptr_t)sobol_qrng_directions, n_dimensions*n_directions, sizeof(unsigned));
  114. unsigned *cnt_array = malloc(ntasks*sizeof(unsigned));
  115. STARPU_ASSERT(cnt_array);
  116. starpu_data_handle_t cnt_array_handle;
  117. starpu_vector_data_register(&cnt_array_handle, STARPU_MAIN_RAM, (uintptr_t)cnt_array, ntasks, sizeof(unsigned));
  118. /* Use a write-through policy : when the data is modified on an
  119. * accelerator, we know that it will only be modified once and be
  120. * accessed by the CPU later on */
  121. starpu_data_set_wt_mask(cnt_array_handle, (1<<0));
  122. struct starpu_data_filter f =
  123. {
  124. .filter_func = starpu_vector_filter_block,
  125. .nchildren = ntasks
  126. };
  127. starpu_data_partition(cnt_array_handle, &f);
  128. double start;
  129. double end;
  130. start = starpu_timing_now();
  131. for (i = 0; i < ntasks; i++)
  132. {
  133. struct starpu_task *task = starpu_task_create();
  134. task->cl = &pi_cl;
  135. STARPU_ASSERT(starpu_data_get_sub_data(cnt_array_handle, 1, i));
  136. task->handles[0] = sobol_qrng_direction_handle;
  137. task->handles[1] = starpu_data_get_sub_data(cnt_array_handle, 1, i);
  138. ret = starpu_task_submit(task);
  139. STARPU_ASSERT(!ret);
  140. }
  141. starpu_task_wait_for_all();
  142. /* Get the cnt_array back in main memory */
  143. starpu_data_unpartition(cnt_array_handle, STARPU_MAIN_RAM);
  144. starpu_data_unregister(cnt_array_handle);
  145. starpu_data_unregister(sobol_qrng_direction_handle);
  146. /* Count the total number of entries */
  147. unsigned long total_cnt = 0;
  148. for (i = 0; i < ntasks; i++)
  149. total_cnt += cnt_array[i];
  150. end = starpu_timing_now();
  151. double timing = end - start;
  152. unsigned long total_shot_cnt = ntasks * nshot_per_task;
  153. /* Total surface : Pi * r^ 2 = Pi*1^2, total square surface : 2^2 = 4, probability to impact the disk: pi/4 */
  154. FPRINTF(stderr, "Pi approximation : %f (%ld / %ld)\n", ((TYPE)total_cnt*4)/(total_shot_cnt), total_cnt, total_shot_cnt);
  155. FPRINTF(stderr, "Total time : %f ms\n", timing/1000.0);
  156. FPRINTF(stderr, "Speed : %f GShot/s\n", total_shot_cnt/(1e3*timing));
  157. if (!getenv("STARPU_SSILENT")) starpu_codelet_display_stats(&pi_cl);
  158. starpu_shutdown();
  159. return 0;
  160. }