vector_scale.c 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290
  1. /* StarPURM --- StarPU Resource Management Layer.
  2. *
  3. * Copyright (C) 2017 Inria
  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. /* This example shows a basic StarPU vector scale app on top of StarPURM with a nVidia CUDA kernel */
  17. #include <stdlib.h>
  18. #include <stdio.h>
  19. #include <assert.h>
  20. #include <starpu.h>
  21. #include <starpurm.h>
  22. static int rm_cpu_type_id = -1;
  23. static int rm_cuda_type_id = -1;
  24. static int rm_nb_cpu_units = 0;
  25. static int rm_nb_cuda_units = 0;
  26. static void usage(void);
  27. static void test1(const int N);
  28. static void test2(const int N, const int task_mult);
  29. static void init_rm_infos(void);
  30. /* vector scale codelet */
  31. static void vector_scale_func(void *cl_buffers[], void *cl_arg)
  32. {
  33. float scalar = -1.0;
  34. int n = STARPU_VECTOR_GET_NX(cl_buffers[0]);
  35. float *vector = (float *)STARPU_VECTOR_GET_PTR(cl_buffers[0]);
  36. int i;
  37. starpu_codelet_unpack_args(cl_arg, &scalar);
  38. {
  39. int workerid = starpu_worker_get_id();
  40. hwloc_cpuset_t worker_cpuset = starpu_worker_get_hwloc_cpuset(workerid);
  41. hwloc_cpuset_t check_cpuset = starpurm_get_selected_cpuset();
  42. #if 0
  43. {
  44. int strl1 = hwloc_bitmap_snprintf(NULL, 0, worker_cpuset);
  45. char str1[strl1+1];
  46. hwloc_bitmap_snprintf(str1, strl1+1, worker_cpuset);
  47. int strl2 = hwloc_bitmap_snprintf(NULL, 0, check_cpuset);
  48. char str2[strl2+1];
  49. hwloc_bitmap_snprintf(str2, strl2+1, check_cpuset);
  50. printf("worker[%03d] - task: vector=%p, n=%d, scalar=%lf, worker cpuset = %s, selected cpuset = %s\n", workerid, vector, n, scalar, str1, str2);
  51. }
  52. #endif
  53. hwloc_bitmap_and(check_cpuset, check_cpuset, worker_cpuset);
  54. assert(!hwloc_bitmap_iszero(check_cpuset));
  55. hwloc_bitmap_free(check_cpuset);
  56. hwloc_bitmap_free(worker_cpuset);
  57. }
  58. for (i = 0; i < n; i++)
  59. {
  60. vector[i] *= scalar;
  61. }
  62. }
  63. extern void vector_scale_cuda_func(void *cl_buffers[], void *cl_arg);
  64. static struct starpu_codelet vector_scale_cl =
  65. {
  66. .cpu_funcs = {vector_scale_func},
  67. .cuda_funcs = {vector_scale_cuda_func},
  68. .cuda_flags = {STARPU_CUDA_ASYNC},
  69. .nbuffers = 1
  70. };
  71. /* main routines */
  72. static void usage(void)
  73. {
  74. fprintf(stderr, "usage: 05_vector_scale [VECTOR_SIZE]\n");
  75. exit(1);
  76. }
  77. static void test1(const int N)
  78. {
  79. float *vector = NULL;
  80. const float scalar = 2.0;
  81. starpu_data_handle_t vector_handle;
  82. int ret;
  83. vector = malloc(N * sizeof(*vector));
  84. {
  85. int i;
  86. for (i = 0; i < N; i++)
  87. {
  88. vector[i] = i;
  89. }
  90. }
  91. starpu_vector_data_register(&vector_handle, STARPU_MAIN_RAM, (uintptr_t)vector, N, sizeof(*vector));
  92. ret = starpu_task_insert(&vector_scale_cl,
  93. STARPU_RW, vector_handle,
  94. STARPU_VALUE, &scalar, sizeof(scalar),
  95. 0);
  96. assert(ret == 0);
  97. starpu_task_wait_for_all();
  98. starpu_data_unregister(vector_handle);
  99. {
  100. int i;
  101. for (i = 0; i < N; i++)
  102. {
  103. float d_i = i;
  104. if (vector[i] != d_i*scalar)
  105. {
  106. fprintf(stderr, "%s: check_failed, vector[%d]: %f != %f\n", __func__, i, vector[i], d_i*scalar);
  107. exit(1);
  108. }
  109. }
  110. }
  111. free(vector);
  112. }
  113. static void test2(const int N, const int task_mult)
  114. {
  115. float *vector = NULL;
  116. const float scalar = 3.0;
  117. starpu_data_handle_t vector_handle;
  118. int ret;
  119. vector = malloc(N * sizeof(*vector));
  120. {
  121. int i;
  122. for (i = 0; i < N; i++)
  123. {
  124. vector[i] = i;
  125. }
  126. }
  127. const int nparts = (rm_nb_cpu_units+rm_nb_cuda_units) * task_mult;
  128. starpu_vector_data_register(&vector_handle, STARPU_MAIN_RAM, (uintptr_t)vector, N, sizeof(*vector));
  129. struct starpu_data_filter partition_filter =
  130. {
  131. .filter_func = starpu_vector_filter_block,
  132. .nchildren = nparts
  133. };
  134. starpu_data_partition(vector_handle, &partition_filter);
  135. {
  136. int i;
  137. for (i = 0; i < nparts; i++)
  138. {
  139. starpu_data_handle_t sub_vector_handle = starpu_data_get_sub_data(vector_handle, 1, i);
  140. ret = starpu_task_insert(&vector_scale_cl,
  141. STARPU_RW, sub_vector_handle,
  142. STARPU_VALUE, &scalar, sizeof(scalar),
  143. 0);
  144. assert(ret == 0);
  145. }
  146. }
  147. starpu_task_wait_for_all();
  148. starpu_data_unpartition(vector_handle, STARPU_MAIN_RAM);
  149. starpu_data_unregister(vector_handle);
  150. {
  151. int i;
  152. for (i = 0; i < N; i++)
  153. {
  154. float d_i = i;
  155. if (vector[i] != d_i*scalar)
  156. {
  157. fprintf(stderr, "%s: check_failed, vector[%d]: %f != %f\n", __func__, i, vector[i], d_i*scalar);
  158. exit(1);
  159. }
  160. }
  161. }
  162. free(vector);
  163. }
  164. static void init_rm_infos(void)
  165. {
  166. int cpu_type = starpurm_get_device_type_id("cpu");
  167. int nb_cpu_units = starpurm_get_nb_devices_by_type(cpu_type);
  168. if (nb_cpu_units < 1)
  169. {
  170. /* No CPU unit available. */
  171. exit(77);
  172. }
  173. int cuda_type = starpurm_get_device_type_id("cuda");
  174. int nb_cuda_units = starpurm_get_nb_devices_by_type(cuda_type);
  175. if (nb_cuda_units < 1)
  176. {
  177. /* No CUDA unit available. */
  178. exit(77);
  179. }
  180. rm_cpu_type_id = cpu_type;
  181. rm_cuda_type_id = cuda_type;
  182. rm_nb_cpu_units = nb_cpu_units;
  183. rm_nb_cuda_units = nb_cuda_units;
  184. }
  185. static void disp_selected_cpuset(void)
  186. {
  187. hwloc_cpuset_t selected_cpuset = starpurm_get_selected_cpuset();
  188. int strl = hwloc_bitmap_snprintf(NULL, 0, selected_cpuset);
  189. char str[strl+1];
  190. hwloc_bitmap_snprintf(str, strl+1, selected_cpuset);
  191. printf("selected cpuset = %s\n", str);
  192. }
  193. int main(int argc, char *argv[])
  194. {
  195. int param_N = 1000000;
  196. int drs_enabled;
  197. if (argc > 1)
  198. {
  199. param_N = atoi(argv[1]);
  200. if (param_N < 1)
  201. {
  202. usage();
  203. }
  204. }
  205. starpurm_initialize();
  206. init_rm_infos();
  207. printf("using default units\n");
  208. disp_selected_cpuset();
  209. test1(param_N);
  210. test2(param_N, 1);
  211. test2(param_N, 10);
  212. test2(param_N, 100);
  213. if (rm_nb_cpu_units > 1 && rm_nb_cuda_units > 1)
  214. {
  215. int nb_cpus = rm_nb_cpu_units;
  216. const int nb_cudas = rm_nb_cuda_units;
  217. const int cuda_type = rm_cuda_type_id;
  218. printf("nb_cpu_units = %d\n", nb_cpus);
  219. printf("nb_cuda_units = %d\n", nb_cudas);
  220. /* Keep at least one CPU core */
  221. nb_cpus--;
  222. starpurm_set_drs_enable(NULL);
  223. drs_enabled = starpurm_drs_enabled_p();
  224. assert(drs_enabled != 0);
  225. printf("withdrawing %d cpus from StarPU\n", nb_cpus);
  226. starpurm_withdraw_cpus_from_starpu(NULL, nb_cpus);
  227. disp_selected_cpuset();
  228. test2(param_N, 1);
  229. test2(param_N, 10);
  230. test2(param_N, 100);
  231. printf("assigning %d cpus to StarPU\n", nb_cpus);
  232. starpurm_assign_cpus_to_starpu(NULL, nb_cpus);
  233. disp_selected_cpuset();
  234. test2(param_N, 1);
  235. test2(param_N, 10);
  236. test2(param_N, 100);
  237. printf("withdrawing %d cuda devices from StarPU\n", nb_cudas);
  238. starpurm_withdraw_devices_from_starpu(NULL, cuda_type, nb_cudas);
  239. disp_selected_cpuset();
  240. test2(param_N, 1);
  241. test2(param_N, 10);
  242. test2(param_N, 100);
  243. printf("lending %d cuda devices to StarPU\n", nb_cudas);
  244. starpurm_assign_devices_to_starpu(NULL, cuda_type, nb_cudas);
  245. disp_selected_cpuset();
  246. test2(param_N, 1);
  247. test2(param_N, 10);
  248. test2(param_N, 100);
  249. starpurm_set_drs_disable(NULL);
  250. drs_enabled = starpurm_drs_enabled_p();
  251. assert(drs_enabled == 0);
  252. }
  253. starpurm_shutdown();
  254. return 0;
  255. }