1 /* Copyright (c) 2013-2022. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
7 #include "../colls_private.hpp"
8 /*****************************************************************************
10 * Function: alltoall_ring_light_barrier
15 send_buff: send input buffer
16 send_count: number of elements to send
17 send_type: data type of elements being sent
18 recv_buff: receive output buffer
19 recv_count: number of elements to received
20 recv_type: data type of elements being received
23 * Descrp: Function works in P - 1 steps. In step i, node j - i -> j -> j + i.
24 Light barriers are inserted between communications in different
29 ****************************************************************************/
33 alltoall__ring_light_barrier(const void *send_buff, int send_count,
34 MPI_Datatype send_type,
35 void *recv_buff, int recv_count,
36 MPI_Datatype recv_type,
39 MPI_Aint send_chunk, recv_chunk;
41 int i, src, dst, rank, num_procs, next_dst, next_src;
42 int tag = COLL_TAG_ALLTOALL;
44 char send_sync = 'a', recv_sync = 'b';
45 char *send_ptr = (char *) send_buff;
46 char *recv_ptr = (char *) recv_buff;
49 num_procs = comm->size();
50 send_chunk = send_type->get_extent();
51 recv_chunk = recv_type->get_extent();
53 send_chunk *= send_count;
54 recv_chunk *= recv_count;
56 Request::sendrecv(send_ptr + rank * send_chunk, send_count, send_type, rank, tag,
57 recv_ptr + rank * recv_chunk, recv_count, recv_type, rank, tag,
60 for (i = 1; i < num_procs; i++) {
61 src = (rank - i + num_procs) % num_procs;
62 dst = (rank + i) % num_procs;
64 Request::sendrecv(send_ptr + dst * send_chunk, send_count, send_type,
65 dst, tag, recv_ptr + src * recv_chunk, recv_count,
66 recv_type, src, tag, comm, &s);
68 if ((i + 1) < num_procs) {
69 next_src = (rank - (i + 1) + num_procs) % num_procs;
70 next_dst = (rank + (i + 1) + num_procs) % num_procs;
71 Request::sendrecv(&send_sync, 1, MPI_CHAR, next_src, tag,
72 &recv_sync, 1, MPI_CHAR, next_dst, tag, comm, &s);