Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Use simgrid function instead of MPI in collectives
[simgrid.git] / src / smpi / colls / alltoall-simple.c
1 #include "colls_private.h"
2
3 /*****************************************************************************
4
5  * Function: alltoall_spreading_simple
6
7  * Return: int
8
9  *  Inputs:
10     send_buff: send input buffer
11     send_count: number of elements to send
12     send_type: data type of elements being sent
13     recv_buff: receive output buffer
14     recv_count: number of elements to received
15     recv_type: data type of elements being received
16     comm: communicator
17
18  * Descrp: Let i -> j denote the communication from node i to node j. The
19           order of communications for node i is i -> i + 1, i -> i + 2, ...,
20           i -> (i + p -1) % P.
21  
22  * Auther: Ahmad Faraj
23
24  ****************************************************************************/
25 int smpi_coll_tuned_alltoall_simple(void *send_buff, int send_count,
26                                     MPI_Datatype send_type,
27                                     void *recv_buff, int recv_count,
28                                     MPI_Datatype recv_type, MPI_Comm comm)
29 {
30   int i, rank, size, nreqs, src, dst, tag = 101;
31   char *psnd;
32   char *prcv;
33   MPI_Aint sndinc;
34   MPI_Aint rcvinc;
35   MPI_Request *req;
36   MPI_Request *preq;
37   MPI_Request *qreq;
38   MPI_Status s, *statuses;
39
40
41   size = smpi_comm_size(comm);
42   rank = smpi_comm_rank(comm);
43   sndinc = smpi_datatype_get_extent(send_type);
44   rcvinc = smpi_datatype_get_extent(recv_type);
45   sndinc *= send_count;
46   rcvinc *= recv_count;
47
48   /* Allocate arrays of requests. */
49
50   nreqs = 2 * (size - 1);
51   if (nreqs > 0) {
52     req = (MPI_Request *) xbt_malloc(nreqs * sizeof(MPI_Request));
53     statuses = (MPI_Status *) xbt_malloc(nreqs * sizeof(MPI_Status));
54   } else {
55     req = NULL;
56     statuses = NULL;
57   }
58
59   /* simple optimization */
60
61   psnd = ((char *) send_buff) + (rank * sndinc);
62   prcv = ((char *) recv_buff) + (rank * rcvinc);
63   smpi_mpi_sendrecv(psnd, send_count, send_type, rank, tag,
64                prcv, recv_count, recv_type, rank, tag, comm, &s);
65
66
67   /* Initiate all send/recv to/from others. */
68
69   preq = req;
70   qreq = req + size - 1;
71   prcv = (char *) recv_buff;
72   psnd = (char *) send_buff;
73   for (i = 0; i < size; i++) {
74     src = dst = (rank + i) % size;
75     if (src == rank)
76       continue;
77     if (dst == rank)
78       continue;
79     *(preq++) = smpi_mpi_recv_init(prcv + (src * rcvinc), recv_count, recv_type, src,
80                   tag, comm);
81     *(qreq++) = smpi_mpi_send_init(psnd + (dst * sndinc), send_count, send_type, dst,
82                   tag, comm);
83   }
84
85   /* Start all the requests. */
86
87   smpi_mpi_startall(nreqs, req);
88
89   /* Wait for them all. */
90
91   smpi_mpi_waitall(nreqs, req, statuses);
92
93   for (i = 0, preq = req; i < nreqs; ++i, ++preq) {
94     smpi_mpi_request_free(preq);
95   }
96
97   /* All done */
98
99   if (req)
100     free((char *) req);
101   if (statuses)
102     free(statuses);
103   return MPI_SUCCESS;
104 }