Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
6455db135aebb0a8a08af9b98c1e3facf261619f
[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 = 50;
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   rank = smpi_comm_rank(comm);
20   MPI_Aint rextent, sextent;
21   rextent = smpi_datatype_get_extent(rtype);
22   sextent = smpi_datatype_get_extent(stype);
23   MPI_Request inter_rrequest;
24   MPI_Request rrequest_array[128];
25   MPI_Request srequest_array[128];
26   MPI_Request inter_srequest_array[128];
27
28
29   int rrequest_count = 0;
30   int srequest_count = 0;
31   int inter_srequest_count = 0;
32
33   MPI_Status status;
34
35   intra_rank = rank % NUM_CORE;
36   inter_rank = rank / NUM_CORE;
37   inter_comm_size = (comm_size + NUM_CORE - 1) / NUM_CORE;
38   intra_comm_size = NUM_CORE;
39
40   int src_seg, dst_seg;
41
42   //copy corresponding message from sbuf to rbuf
43   recv_offset = rank * rextent * rcount;
44   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
45                (char *)rbuf + recv_offset, rcount, rtype, rank, tag, comm, &status);
46
47   int dst, src;
48   int inter_send_offset, inter_recv_offset;
49
50   rrequest_count = 0;
51   srequest_count = 0;
52   inter_srequest_count = 0;
53
54   for (i = 0; i < inter_comm_size; i++) {
55
56     // inter_communication
57
58     inter_dst = (rank + intra_comm_size) % comm_size;
59     inter_src = (rank - intra_comm_size + comm_size) % comm_size;
60
61     src_seg =
62         ((inter_rank - 1 - i +
63           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
64     dst_seg =
65         ((inter_rank - i +
66           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
67
68     inter_send_offset = dst_seg * sextent * scount;
69     inter_recv_offset = src_seg * rextent * rcount;
70
71     for (j = 0; j < intra_comm_size; j++) {
72
73       // inter communication
74       if (intra_rank == j) {
75         if (i != inter_comm_size - 1) {
76
77           inter_rrequest = smpi_mpi_irecv((char *)rbuf + inter_recv_offset, rcount, rtype,
78                                           inter_src, tag, comm);
79           inter_srequest_array[inter_srequest_count++] = smpi_mpi_isend((char *)rbuf + inter_send_offset, scount, stype,
80                                                                         inter_dst, tag, comm);
81         }
82       }
83       //intra_communication
84       src = inter_rank * intra_comm_size + j;
85       dst = inter_rank * intra_comm_size + j;
86
87       src_seg =
88           ((inter_rank - i +
89             inter_comm_size) % inter_comm_size) * intra_comm_size + j;
90       dst_seg =
91           ((inter_rank - i +
92             inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
93
94       send_offset = dst_seg * sextent * scount;
95       recv_offset = src_seg * rextent * rcount;
96
97
98       if (j != intra_rank) {
99
100         rrequest_array[rrequest_count++] = smpi_mpi_irecv((char *)rbuf + recv_offset, rcount, rtype, src, tag, comm);
101         srequest_array[srequest_count++] = smpi_mpi_isend((char *)rbuf + send_offset, scount, stype, dst, tag, comm);
102
103       }
104     }                           // intra loop
105
106
107     // wait for inter communication to finish for these rounds (# of round equals NUM_CORE)
108     if (i != inter_comm_size - 1) {
109       smpi_mpi_wait(&inter_rrequest, &status);
110     }
111
112   }                             //inter loop
113
114   smpi_mpi_waitall(rrequest_count, rrequest_array, MPI_STATUSES_IGNORE);
115   smpi_mpi_waitall(srequest_count, srequest_array, MPI_STATUSES_IGNORE);
116   smpi_mpi_waitall(inter_srequest_count, inter_srequest_array, MPI_STATUSES_IGNORE);
117
118   return MPI_SUCCESS;
119 }