pi.c 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010-2011 Université de Bordeaux 1
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011 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, args ...) do { if (!getenv("STARPU_SSILENT")) {fprintf(ofile, fmt, ##args); }} while(0)
  26. /* default value */
  27. static unsigned ntasks = 1024;
  28. static 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. static void parse_args(int argc, char **argv)
  51. {
  52. int i;
  53. for (i = 1; i < argc; i++) {
  54. if (strcmp(argv[i], "-ntasks") == 0) {
  55. char *argptr;
  56. ntasks = strtol(argv[++i], &argptr, 10);
  57. }
  58. }
  59. }
  60. int main(int argc, char **argv)
  61. {
  62. unsigned i;
  63. parse_args(argc, argv);
  64. starpu_init(NULL);
  65. /* Initialize the random number generator */
  66. unsigned *sobol_qrng_directions = malloc(n_dimensions*n_directions*sizeof(unsigned));
  67. STARPU_ASSERT(sobol_qrng_directions);
  68. initSobolDirectionVectors(n_dimensions, sobol_qrng_directions);
  69. /* Any worker may use that array now */
  70. starpu_data_handle_t sobol_qrng_direction_handle;
  71. starpu_vector_data_register(&sobol_qrng_direction_handle, 0,
  72. (uintptr_t)sobol_qrng_directions, n_dimensions*n_directions, sizeof(unsigned));
  73. unsigned *cnt_array = malloc(ntasks*sizeof(unsigned));
  74. STARPU_ASSERT(cnt_array);
  75. starpu_data_handle_t cnt_array_handle;
  76. starpu_vector_data_register(&cnt_array_handle, 0, (uintptr_t)cnt_array, ntasks, sizeof(unsigned));
  77. /* Use a write-through policy : when the data is modified on an
  78. * accelerator, we know that it will only be modified once and be
  79. * accessed by the CPU later on */
  80. starpu_data_set_wt_mask(cnt_array_handle, (1<<0));
  81. struct starpu_data_filter f = {
  82. .filter_func = starpu_block_filter_func_vector,
  83. .nchildren = ntasks
  84. };
  85. starpu_data_partition(cnt_array_handle, &f);
  86. static struct starpu_perfmodel model = {
  87. .type = STARPU_HISTORY_BASED,
  88. .symbol = "monte_carlo_pi"
  89. };
  90. struct starpu_codelet cl = {
  91. .where = STARPU_CPU|STARPU_CUDA,
  92. .cpu_funcs = {cpu_kernel, NULL},
  93. #ifdef STARPU_USE_CUDA
  94. .cuda_funcs = {cuda_kernel, NULL},
  95. #endif
  96. .nbuffers = 2,
  97. .model = &model
  98. };
  99. struct timeval start;
  100. struct timeval end;
  101. gettimeofday(&start, NULL);
  102. for (i = 0; i < ntasks; i++)
  103. {
  104. struct starpu_task *task = starpu_task_create();
  105. task->cl = &cl;
  106. STARPU_ASSERT(starpu_data_get_sub_data(cnt_array_handle, 1, i));
  107. task->buffers[0].handle = sobol_qrng_direction_handle;
  108. task->buffers[0].mode = STARPU_R;
  109. task->buffers[1].handle = starpu_data_get_sub_data(cnt_array_handle, 1, i);
  110. task->buffers[1].mode = STARPU_W;
  111. int ret = starpu_task_submit(task);
  112. STARPU_ASSERT(!ret);
  113. }
  114. starpu_task_wait_for_all();
  115. /* Get the cnt_array back in main memory */
  116. starpu_data_unpartition(cnt_array_handle, 0);
  117. starpu_data_unregister(cnt_array_handle);
  118. starpu_data_unregister(sobol_qrng_direction_handle);
  119. /* Count the total number of entries */
  120. unsigned long total_cnt = 0;
  121. for (i = 0; i < ntasks; i++)
  122. total_cnt += cnt_array[i];
  123. gettimeofday(&end, NULL);
  124. double timing = (double)((end.tv_sec - start.tv_sec)*1000000 + (end.tv_usec - start.tv_usec));
  125. unsigned long total_shot_cnt = ntasks * NSHOT_PER_TASK;
  126. /* Total surface : Pi * r^ 2 = Pi*1^2, total square surface : 2^2 = 4, probability to impact the disk: pi/4 */
  127. FPRINTF(stderr, "Pi approximation : %f (%ld / %ld)\n", ((TYPE)total_cnt*4)/(total_shot_cnt), total_cnt, total_shot_cnt);
  128. FPRINTF(stderr, "Total time : %f ms\n", timing/1000.0);
  129. FPRINTF(stderr, "Speed : %f GShot/s\n", total_shot_cnt/(1e3*timing));
  130. if (!getenv("STARPU_SSILENT")) starpu_display_codelet_stats(&cl);
  131. starpu_shutdown();
  132. return 0;
  133. }