Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
have smp-aware algorithms use number of cores on the node as basis for their computat...
[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   int num_core = simcall_host_get_core(SIMIX_host_self());
21   // do we use the default one or the number of cores in the platform ?
22   // if the number of cores is one, the platform may be simulated with 1 node = 1 core
23   if (num_core == 1) num_core = NUM_CORE;
24
25   if(comm_size%num_core)
26     THROWF(arg_error,0, "allgather loosely lr algorithm can't be used with non multiple of NUM_CORE=%d number of processes ! ",num_core);
27
28   rank = smpi_comm_rank(comm);
29   MPI_Aint rextent, sextent;
30   rextent = smpi_datatype_get_extent(rtype);
31   sextent = smpi_datatype_get_extent(stype);
32   MPI_Request inter_rrequest;
33   MPI_Request rrequest_array[128];
34   MPI_Request srequest_array[128];
35   MPI_Request inter_srequest_array[128];
36
37
38   int rrequest_count = 0;
39   int srequest_count = 0;
40   int inter_srequest_count = 0;
41
42   MPI_Status status;
43
44   intra_rank = rank % num_core;
45   inter_rank = rank / num_core;
46   inter_comm_size = (comm_size + num_core - 1) / num_core;
47   intra_comm_size = num_core;
48
49   int src_seg, dst_seg;
50
51   //copy corresponding message from sbuf to rbuf
52   recv_offset = rank * rextent * rcount;
53   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
54                (char *)rbuf + recv_offset, rcount, rtype, rank, tag, comm, &status);
55
56   int dst, src;
57   int inter_send_offset, inter_recv_offset;
58
59   rrequest_count = 0;
60   srequest_count = 0;
61   inter_srequest_count = 0;
62
63   for (i = 0; i < inter_comm_size; i++) {
64
65     // inter_communication
66
67     inter_dst = (rank + intra_comm_size) % comm_size;
68     inter_src = (rank - intra_comm_size + comm_size) % comm_size;
69
70     src_seg =
71         ((inter_rank - 1 - i +
72           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
73     dst_seg =
74         ((inter_rank - i +
75           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
76
77     inter_send_offset = dst_seg * sextent * scount;
78     inter_recv_offset = src_seg * rextent * rcount;
79
80     for (j = 0; j < intra_comm_size; j++) {
81
82       // inter communication
83       if (intra_rank == j) {
84         if (i != inter_comm_size - 1) {
85
86           inter_rrequest = smpi_mpi_irecv((char *)rbuf + inter_recv_offset, rcount, rtype,
87                                           inter_src, tag, comm);
88           inter_srequest_array[inter_srequest_count++] = smpi_mpi_isend((char *)rbuf + inter_send_offset, scount, stype,
89                                                                         inter_dst, tag, comm);
90         }
91       }
92       //intra_communication
93       src = inter_rank * intra_comm_size + j;
94       dst = inter_rank * intra_comm_size + j;
95
96       src_seg =
97           ((inter_rank - i +
98             inter_comm_size) % inter_comm_size) * intra_comm_size + j;
99       dst_seg =
100           ((inter_rank - i +
101             inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
102
103       send_offset = dst_seg * sextent * scount;
104       recv_offset = src_seg * rextent * rcount;
105
106
107       if (j != intra_rank) {
108
109         rrequest_array[rrequest_count++] = smpi_mpi_irecv((char *)rbuf + recv_offset, rcount, rtype, src, tag, comm);
110         srequest_array[srequest_count++] = smpi_mpi_isend((char *)rbuf + send_offset, scount, stype, dst, tag, comm);
111
112       }
113     }                           // intra loop
114
115
116     // wait for inter communication to finish for these rounds (# of round equals num_core)
117     if (i != inter_comm_size - 1) {
118       smpi_mpi_wait(&inter_rrequest, &status);
119     }
120
121   }                             //inter loop
122
123   smpi_mpi_waitall(rrequest_count, rrequest_array, MPI_STATUSES_IGNORE);
124   smpi_mpi_waitall(srequest_count, srequest_array, MPI_STATUSES_IGNORE);
125   smpi_mpi_waitall(inter_srequest_count, inter_srequest_array, MPI_STATUSES_IGNORE);
126
127   return MPI_SUCCESS;
128 }