black_scholes_with_generated.c 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2019 Mael Keryell
  4. *
  5. * StarPU is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU Lesser General Public License as published by
  7. * the Free Software Foundation; either version 2.1 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * StarPU is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  13. *
  14. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  15. */
  16. #include <stdio.h>
  17. #include <stdlib.h>
  18. #include <starpu.h>
  19. #include "../includes/sorting.h"
  20. void black_scholes(void **, void *);
  21. void CUDA_black_scholes(void **, void*);
  22. static struct starpu_perfmodel model =
  23. {
  24. .type = STARPU_HISTORY_BASED,
  25. .symbol = "history_perf"
  26. };
  27. struct starpu_codelet cl =
  28. {
  29. .cpu_funcs = {black_scholes},
  30. .cuda_funcs = {CUDA_black_scholes},
  31. .nbuffers = 2,
  32. .modes = {STARPU_R, STARPU_RW},
  33. .model = &model
  34. };
  35. void black_scholes_with_starpu(double *data, double *res, unsigned nslices, unsigned nbr_data)
  36. {
  37. starpu_data_handle_t D_handle, RES_handle;
  38. starpu_matrix_data_register(&D_handle, STARPU_MAIN_RAM, (uintptr_t)data, 5, 5, nbr_data, sizeof(double));
  39. starpu_matrix_data_register(&RES_handle, STARPU_MAIN_RAM, (uintptr_t)res, 2, 2, nbr_data, sizeof(double));
  40. struct starpu_data_filter vert =
  41. {
  42. .filter_func = starpu_matrix_filter_vertical_block,
  43. .nchildren = nslices
  44. };
  45. starpu_data_partition(D_handle, &vert);
  46. starpu_data_partition(RES_handle, &vert);
  47. unsigned taskx;
  48. for (taskx = 0; taskx < nslices; taskx++){
  49. struct starpu_task *task = starpu_task_create();
  50. task->cl = &cl;
  51. task->handles[0] = starpu_data_get_sub_data(D_handle, 1, taskx);
  52. task->handles[1] = starpu_data_get_sub_data(RES_handle, 1, taskx);
  53. starpu_task_submit(task);
  54. }
  55. starpu_task_wait_for_all();
  56. starpu_data_unpartition(D_handle, STARPU_MAIN_RAM);
  57. starpu_data_unpartition(RES_handle, STARPU_MAIN_RAM);
  58. starpu_data_unregister(D_handle);
  59. starpu_data_unregister(RES_handle);
  60. }
  61. void init_data(double *data, unsigned nbr_data)
  62. {
  63. unsigned i;
  64. for (i = 0; i < nbr_data; i++){
  65. data[5*i] = 100. * rand() / (double) RAND_MAX;
  66. data[5*i + 1] = 100. * rand() / (double) RAND_MAX;
  67. data[5*i + 2] = rand() / (double) RAND_MAX;
  68. data[5*i + 3] = 10. * rand() / (double) RAND_MAX;
  69. data[5*i + 4] = 10. * rand() / (double) RAND_MAX;
  70. }
  71. }
  72. double median_time(unsigned nbr_data, unsigned nslices, unsigned nbr_tests)
  73. {
  74. double *data = malloc(5 * nbr_data * sizeof(double));
  75. double *res = calloc(2 * nbr_data, sizeof(double));
  76. double exec_times[nbr_tests];
  77. /* printf("nbr_data: %u\n", nbr_data); */
  78. unsigned i;
  79. for (i = 0; i < nbr_tests; i++){
  80. init_data(data, nbr_data);
  81. /* data[0] = 100.0; */
  82. /* data[1] = 100.0; */
  83. /* data[2] = 0.05; */
  84. /* data[3] = 1.0; */
  85. /* data[4] = 0.2; */
  86. double start = starpu_timing_now();
  87. black_scholes_with_starpu(data, res, nslices, nbr_data);
  88. double stop = starpu_timing_now();
  89. exec_times[i] = (stop-start)/1.e6;
  90. }
  91. /* printf("RES:\n%f\n%f\n", res[0], res[1]); */
  92. free(data);
  93. free(res);
  94. quicksort(exec_times, 0, nbr_tests - 1);
  95. return exec_times[nbr_tests/2];
  96. }
  97. void display_times(unsigned start_nbr, unsigned step_nbr, unsigned stop_nbr, unsigned nslices, unsigned nbr_tests){
  98. double t;
  99. unsigned nbr_data;
  100. FILE *myfile;
  101. myfile = fopen("DAT/black_scholes_c_generated_times.dat", "w");
  102. for (nbr_data = start_nbr; nbr_data <= stop_nbr; nbr_data+=step_nbr){
  103. t = median_time(nbr_data, nslices, nbr_tests);
  104. printf("Number of data: %u\nTime: %f\n", nbr_data, t);
  105. fprintf(myfile, "%f\n", t);
  106. }
  107. fclose(myfile);
  108. }
  109. int main(int argc, char *argv[])
  110. {
  111. if (argc != 6){
  112. printf("Usage: %s start_nbr step_nbr stop_nbr nslices nbr_tests\n", argv[0]);
  113. return 1;
  114. }
  115. if (starpu_init(NULL) != EXIT_SUCCESS){
  116. fprintf(stderr, "ERROR\n");
  117. return 77;
  118. }
  119. unsigned start_nbr = (unsigned) atoi(argv[1]);
  120. unsigned step_nbr = (unsigned) atoi(argv[2]);
  121. unsigned stop_nbr = (unsigned) atoi(argv[3]);
  122. unsigned nslices = (unsigned) atoi(argv[4]);
  123. unsigned nbr_tests = (unsigned) atoi(argv[5]);
  124. display_times(start_nbr, step_nbr, stop_nbr, nslices, nbr_tests);
  125. starpu_shutdown();
  126. return 0;
  127. }