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