16mpi_support.doxy 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540
  1. /*
  2. * This file is part of the StarPU Handbook.
  3. * Copyright (C) 2009--2011 Universit@'e de Bordeaux
  4. * Copyright (C) 2010, 2011, 2012, 2013, 2014 Centre National de la Recherche Scientifique
  5. * Copyright (C) 2011, 2012 Institut National de Recherche en Informatique et Automatique
  6. * See the file version.doxy for copying conditions.
  7. */
  8. /*! \page MPISupport MPI Support
  9. The integration of MPI transfers within task parallelism is done in a
  10. very natural way by the means of asynchronous interactions between the
  11. application and StarPU. This is implemented in a separate libstarpumpi library
  12. which basically provides "StarPU" equivalents of <c>MPI_*</c> functions, where
  13. <c>void *</c> buffers are replaced with ::starpu_data_handle_t, and all
  14. GPU-RAM-NIC transfers are handled efficiently by StarPU-MPI. The user has to
  15. use the usual <c>mpirun</c> command of the MPI implementation to start StarPU on
  16. the different MPI nodes.
  17. An MPI Insert Task function provides an even more seamless transition to a
  18. distributed application, by automatically issuing all required data transfers
  19. according to the task graph and an application-provided distribution.
  20. \section SimpleExample Simple Example
  21. The flags required to compile or link against the MPI layer are
  22. accessible with the following commands:
  23. \verbatim
  24. $ pkg-config --cflags starpumpi-1.2 # options for the compiler
  25. $ pkg-config --libs starpumpi-1.2 # options for the linker
  26. \endverbatim
  27. You also need pass the option <c>--static</c> if the application is to
  28. be linked statically.
  29. \code{.c}
  30. void increment_token(void)
  31. {
  32. struct starpu_task *task = starpu_task_create();
  33. task->cl = &increment_cl;
  34. task->handles[0] = token_handle;
  35. starpu_task_submit(task);
  36. }
  37. int main(int argc, char **argv)
  38. {
  39. int rank, size;
  40. starpu_init(NULL);
  41. starpu_mpi_initialize_extended(&rank, &size);
  42. starpu_vector_data_register(&token_handle, STARPU_MAIN_RAM, (uintptr_t)&token, 1, sizeof(unsigned));
  43. unsigned nloops = NITER;
  44. unsigned loop;
  45. unsigned last_loop = nloops - 1;
  46. unsigned last_rank = size - 1;
  47. for (loop = 0; loop < nloops; loop++) {
  48. int tag = loop*size + rank;
  49. if (loop == 0 && rank == 0)
  50. {
  51. token = 0;
  52. fprintf(stdout, "Start with token value %d\n", token);
  53. }
  54. else
  55. {
  56. starpu_mpi_irecv_detached(token_handle, (rank+size-1)%size, tag,
  57. MPI_COMM_WORLD, NULL, NULL);
  58. }
  59. increment_token();
  60. if (loop == last_loop && rank == last_rank)
  61. {
  62. starpu_data_acquire(token_handle, STARPU_R);
  63. fprintf(stdout, "Finished: token value %d\n", token);
  64. starpu_data_release(token_handle);
  65. }
  66. else
  67. {
  68. starpu_mpi_isend_detached(token_handle, (rank+1)%size, tag+1,
  69. MPI_COMM_WORLD, NULL, NULL);
  70. }
  71. }
  72. starpu_task_wait_for_all();
  73. starpu_mpi_shutdown();
  74. starpu_shutdown();
  75. if (rank == last_rank)
  76. {
  77. fprintf(stderr, "[%d] token = %d == %d * %d ?\n", rank, token, nloops, size);
  78. STARPU_ASSERT(token == nloops*size);
  79. }
  80. \endcode
  81. \section PointToPointCommunication Point To Point Communication
  82. The standard point to point communications of MPI have been
  83. implemented. The semantic is similar to the MPI one, but adapted to
  84. the DSM provided by StarPU. A MPI request will only be submitted when
  85. the data is available in the main memory of the node submitting the
  86. request.
  87. There is two types of asynchronous communications: the classic
  88. asynchronous communications and the detached communications. The
  89. classic asynchronous communications (starpu_mpi_isend() and
  90. starpu_mpi_irecv()) need to be followed by a call to
  91. starpu_mpi_wait() or to starpu_mpi_test() to wait for or to
  92. test the completion of the communication. Waiting for or testing the
  93. completion of detached communications is not possible, this is done
  94. internally by StarPU-MPI, on completion, the resources are
  95. automatically released. This mechanism is similar to the pthread
  96. detach state attribute which determines whether a thread will be
  97. created in a joinable or a detached state.
  98. Internally, all communication are divided in 2 communications, a first
  99. message is used to exchange an envelope describing the data (i.e its
  100. tag and its size), the data itself is sent in a second message. All
  101. MPI communications submitted by StarPU uses a unique tag which has a
  102. default value, and can be accessed with the functions
  103. starpu_mpi_get_communication_tag() and
  104. starpu_mpi_set_communication_tag(). The matching of tags with
  105. corresponding requests is done within StarPU-MPI.
  106. For any userland communication, the call of the corresponding function
  107. (e.g starpu_mpi_isend()) will result in the creation of a StarPU-MPI
  108. request, the function starpu_data_acquire_cb() is then called to
  109. asynchronously request StarPU to fetch the data in main memory; when
  110. the data is ready and the corresponding buffer has already been
  111. received by MPI, it will be copied in the memory of the data,
  112. otherwise the request is stored in the <em>early requests list</em>. Sending
  113. requests are stored in the <em>ready requests list</em>.
  114. While requests need to be processed, the StarPU-MPI progression thread
  115. does the following:
  116. <ol>
  117. <li> it polls the <em>ready requests list</em>. For all the ready
  118. requests, the appropriate function is called to post the corresponding
  119. MPI call. For example, an initial call to starpu_mpi_isend() will
  120. result in a call to <c>MPI_Isend</c>. If the request is marked as
  121. detached, the request will then be added in the <em>detached requests
  122. list</em>.
  123. </li>
  124. <li> it posts a <c>MPI_Irecv()</c> to retrieve a data envelope.
  125. </li>
  126. <li> it polls the <em>detached requests list</em>. For all the detached
  127. requests, it tests its completion of the MPI request by calling
  128. <c>MPI_Test</c>. On completion, the data handle is released, and if a
  129. callback was defined, it is called.
  130. </li>
  131. <li> finally, it checks if a data envelope has been received. If so,
  132. if the data envelope matches a request in the <em>early requests list</em> (i.e
  133. the request has already been posted by the application), the
  134. corresponding MPI call is posted (similarly to the first step above).
  135. If the data envelope does not match any application request, a
  136. temporary handle is created to receive the data, a StarPU-MPI request
  137. is created and added into the <em>ready requests list</em>, and thus will be
  138. processed in the first step of the next loop.
  139. </li>
  140. </ol>
  141. \ref MPIPtpCommunication "Communication" gives the list of all the
  142. point to point communications defined in StarPU-MPI.
  143. \section ExchangingUserDefinedDataInterface Exchanging User Defined Data Interface
  144. New data interfaces defined as explained in \ref
  145. DefiningANewDataInterface can also be used within StarPU-MPI and
  146. exchanged between nodes. Two functions needs to be defined through the
  147. type starpu_data_interface_ops. The function
  148. starpu_data_interface_ops::pack_data takes a handle and returns a
  149. contiguous memory buffer along with its size where data to be conveyed
  150. to another node should be copied. The reversed operation is
  151. implemented in the function starpu_data_interface_ops::unpack_data which
  152. takes a contiguous memory buffer and recreates the data handle.
  153. \code{.c}
  154. static int complex_pack_data(starpu_data_handle_t handle, unsigned node, void **ptr, ssize_t *count)
  155. {
  156. STARPU_ASSERT(starpu_data_test_if_allocated_on_node(handle, node));
  157. struct starpu_complex_interface *complex_interface =
  158. (struct starpu_complex_interface *) starpu_data_get_interface_on_node(handle, node);
  159. *count = complex_get_size(handle);
  160. *ptr = malloc(*count);
  161. memcpy(*ptr, complex_interface->real, complex_interface->nx*sizeof(double));
  162. memcpy(*ptr+complex_interface->nx*sizeof(double), complex_interface->imaginary,
  163. complex_interface->nx*sizeof(double));
  164. return 0;
  165. }
  166. static int complex_unpack_data(starpu_data_handle_t handle, unsigned node, void *ptr, size_t count)
  167. {
  168. STARPU_ASSERT(starpu_data_test_if_allocated_on_node(handle, node));
  169. struct starpu_complex_interface *complex_interface =
  170. (struct starpu_complex_interface *) starpu_data_get_interface_on_node(handle, node);
  171. memcpy(complex_interface->real, ptr, complex_interface->nx*sizeof(double));
  172. memcpy(complex_interface->imaginary, ptr+complex_interface->nx*sizeof(double),
  173. complex_interface->nx*sizeof(double));
  174. return 0;
  175. }
  176. static struct starpu_data_interface_ops interface_complex_ops =
  177. {
  178. ...
  179. .pack_data = complex_pack_data,
  180. .unpack_data = complex_unpack_data
  181. };
  182. \endcode
  183. \section MPIInsertTaskUtility MPI Insert Task Utility
  184. To save the programmer from having to explicit all communications, StarPU
  185. provides an "MPI Insert Task Utility". The principe is that the application
  186. decides a distribution of the data over the MPI nodes by allocating it and
  187. notifying StarPU of that decision, i.e. tell StarPU which MPI node "owns"
  188. which data. It also decides, for each handle, an MPI tag which will be used to
  189. exchange the content of the handle. All MPI nodes then process the whole task
  190. graph, and StarPU automatically determines which node actually execute which
  191. task, and trigger the required MPI transfers.
  192. The list of functions is described in \ref MPIInsertTask "MPI Insert Task".
  193. Here an stencil example showing how to use starpu_mpi_task_insert(). One
  194. first needs to define a distribution function which specifies the
  195. locality of the data. Note that the data needs to be registered to MPI
  196. by calling starpu_mpi_data_register(). This function allows to set
  197. the distribution information and the MPI tag which should be used when
  198. communicating the data. The function starpu_mpi_data_register() should
  199. be prefered to starpu_data_set_rank() and starpu_data_set_tag() as
  200. it also allows to automatically clear the MPI communication cache
  201. when unregistering the data.
  202. \code{.c}
  203. /* Returns the MPI node number where data is */
  204. int my_distrib(int x, int y, int nb_nodes) {
  205. /* Block distrib */
  206. return ((int)(x / sqrt(nb_nodes) + (y / sqrt(nb_nodes)) * sqrt(nb_nodes))) % nb_nodes;
  207. // /* Other examples useful for other kinds of computations */
  208. // /* / distrib */
  209. // return (x+y) % nb_nodes;
  210. // /* Block cyclic distrib */
  211. // unsigned side = sqrt(nb_nodes);
  212. // return x % side + (y % side) * size;
  213. }
  214. \endcode
  215. Now the data can be registered within StarPU. Data which are not
  216. owned but will be needed for computations can be registered through
  217. the lazy allocation mechanism, i.e. with a <c>home_node</c> set to <c>-1</c>.
  218. StarPU will automatically allocate the memory when it is used for the
  219. first time.
  220. One can note an optimization here (the <c>else if</c> test): we only register
  221. data which will be needed by the tasks that we will execute.
  222. \code{.c}
  223. unsigned matrix[X][Y];
  224. starpu_data_handle_t data_handles[X][Y];
  225. for(x = 0; x < X; x++) {
  226. for (y = 0; y < Y; y++) {
  227. int mpi_rank = my_distrib(x, y, size);
  228. if (mpi_rank == my_rank)
  229. /* Owning data */
  230. starpu_variable_data_register(&data_handles[x][y], STARPU_MAIN_RAM,
  231. (uintptr_t)&(matrix[x][y]), sizeof(unsigned));
  232. else if (my_rank == my_distrib(x+1, y, size) || my_rank == my_distrib(x-1, y, size)
  233. || my_rank == my_distrib(x, y+1, size) || my_rank == my_distrib(x, y-1, size))
  234. /* I don't own that index, but will need it for my computations */
  235. starpu_variable_data_register(&data_handles[x][y], -1,
  236. (uintptr_t)NULL, sizeof(unsigned));
  237. else
  238. /* I know it's useless to allocate anything for this */
  239. data_handles[x][y] = NULL;
  240. if (data_handles[x][y]) {
  241. starpu_mpi_data_register(data_handles[x][y], x*X+y, mpi_rank);
  242. }
  243. }
  244. }
  245. \endcode
  246. Now starpu_mpi_task_insert() can be called for the different
  247. steps of the application.
  248. \code{.c}
  249. for(loop=0 ; loop<niter; loop++)
  250. for (x = 1; x < X-1; x++)
  251. for (y = 1; y < Y-1; y++)
  252. starpu_mpi_task_insert(MPI_COMM_WORLD, &stencil5_cl,
  253. STARPU_RW, data_handles[x][y],
  254. STARPU_R, data_handles[x-1][y],
  255. STARPU_R, data_handles[x+1][y],
  256. STARPU_R, data_handles[x][y-1],
  257. STARPU_R, data_handles[x][y+1],
  258. 0);
  259. starpu_task_wait_for_all();
  260. \endcode
  261. I.e. all MPI nodes process the whole task graph, but as mentioned above, for
  262. each task, only the MPI node which owns the data being written to (here,
  263. <c>data_handles[x][y]</c>) will actually run the task. The other MPI nodes will
  264. automatically send the required data.
  265. This can be a concern with a growing number of nodes. To avoid this, the
  266. application can prune the task for loops according to the data distribution,
  267. so as to only submit tasks on nodes which have to care about them (either to
  268. execute them, or to send the required data).
  269. A way to do some of this quite easily can be to just add an <c>if</c> like this:
  270. \code{.c}
  271. for(loop=0 ; loop<niter; loop++)
  272. for (x = 1; x < X-1; x++)
  273. for (y = 1; y < Y-1; y++)
  274. if (my_distrib(x,y,size) == my_rank
  275. || my_distrib(x-1,y,size) == my_rank
  276. || my_distrib(x+1,y,size) == my_rank
  277. || my_distrib(x,y-1,size) == my_rank
  278. || my_distrib(x,y+1,size) == my_rank)
  279. starpu_mpi_insert_task(MPI_COMM_WORLD, &stencil5_cl,
  280. STARPU_RW, data_handles[x][y],
  281. STARPU_R, data_handles[x-1][y],
  282. STARPU_R, data_handles[x+1][y],
  283. STARPU_R, data_handles[x][y-1],
  284. STARPU_R, data_handles[x][y+1],
  285. 0);
  286. starpu_task_wait_for_all();
  287. \endcode
  288. This permits to drop the cost of function call argument passing and parsing.
  289. If the <c>my_distrib</c> function can be inlined by the compiler, the latter can
  290. improve the test.
  291. If the <c>size</c> can be made a compile-time constant, the compiler can
  292. considerably improve the test further.
  293. If the distribution function is not too complex and the compiler is very good,
  294. the latter can even optimize the <c>for</c> loops, thus dramatically reducing
  295. the cost of task submission.
  296. A function starpu_mpi_task_build() is also provided with the aim to
  297. only construct the task structure. All MPI nodes need to call the
  298. function, only the node which is to execute the task will return a
  299. valid task structure. Following the execution of the task, all nodes
  300. need to call the function starpu_mpi_task_post_build() -- with the same
  301. list of arguments as starpu_mpi_task_build() -- to post all the
  302. necessary data communications.
  303. \code{.c}
  304. struct starpu_task *task;
  305. task = starpu_mpi_task_build(MPI_COMM_WORLD, &cl,
  306. STARPU_RW, data_handles[0],
  307. STARPU_R, data_handles[1],
  308. 0);
  309. if (task) starpu_task_submit(task);
  310. starpu_mpi_task_post_build(MPI_COMM_WORLD, &cl,
  311. STARPU_RW, data_handles[0],
  312. STARPU_R, data_handles[1],
  313. 0);
  314. \endcode
  315. \section MPICache MPI cache support
  316. StarPU-MPI automatically optimizes duplicate data transmissions: if an MPI
  317. node B needs a piece of data D from MPI node A for several tasks, only one
  318. transmission of D will take place from A to B, and the value of D will be kept
  319. on B as long as no task modifies D.
  320. If a task modifies D, B will wait for all tasks which need the previous value of
  321. D, before invalidating the value of D. As a consequence, it releases the memory
  322. occupied by D. Whenever a task running on B needs the new value of D, allocation
  323. will take place again to receive it.
  324. Since tasks can be submitted dynamically, StarPU-MPI can not know whether the
  325. current value of data D will again be used by a newly-submitted task before
  326. being modified by another newly-submitted task, so until a task is submitted to
  327. modify the current value, it can not decide by itself whether to flush the cache
  328. or not. The application can however explicitly tell StarPU-MPI to flush the
  329. cache by calling starpu_mpi_cache_flush() or starpu_mpi_cache_flush_all_data(),
  330. for instance in case the data will not be used at all any more (see for instance
  331. the cholesky example in mpi/examples/matrix_decomposition), or at least not in
  332. the close future. If a newly-submitted task actually needs the value again,
  333. another transmission of D will be initiated from A to B.
  334. The whole caching behavior can be disabled thanks to the ::STARPU_MPI_CACHE
  335. environment variable. The variable ::STARPU_MPI_CACHE_STATS can be set to 1
  336. to enable the runtime to display messages when data are added or removed
  337. from the cache holding the received data.
  338. \section MPIMigration MPI Data migration
  339. The application can dynamically change its mind about the data distribution, to
  340. balance the load over MPI nodes for instance. This can be done very simply by
  341. requesting an explicit move and then change the registered rank. For instance,
  342. we here switch to a new distribution function <c>my_distrib2</c>: we first
  343. register any data that wasn't registered already and will be needed, then
  344. migrate the data, and register the new location.
  345. \code{.c}
  346. for(x = 0; x < X; x++) {
  347. for (y = 0; y < Y; y++) {
  348. int mpi_rank = my_distrib2(x, y, size);
  349. if (!data_handles[x][y] && (mpi_rank == my_rank
  350. || my_rank == my_distrib(x+1, y, size) || my_rank == my_distrib(x-1, y, size)
  351. || my_rank == my_distrib(x, y+1, size) || my_rank == my_distrib(x, y-1, size)))
  352. /* Register newly-needed data */
  353. starpu_variable_data_register(&data_handles[x][y], -1,
  354. (uintptr_t)NULL, sizeof(unsigned));
  355. if (data_handles[x][y]) {
  356. /* Migrate the data */
  357. starpu_mpi_get_data_on_node_detached(MPI_COMM_WORLD, data_handles[x][y], mpi_rank, NULL, NULL);
  358. /* And register the new rank of the matrix */
  359. starpu_data_set_rank(data_handles[x][y], mpi_rank);
  360. }
  361. }
  362. }
  363. \endcode
  364. From then on, further tasks submissions will use the new data distribution,
  365. which will thus change both MPI communications and task assignments.
  366. Very importantly, since all nodes have to agree on which node owns which data
  367. so as to determine MPI communications and task assignments the same way, all
  368. nodes have to perform the same data migration, and at the same point among task
  369. submissions. It thus does not require a strict synchronization, just a clear
  370. separation of task submissions before and after the data redistribution.
  371. Before data unregistration, it has to be migrated back to its original home
  372. node (the value, at least), since that is where the user-provided buffer
  373. resides. Otherwise the unregistration will complain that it does not have the
  374. latest value on the original home node.
  375. \code{.c}
  376. for(x = 0; x < X; x++) {
  377. for (y = 0; y < Y; y++) {
  378. if (data_handles[x][y]) {
  379. int mpi_rank = my_distrib(x, y, size);
  380. /* Get back data to original place where the user-provided buffer is. */
  381. starpu_mpi_get_data_on_node_detached(MPI_COMM_WORLD, data_handles[x][y], mpi_rank, NULL, NULL);
  382. /* And unregister it */
  383. starpu_data_unregister(data_handles[x][y]);
  384. }
  385. }
  386. }
  387. \endcode
  388. \section MPICollective MPI Collective Operations
  389. The functions are described in \ref MPICollectiveOperations "MPI Collective Operations".
  390. \code{.c}
  391. if (rank == root)
  392. {
  393. /* Allocate the vector */
  394. vector = malloc(nblocks * sizeof(float *));
  395. for(x=0 ; x<nblocks ; x++)
  396. {
  397. starpu_malloc((void **)&vector[x], block_size*sizeof(float));
  398. }
  399. }
  400. /* Allocate data handles and register data to StarPU */
  401. data_handles = malloc(nblocks*sizeof(starpu_data_handle_t *));
  402. for(x = 0; x < nblocks ; x++)
  403. {
  404. int mpi_rank = my_distrib(x, nodes);
  405. if (rank == root) {
  406. starpu_vector_data_register(&data_handles[x], STARPU_MAIN_RAM, (uintptr_t)vector[x],
  407. blocks_size, sizeof(float));
  408. }
  409. else if ((mpi_rank == rank) || ((rank == mpi_rank+1 || rank == mpi_rank-1))) {
  410. /* I own that index, or i will need it for my computations */
  411. starpu_vector_data_register(&data_handles[x], -1, (uintptr_t)NULL,
  412. block_size, sizeof(float));
  413. }
  414. else {
  415. /* I know it's useless to allocate anything for this */
  416. data_handles[x] = NULL;
  417. }
  418. if (data_handles[x]) {
  419. starpu_mpi_data_register(data_handles[x], x*nblocks+y, mpi_rank);
  420. }
  421. }
  422. /* Scatter the matrix among the nodes */
  423. starpu_mpi_scatter_detached(data_handles, nblocks, root, MPI_COMM_WORLD);
  424. /* Calculation */
  425. for(x = 0; x < nblocks ; x++) {
  426. if (data_handles[x]) {
  427. int owner = starpu_data_get_rank(data_handles[x]);
  428. if (owner == rank) {
  429. starpu_task_insert(&cl, STARPU_RW, data_handles[x], 0);
  430. }
  431. }
  432. }
  433. /* Gather the matrix on main node */
  434. starpu_mpi_gather_detached(data_handles, nblocks, 0, MPI_COMM_WORLD);
  435. \endcode
  436. */
  437. \section MPIExamples More MPI examples
  438. MPI examples are available in the StarPU source code in mpi/examples:
  439. <ul>
  440. <li><c>complex</c> is a simple example using a user-define data interface over
  441. MPI (complex numbers),
  442. <li><c>stencil5</c> is a simple stencil example using starpu_mpi_task_insert(),
  443. <li><c>matrix_decomposition</c> is a cholesky decomposition example using
  444. starpu_mpi_task_insert(). The non-distributed version can check for
  445. <algorithm correctness in 1-node configuration, the distributed version uses
  446. exactly the same source code, to be used over MPI,
  447. <li><c>mpi_lu</c> is an LU decomposition example, provided in three versions:
  448. <c>plu_example</c> uses explicit MPI data transfers, <c>plu_implicit_example</c>
  449. uses implicit MPI data transfers, <c>plu_outofcore_example</c> uses implicit MPI
  450. data transfers and supports data matrices which do not fit in memory (out-of-core).
  451. </ul>