Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
another allreduce ported
[simgrid.git] / src / smpi / colls / alltoall-simple.c
1 #include "colls.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, err, 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   MPI_Comm_size(comm, &size);
42   MPI_Comm_rank(comm, &rank);
43   MPI_Type_extent(send_type, &sndinc);
44   MPI_Type_extent(recv_type, &rcvinc);
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 *) malloc(nreqs * sizeof(MPI_Request));
53     statuses = (MPI_Status *) malloc(nreqs * sizeof(MPI_Status));
54     if (!req || !statuses) {
55       free(req);
56       free(statuses);
57       return 0;
58     }
59   } else {
60     req = NULL;
61     statuses = NULL;
62   }
63
64   /* simple optimization */
65
66   psnd = ((char *) send_buff) + (rank * sndinc);
67   prcv = ((char *) recv_buff) + (rank * rcvinc);
68   MPI_Sendrecv(psnd, send_count, send_type, rank, tag,
69                prcv, recv_count, recv_type, rank, tag, comm, &s);
70
71
72   /* Initiate all send/recv to/from others. */
73
74   preq = req;
75   qreq = req + size - 1;
76   prcv = (char *) recv_buff;
77   psnd = (char *) send_buff;
78   for (i = 0; i < size; i++) {
79     src = dst = (rank + i) % size;
80     if (src == rank)
81       continue;
82     if (dst == rank)
83       continue;
84     MPI_Recv_init(prcv + (src * rcvinc), recv_count, recv_type, src,
85                   tag, comm, preq++);
86     MPI_Send_init(psnd + (dst * sndinc), send_count, send_type, dst,
87                   tag, comm, qreq++);
88   }
89
90   /* Start all the requests. */
91
92   err = MPI_Startall(nreqs, req);
93
94   /* Wait for them all. */
95
96   err = MPI_Waitall(nreqs, req, statuses);
97
98   if (err != MPI_SUCCESS) {
99     if (req)
100       free((char *) req);
101     return err;
102   }
103
104   for (i = 0, preq = req; i < nreqs; ++i, ++preq) {
105     err = MPI_Request_free(preq);
106     if (err != MPI_SUCCESS) {
107       if (req)
108         free((char *) req);
109       if (statuses)
110         free(statuses);
111       return err;
112     }
113   }
114
115   /* All done */
116
117   if (req)
118     free((char *) req);
119   if (statuses)
120     free(statuses);
121   return (1);
122 }