Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' of git+ssh://scm.gforge.inria.fr//gitroot/simgrid/simgrid
[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
26 smpi_coll_tuned_alltoall_simple(void * send_buff, int send_count,
27                 MPI_Datatype send_type, void * recv_buff,
28                 int recv_count, MPI_Datatype recv_type,
29                 MPI_Comm comm)
30 {
31   int i, rank, size, nreqs, err, src, dst, tag = 101;
32   char *psnd;
33   char *prcv;
34   MPI_Aint sndinc;
35   MPI_Aint rcvinc;
36   MPI_Request *req;
37   MPI_Request *preq;
38   MPI_Request *qreq;
39   MPI_Status s, * statuses;
40
41
42   MPI_Comm_size(comm, &size);
43   MPI_Comm_rank(comm, &rank);
44   MPI_Type_extent(send_type, &sndinc);
45   MPI_Type_extent(recv_type, &rcvinc);
46   sndinc *= send_count;
47   rcvinc *= recv_count;
48
49   /* Allocate arrays of requests. */
50
51   nreqs = 2 * (size - 1);
52   if (nreqs > 0)
53     {
54       req = (MPI_Request *) malloc(nreqs * sizeof(MPI_Request));
55       statuses = (MPI_Status *) malloc(nreqs * sizeof(MPI_Status));
56       if (!req || !statuses)
57         {
58           free(req);
59           free(statuses);
60           return 0;
61         }
62     }
63   else
64     req = 0;
65
66   /* simple optimization */
67
68   psnd = ((char *) send_buff) + (rank * sndinc);
69   prcv = ((char *) recv_buff) + (rank * rcvinc);
70   MPI_Sendrecv (psnd, send_count, send_type, rank, tag,
71                 prcv, recv_count, recv_type,
72                 rank, tag, comm, &s);
73
74
75   /* Initiate all send/recv to/from others. */
76
77   preq = req;
78   qreq = req + size - 1;
79   prcv = (char*) recv_buff;
80   psnd = (char*) send_buff;
81   for (i = 0; i < size; i++)
82     {
83       src = dst = (rank + i) % size;
84       if (src == rank) continue;
85       if (dst == rank) continue;      
86       MPI_Recv_init(prcv + (src * rcvinc), recv_count, recv_type, src,
87                     tag, comm, preq++);
88       MPI_Send_init(psnd + (dst * sndinc), send_count, send_type, dst,
89                     tag, comm, qreq++);
90     }
91
92   /* Start all the requests. */
93
94   err = MPI_Startall(nreqs, req);
95
96   /* Wait for them all. */
97
98   err = MPI_Waitall(nreqs, req, statuses);
99
100   if (err != MPI_SUCCESS) {
101     if (req)
102       free((char *) req);
103     return err;
104   }
105
106   for (i = 0, preq = req; i < nreqs; ++i, ++preq) {
107     err = MPI_Request_free(preq);
108     if (err != MPI_SUCCESS) {
109       if (req)
110         free((char *) req);
111       if (statuses)
112         free(statuses);
113       return err;
114     }
115   }
116
117   /* All done */
118
119   if (req)
120     free((char *) req);
121   if (statuses)
122     free(statuses);
123   return (1);
124 }
125
126