Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
7abe6e0d2b63c2f20f86b764975102bad73d92c0
[simgrid.git] / src / smpi / colls / allgather-loosely-lr.c
1 #include "colls_private.h"
2
3 #ifndef NUM_CORE
4 #define NUM_CORE 4
5 #endif
6
7 int smpi_coll_tuned_allgather_loosely_lr(void *sbuf, int scount,
8                                          MPI_Datatype stype, void *rbuf,
9                                          int rcount, MPI_Datatype rtype,
10                                          MPI_Comm comm)
11 {
12   int comm_size, rank;
13   int tag = COLL_TAG_ALLGATHER;
14   int i, j, send_offset, recv_offset;
15   int intra_rank, inter_rank, inter_comm_size, intra_comm_size;
16   int inter_dst, inter_src;
17
18   comm_size = smpi_comm_size(comm);
19
20   if(comm_size%4)
21     THROWF(arg_error,0, "allgather loosely lr algorithm can't be used with non multiple of NUM_CORE=4 number of processes ! ");
22
23   rank = smpi_comm_rank(comm);
24   MPI_Aint rextent, sextent;
25   rextent = smpi_datatype_get_extent(rtype);
26   sextent = smpi_datatype_get_extent(stype);
27   MPI_Request inter_rrequest;
28   MPI_Request rrequest_array[128];
29   MPI_Request srequest_array[128];
30   MPI_Request inter_srequest_array[128];
31
32
33   int rrequest_count = 0;
34   int srequest_count = 0;
35   int inter_srequest_count = 0;
36
37   MPI_Status status;
38
39   intra_rank = rank % NUM_CORE;
40   inter_rank = rank / NUM_CORE;
41   inter_comm_size = (comm_size + NUM_CORE - 1) / NUM_CORE;
42   intra_comm_size = NUM_CORE;
43
44   int src_seg, dst_seg;
45
46   //copy corresponding message from sbuf to rbuf
47   recv_offset = rank * rextent * rcount;
48   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
49                (char *)rbuf + recv_offset, rcount, rtype, rank, tag, comm, &status);
50
51   int dst, src;
52   int inter_send_offset, inter_recv_offset;
53
54   rrequest_count = 0;
55   srequest_count = 0;
56   inter_srequest_count = 0;
57
58   for (i = 0; i < inter_comm_size; i++) {
59
60     // inter_communication
61
62     inter_dst = (rank + intra_comm_size) % comm_size;
63     inter_src = (rank - intra_comm_size + comm_size) % comm_size;
64
65     src_seg =
66         ((inter_rank - 1 - i +
67           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
68     dst_seg =
69         ((inter_rank - i +
70           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
71
72     inter_send_offset = dst_seg * sextent * scount;
73     inter_recv_offset = src_seg * rextent * rcount;
74
75     for (j = 0; j < intra_comm_size; j++) {
76
77       // inter communication
78       if (intra_rank == j) {
79         if (i != inter_comm_size - 1) {
80
81           inter_rrequest = smpi_mpi_irecv((char *)rbuf + inter_recv_offset, rcount, rtype,
82                                           inter_src, tag, comm);
83           inter_srequest_array[inter_srequest_count++] = smpi_mpi_isend((char *)rbuf + inter_send_offset, scount, stype,
84                                                                         inter_dst, tag, comm);
85         }
86       }
87       //intra_communication
88       src = inter_rank * intra_comm_size + j;
89       dst = inter_rank * intra_comm_size + j;
90
91       src_seg =
92           ((inter_rank - i +
93             inter_comm_size) % inter_comm_size) * intra_comm_size + j;
94       dst_seg =
95           ((inter_rank - i +
96             inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
97
98       send_offset = dst_seg * sextent * scount;
99       recv_offset = src_seg * rextent * rcount;
100
101
102       if (j != intra_rank) {
103
104         rrequest_array[rrequest_count++] = smpi_mpi_irecv((char *)rbuf + recv_offset, rcount, rtype, src, tag, comm);
105         srequest_array[srequest_count++] = smpi_mpi_isend((char *)rbuf + send_offset, scount, stype, dst, tag, comm);
106
107       }
108     }                           // intra loop
109
110
111     // wait for inter communication to finish for these rounds (# of round equals NUM_CORE)
112     if (i != inter_comm_size - 1) {
113       smpi_mpi_wait(&inter_rrequest, &status);
114     }
115
116   }                             //inter loop
117
118   smpi_mpi_waitall(rrequest_count, rrequest_array, MPI_STATUSES_IGNORE);
119   smpi_mpi_waitall(srequest_count, srequest_array, MPI_STATUSES_IGNORE);
120   smpi_mpi_waitall(inter_srequest_count, inter_srequest_array, MPI_STATUSES_IGNORE);
121
122   return MPI_SUCCESS;
123 }