Logo AND Algorithmique Numérique Distribuée

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