Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Reindent files before changes.
[simgrid.git] / src / smpi / colls / alltoall-pair.c
1 #include "smpi/mpi.h"
2
3 /*****************************************************************************
4
5  * Function: alltoall_pair
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: Function works when P is power of two. In each phase of P - 1
19            phases, nodes in pair communicate their data.
20
21  * Auther: Ahmad Faraj
22
23  ****************************************************************************/
24 /*
25 int alltoall_pair(void *send_buff, int send_count, MPI_Datatype send_type,
26                   void *recv_buff, int recv_count, MPI_Datatype recv_type,
27                   MPI_Comm comm)
28 {
29
30   MPI_Aint send_chunk, recv_chunk;
31   MPI_Status s;
32   MPI_Win win;
33   int assert = 0;
34   int i, src, dst, rank, num_procs;
35   int tag = 1, success = 1, failure = 0, pof2 = 1;
36
37   char *send_ptr = (char *) send_buff;
38   char *recv_ptr = (char *) recv_buff;
39
40   MPI_Comm_rank(comm, &rank);
41   MPI_Comm_size(comm, &num_procs);
42   MPI_Type_extent(send_type, &send_chunk);
43   MPI_Type_extent(recv_type, &recv_chunk);
44
45   MPI_Win_create(recv_buff, num_procs * recv_chunk * send_count, recv_chunk, 0,
46                  comm, &win);
47   send_chunk *= send_count;
48   recv_chunk *= recv_count;
49
50   MPI_Win_fence(assert, win);
51   for (i = 0; i < num_procs; i++) {
52     src = dst = rank ^ i;
53     MPI_Put(send_ptr + dst * send_chunk, send_count, send_type, dst,
54             rank * send_chunk, send_count, send_type, win);
55   }
56   MPI_Win_fence(assert, win);
57   MPI_Win_free(&win);
58   return 0;
59 }
60 */
61
62 int smpi_coll_tuned_alltoall_pair(void *send_buff, int send_count,
63                                   MPI_Datatype send_type,
64                                   void *recv_buff, int recv_count,
65                                   MPI_Datatype recv_type,
66                                   MPI_Comm comm)
67 {
68
69   MPI_Aint send_chunk, recv_chunk;
70   MPI_Status s;
71   int i, src, dst, rank, num_procs;
72   int tag = 1, success = 1;
73
74   char *send_ptr = (char *) send_buff;
75   char *recv_ptr = (char *) recv_buff;
76
77   MPI_Comm_rank(comm, &rank);
78   MPI_Comm_size(comm, &num_procs);
79   MPI_Type_extent(send_type, &send_chunk);
80   MPI_Type_extent(recv_type, &recv_chunk);
81
82   send_chunk *= send_count;
83   recv_chunk *= recv_count;
84
85   for (i = 0; i < num_procs; i++) {
86     src = dst = rank ^ i;
87     MPI_Sendrecv(send_ptr + dst * send_chunk, send_count, send_type, dst,
88                  tag, recv_ptr + src * recv_chunk, recv_count, recv_type,
89                  src, tag, comm, &s);
90   }
91
92   return success;
93 }