Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add missing file in DefinePackage
[simgrid.git] / src / smpi / colls / allgather-SMP-NTS.c
1 #include "colls_private.h"
2 #ifndef NUM_CORE
3 #define NUM_CORE 8
4 #endif
5
6 int smpi_coll_tuned_allgather_SMP_NTS(void *sbuf, int scount,
7                                       MPI_Datatype stype, void *rbuf,
8                                       int rcount, MPI_Datatype rtype,
9                                       MPI_Comm comm)
10 {
11   int src, dst, comm_size, rank;
12   comm_size = smpi_comm_size(comm);
13   rank = smpi_comm_rank(comm);
14   MPI_Aint rextent, sextent;
15   rextent = smpi_datatype_get_extent(rtype);
16   sextent = smpi_datatype_get_extent(stype);
17   int tag = 50;
18   MPI_Request request;
19   MPI_Request rrequest_array[128];
20   MPI_Request srequest_array[128];
21
22   MPI_Status status;
23   int i, send_offset, recv_offset;
24   int intra_rank, inter_rank;
25   intra_rank = rank % NUM_CORE;
26   inter_rank = rank / NUM_CORE;
27   int inter_comm_size = (comm_size + NUM_CORE - 1) / NUM_CORE;
28   int num_core_in_current_smp = NUM_CORE;
29
30   /* for too small number of processes, use default implementation */
31   if (comm_size <= NUM_CORE) {
32     XBT_WARN("MPI_allgather_SMP_NTS use default MPI_allgather.");         
33     smpi_mpi_allgather(sbuf, scount, stype, rbuf, rcount, rtype, comm);
34     return MPI_SUCCESS;    
35   }
36
37   // the last SMP node may have fewer number of running processes than all others
38   if (inter_rank == (inter_comm_size - 1)) {
39     num_core_in_current_smp = comm_size - (inter_rank * NUM_CORE);
40   }
41   //copy corresponding message from sbuf to rbuf
42   recv_offset = rank * rextent * rcount;
43   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
44                ((char *) rbuf + recv_offset), rcount, rtype, rank, tag, comm,
45                &status);
46
47   //gather to root of each SMP
48
49   for (i = 1; i < num_core_in_current_smp; i++) {
50
51     dst =
52         (inter_rank * NUM_CORE) + (intra_rank + i) % (num_core_in_current_smp);
53     src =
54         (inter_rank * NUM_CORE) + (intra_rank - i +
55                                    num_core_in_current_smp) %
56         (num_core_in_current_smp);
57     recv_offset = src * rextent * rcount;
58
59     smpi_mpi_sendrecv(sbuf, scount, stype, dst, tag,
60                  ((char *) rbuf + recv_offset), rcount, rtype, src, tag, comm,
61                  &status);
62
63   }
64
65   // INTER-SMP-ALLGATHER 
66   // Every root of each SMP node post INTER-Sendrecv, then do INTRA-Bcast for each receiving message
67   // Use logical ring algorithm
68
69   // root of each SMP
70   if (intra_rank == 0) {
71     src = ((inter_rank - 1 + inter_comm_size) % inter_comm_size) * NUM_CORE;
72     dst = ((inter_rank + 1) % inter_comm_size) * NUM_CORE;
73
74     // post all inter Irecv
75     for (i = 0; i < inter_comm_size - 1; i++) {
76       recv_offset =
77           ((inter_rank - i - 1 +
78             inter_comm_size) % inter_comm_size) * NUM_CORE * sextent * scount;
79       rrequest_array[i] = smpi_mpi_irecv((char *)rbuf+recv_offset, rcount * NUM_CORE, rtype, src, tag+i, comm);
80     }
81
82     // send first message
83     send_offset =
84         ((inter_rank +
85           inter_comm_size) % inter_comm_size) * NUM_CORE * sextent * scount;
86     srequest_array[0] = smpi_mpi_isend((char *) rbuf + send_offset, scount * NUM_CORE, stype, dst, tag,
87               comm);
88
89     // loop : recv-inter , send-inter, send-intra (linear-bcast)
90     for (i = 0; i < inter_comm_size - 2; i++) {
91       recv_offset =
92           ((inter_rank - i - 1 +
93             inter_comm_size) % inter_comm_size) * NUM_CORE * sextent * scount;
94       smpi_mpi_wait(&rrequest_array[i], &status);
95       srequest_array[i + 1] = smpi_mpi_isend((char *) rbuf + recv_offset, scount * NUM_CORE, stype, dst,
96                 tag + i + 1, comm);
97       if (num_core_in_current_smp > 1) {
98         request = smpi_mpi_isend((char *) rbuf + recv_offset, scount * NUM_CORE, stype,
99                   (rank + 1), tag + i + 1, comm);
100       }
101     }
102
103     // recv last message and send_intra
104     recv_offset =
105         ((inter_rank - i - 1 +
106           inter_comm_size) % inter_comm_size) * NUM_CORE * sextent * scount;
107     //recv_offset = ((inter_rank + 1) % inter_comm_size) * NUM_CORE * sextent * scount;
108     //i=inter_comm_size-2;
109     smpi_mpi_wait(&rrequest_array[i], &status);
110     if (num_core_in_current_smp > 1) {
111       request = smpi_mpi_isend((char *) rbuf + recv_offset, scount * NUM_CORE, stype,
112                 (rank + 1), tag + i + 1, comm);
113     }
114   }
115   // last rank of each SMP
116   else if (intra_rank == (num_core_in_current_smp - 1)) {
117     for (i = 0; i < inter_comm_size - 1; i++) {
118       recv_offset =
119           ((inter_rank - i - 1 +
120             inter_comm_size) % inter_comm_size) * NUM_CORE * sextent * scount;
121       request = smpi_mpi_irecv((char *) rbuf + recv_offset, (rcount * NUM_CORE), rtype,
122                 rank - 1, tag + i + 1, comm);
123       smpi_mpi_wait(&request, &status);
124     }
125   }
126   // intermediate rank of each SMP
127   else {
128     for (i = 0; i < inter_comm_size - 1; i++) {
129       recv_offset =
130           ((inter_rank - i - 1 +
131             inter_comm_size) % inter_comm_size) * NUM_CORE * sextent * scount;
132       request = smpi_mpi_irecv((char *) rbuf + recv_offset, (rcount * NUM_CORE), rtype,
133                 rank - 1, tag + i + 1, comm);
134       smpi_mpi_wait(&request, &status);
135       request = smpi_mpi_isend((char *) rbuf + recv_offset, (scount * NUM_CORE), stype,
136                 (rank + 1), tag + i + 1, comm);
137     }
138   }
139
140   return MPI_SUCCESS;
141 }