Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
protect (hopefully) collective communication algorithms from abuse.
[simgrid.git] / src / smpi / colls / allgather-smp-simple.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_simple(void *send_buf, int scount,
7                                          MPI_Datatype stype, void *recv_buf,
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
14   if(comm_size%NUM_CORE)
15      THROWF(arg_error,0, "allgather SMP simple algorithm can't be used with non multiple of NUM_CORE=%d number of processes ! ", NUM_CORE);
16
17   rank = smpi_comm_rank(comm);
18   MPI_Aint rextent, sextent;
19   rextent = smpi_datatype_get_extent(rtype);
20   sextent = smpi_datatype_get_extent(stype);
21   int tag = COLL_TAG_ALLGATHER;
22   MPI_Status status;
23   int i, send_offset, recv_offset;
24   int intra_rank, inter_rank;
25   int num_core = NUM_CORE;
26   intra_rank = rank % num_core;
27   inter_rank = rank / num_core;
28   int inter_comm_size = (comm_size + num_core - 1) / num_core;
29   int num_core_in_current_smp = num_core;
30
31   // the last SMP node may have fewer number of running processes than all others
32   if (inter_rank == (inter_comm_size - 1)) {
33     num_core_in_current_smp = comm_size - (inter_rank * num_core);
34   }
35   //INTRA-SMP-ALLGATHER
36   recv_offset = rank * rextent * rcount;
37   smpi_mpi_sendrecv(send_buf, scount, stype, rank, tag,
38                ((char *) recv_buf + recv_offset), rcount, rtype, rank, tag,
39                comm, &status);
40   for (i = 1; i < num_core_in_current_smp; i++) {
41
42     dst =
43         (inter_rank * num_core) + (intra_rank + i) % (num_core_in_current_smp);
44     src =
45         (inter_rank * num_core) + (intra_rank - i +
46                                    num_core_in_current_smp) %
47         (num_core_in_current_smp);
48     recv_offset = src * rextent * rcount;
49
50     smpi_mpi_sendrecv(send_buf, scount, stype, dst, tag,
51                  ((char *) recv_buf + recv_offset), rcount, rtype, src, tag,
52                  comm, &status);
53
54   }
55
56   // INTER-SMP-ALLGATHER 
57   // Every root of each SMP node post INTER-Sendrecv, then do INTRA-Bcast for each receiving message
58
59
60
61   if (intra_rank == 0) {
62     MPI_Request *reqs, *req_ptr;
63     int num_req = (inter_comm_size - 1) * 2;
64     reqs = (MPI_Request *) xbt_malloc(num_req * sizeof(MPI_Request));
65     req_ptr = reqs;
66     MPI_Status *stat;
67     stat = (MPI_Status *) xbt_malloc(num_req * sizeof(MPI_Status));
68
69     for (i = 1; i < inter_comm_size; i++) {
70
71       //dst = ((inter_rank+i)%inter_comm_size) * num_core;
72       src = ((inter_rank - i + inter_comm_size) % inter_comm_size) * num_core;
73       //send_offset = (rank * sextent * scount);
74       recv_offset = (src * sextent * scount);
75       //      smpi_mpi_sendrecv((recv_buf+send_offset), (scount * num_core), stype, dst, tag, 
76       //             (recv_buf+recv_offset), (rcount * num_core), rtype, src, tag, comm, &status);
77       //MPIC_Isend((recv_buf+send_offset), (scount * num_core), stype, dst, tag, comm, req_ptr++);
78       *(req_ptr++) = smpi_mpi_irecv(((char *) recv_buf + recv_offset), (rcount * num_core), rtype,
79                 src, tag, comm);
80     }
81     for (i = 1; i < inter_comm_size; i++) {
82
83       dst = ((inter_rank + i) % inter_comm_size) * num_core;
84       //src = ((inter_rank-i+inter_comm_size)%inter_comm_size) * num_core;
85       send_offset = (rank * sextent * scount);
86       //recv_offset = (src * sextent * scount);
87       //      smpi_mpi_sendrecv((recv_buf+send_offset), (scount * num_core), stype, dst, tag, 
88       //             (recv_buf+recv_offset), (rcount * num_core), rtype, src, tag, comm, &status);
89       *(req_ptr++) = smpi_mpi_isend(((char *) recv_buf + send_offset), (scount * num_core), stype,
90                 dst, tag, comm);
91       //MPIC_Irecv((recv_buf+recv_offset), (rcount * num_core), rtype, src, tag, comm, req_ptr++);
92     }
93     smpi_mpi_waitall(num_req, reqs, stat);
94     free(reqs);
95     free(stat);
96
97   }
98   //INTRA-BCAST (use flat tree)
99
100   if (intra_rank == 0) {
101     for (i = 1; i < num_core_in_current_smp; i++) {
102       //printf("rank = %d, num = %d send to %d\n",rank, num_core_in_current_smp, (rank + i));
103       smpi_mpi_send(recv_buf, (scount * comm_size), stype, (rank + i), tag, comm);
104     }
105   } else {
106     //printf("rank = %d recv from %d\n",rank, (inter_rank * num_core));
107     smpi_mpi_recv(recv_buf, (rcount * comm_size), rtype, (inter_rank * num_core),
108              tag, comm, &status);
109   }
110
111
112   return MPI_SUCCESS;
113 }