mandelbrot_between.c 6.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2020 Université de Bordeaux, CNRS (LaBRI UMR 5800), Inria
  4. * Copyright (C) 2019 Mael Keryell
  5. *
  6. * StarPU is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU Lesser General Public License as published by
  8. * the Free Software Foundation; either version 2.1 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * StarPU is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  14. *
  15. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  16. */
  17. #include <stdio.h>
  18. #include <stdlib.h>
  19. #include <starpu.h>
  20. #include <stdint.h>
  21. #include "../includes/display.h"
  22. void mandelbrot(void **, void *);
  23. void CUDA_mandelbrot(void **, void *);
  24. void test(void **, void *); /* Function used to test on my matrix, in the cpu_test_with_generated.c file. */
  25. struct Params
  26. {
  27. float cr;
  28. float ci;
  29. unsigned taskx;
  30. unsigned tasky;
  31. unsigned width;
  32. unsigned height;
  33. };
  34. struct starpu_codelet cl =
  35. {
  36. /* .cpu_funcs = {test}, */
  37. .cpu_funcs = {mandelbrot},
  38. .cuda_funcs = {CUDA_mandelbrot},
  39. //ARRAY PAR
  40. .nbuffers = 3,
  41. .modes = {STARPU_RW, STARPU_R, STARPU_R}
  42. //
  43. //STRUCT PAR
  44. /* .nbuffers = 1, */
  45. /* .modes = {STARPU_RW} */
  46. //
  47. };
  48. void mandelbrot_with_starpu(int64_t *pixels, float cr, float ci, unsigned width, unsigned height, unsigned nslicesx, unsigned nslicesy, double *params)
  49. {
  50. //ARRAY PAR
  51. starpu_data_handle_t p_handle, par_handle, v_handle;
  52. //
  53. //STRUCT PAR
  54. /* starpu_data_handle_t p_handle; */
  55. //
  56. starpu_matrix_data_register(&p_handle, STARPU_MAIN_RAM, (uintptr_t)pixels, width, width, height, sizeof(int64_t));
  57. //ARRAY PAR
  58. starpu_matrix_data_register(&par_handle, STARPU_MAIN_RAM, (uintptr_t)params, 5, 5, 1, sizeof(double));
  59. //
  60. struct starpu_data_filter vert =
  61. {
  62. /* .filter_func = starpu_matrix_filter_block, */
  63. .filter_func = starpu_matrix_filter_vertical_block,
  64. .nchildren = nslicesy
  65. };
  66. struct starpu_data_filter horiz =
  67. {
  68. .filter_func = starpu_matrix_filter_block,
  69. /* .filter_func = starpu_matrix_filter_vertical_block, */
  70. .nchildren = nslicesx
  71. };
  72. starpu_data_map_filters(p_handle, 2, &vert, &horiz);
  73. unsigned taskx, tasky;
  74. //ARRAY PAR
  75. int64_t *V = malloc(2 * nslicesx * nslicesy * sizeof(int64_t));
  76. //
  77. //STRUCT PAR
  78. /* struct Params *parameters = malloc(nslicesx * nslicesy * sizeof(struct Params)); */
  79. //
  80. for (tasky = 0; tasky < nslicesy; tasky++){
  81. for (taskx = 0; taskx < nslicesx; taskx++){
  82. struct starpu_task *task = starpu_task_create();
  83. //ARRAY PAR
  84. V[2 * (taskx + nslicesx * tasky)] = taskx + 1;
  85. V[2 * (taskx + nslicesx * tasky) + 1] = tasky + 1;
  86. starpu_vector_data_register(&v_handle, STARPU_MAIN_RAM, (uintptr_t)&(V[2 * (taskx + nslicesx * tasky)]), 2, sizeof(int64_t));
  87. //
  88. /* printf("Pre-Task%u_%u\n", taskx, tasky); */
  89. task->cl = &cl;
  90. task->handles[0] = starpu_data_get_sub_data(p_handle, 2, tasky, taskx);
  91. //ARRAY PAR
  92. task->handles[1] = par_handle;
  93. task->handles[2] = v_handle;
  94. //
  95. //STRUCT PAR
  96. /* struct Params param = {cr, ci, taskx, tasky, width, height}; */
  97. /* parameters[taskx + tasky * nslicesx] = param; */
  98. /* task->cl_arg = (parameters + taskx + tasky * nslicesx); */
  99. /* task->cl_arg_size = sizeof(struct Params); */
  100. //
  101. starpu_task_submit(task);
  102. //ARRAY PAR
  103. starpu_data_unregister_submit(v_handle);
  104. //
  105. }
  106. }
  107. starpu_task_wait_for_all();
  108. starpu_data_unpartition(p_handle, STARPU_MAIN_RAM);
  109. starpu_data_unregister(p_handle);
  110. //STRUCT PAR
  111. /* free(parameters); */
  112. //
  113. //ARRAY PAR
  114. starpu_data_unregister(par_handle);
  115. /* starpu_data_unregister(v_handle); */
  116. free(V);
  117. //
  118. }
  119. void init_zero(int64_t * pixels, unsigned width, unsigned height)
  120. {
  121. unsigned i,j;
  122. for (i = 0; i < height; i++){
  123. for (j = 0; j < width; j++){
  124. pixels[j + i*width] = 0;
  125. }
  126. }
  127. }
  128. void sort(double *arr, unsigned nbr_tests)
  129. {
  130. unsigned j;
  131. int is_sort = 0;
  132. while (!is_sort){
  133. is_sort = 1;
  134. for (j = 0; j < nbr_tests - 1; j++){
  135. if (arr[j] > arr[j+1]){
  136. is_sort = 0;
  137. double tmp = arr[j];
  138. arr[j] = arr[j+1];
  139. arr[j+1] = tmp;
  140. }
  141. }
  142. }
  143. }
  144. double median_time(float cr, float ci, unsigned width, unsigned height, unsigned nslicesx, unsigned nslicesy, unsigned nbr_tests)
  145. {
  146. int64_t *Pixels = malloc(width*height*sizeof(int64_t));
  147. double *params = malloc(4*sizeof(double));
  148. double max_iterations = (width/2) * 0.049715909 * log10(width * 0.25296875);
  149. params[0] = (double) cr;
  150. params[1] = (double) ci;
  151. params[2] = (double) width;
  152. params[3] = (double) height;
  153. params[4] = (double) max_iterations;
  154. unsigned i;
  155. double exec_times[nbr_tests];
  156. double start, stop, exec_t;
  157. for (i = 0; i < nbr_tests; i++){
  158. init_zero(Pixels, width, height);
  159. start = starpu_timing_now(); // starpu_timing_now() gives the time in microseconds.
  160. mandelbrot_with_starpu(Pixels, cr, ci, width, height, nslicesx, nslicesy, params);
  161. stop = starpu_timing_now();
  162. exec_t = (stop-start)/1.e6;
  163. exec_times[i] = exec_t;
  164. }
  165. char filename[34];
  166. sprintf(filename, "PPM/mandelbrottest%d.ppm", width);
  167. printf("%s\n", filename);
  168. /* Due to Julia registering matrices differently in memory, we need to transpose the matrix we get from the Julia generated kernels */
  169. mandelbrot_graph_transpose(filename, Pixels, width, height);
  170. /* mandelbrot_graph_transpose("PPM/mandelbrottest.ppm", Pixels, width, height); */
  171. free(Pixels);
  172. free(params);
  173. sort(exec_times, nbr_tests);
  174. return exec_times[nbr_tests/2];
  175. }
  176. void display_times(float cr, float ci, unsigned start_dim, unsigned step_dim, unsigned stop_dim, unsigned nslices, unsigned nbr_tests)
  177. {
  178. unsigned dim;
  179. FILE *myfile;
  180. myfile = fopen("DAT/mandelbrot_c_array_times.dat", "w");
  181. for (dim = start_dim; dim <= stop_dim; dim += step_dim){
  182. double t = median_time(cr, ci, dim, dim, nslices, nslices, nbr_tests);
  183. printf("w = %u ; h = %u ; t = %f\n", dim, dim, t);
  184. fprintf(myfile, "%f\n", t);
  185. }
  186. fclose(myfile);
  187. }
  188. int main(int argc, char **argv)
  189. {
  190. if (argc != 8){
  191. printf("Usage: %s cr ci start_dim step_dim stop_dim nslices(must divide dims) nbr_tests\n", argv[0]);
  192. return 1;
  193. }
  194. if (starpu_init(NULL) != EXIT_SUCCESS){
  195. fprintf(stderr, "ERROR\n");
  196. return 77;
  197. }
  198. float cr = (float) atof(argv[1]);
  199. float ci = (float) atof(argv[2]);
  200. unsigned start_dim = (unsigned) atoi(argv[3]);
  201. unsigned step_dim = (unsigned) atoi(argv[4]);
  202. unsigned stop_dim = (unsigned) atoi(argv[5]);
  203. unsigned nslices = (unsigned) atoi(argv[6]);
  204. unsigned nbr_tests = (unsigned) atoi(argv[7]);
  205. display_times(cr, ci, start_dim, step_dim, stop_dim, nslices, nbr_tests);
  206. starpu_shutdown();
  207. return 0;
  208. }