pi.c 5.7 KB

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