Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add/update copyright notices.
[simgrid.git] / src / smpi / colls / allgather-SMP-NTS.c
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "colls_private.h"
8 #ifndef NUM_CORE
9 #define NUM_CORE 8
10 #endif
11
12 int smpi_coll_tuned_allgather_SMP_NTS(void *sbuf, int scount,
13                                       MPI_Datatype stype, void *rbuf,
14                                       int rcount, MPI_Datatype rtype,
15                                       MPI_Comm comm)
16 {
17   int src, dst, comm_size, rank;
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   int tag = COLL_TAG_ALLGATHER;
24
25   int i, send_offset, recv_offset;
26   int intra_rank, inter_rank;
27
28   int num_core = simcall_host_get_core(SIMIX_host_self());
29   // do we use the default one or the number of cores in the platform ?
30   // if the number of cores is one, the platform may be simulated with 1 node = 1 core
31   if (num_core == 1) num_core = NUM_CORE;
32
33
34   intra_rank = rank % num_core;
35   inter_rank = rank / num_core;
36   int inter_comm_size = (comm_size + num_core - 1) / num_core;
37   int num_core_in_current_smp = num_core;
38
39   if(comm_size%num_core)
40     THROWF(arg_error,0, "allgather SMP NTS algorithm can't be used with non multiple of NUM_CORE=%d number of processes ! ", num_core);
41
42   /* for too small number of processes, use default implementation */
43   if (comm_size <= num_core) {
44     XBT_WARN("MPI_allgather_SMP_NTS use default MPI_allgather.");         
45     smpi_mpi_allgather(sbuf, scount, stype, rbuf, rcount, rtype, comm);
46     return MPI_SUCCESS;    
47   }
48
49   // the last SMP node may have fewer number of running processes than all others
50   if (inter_rank == (inter_comm_size - 1)) {
51     num_core_in_current_smp = comm_size - (inter_rank * num_core);
52   }
53   //copy corresponding message from sbuf to rbuf
54   recv_offset = rank * rextent * rcount;
55   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
56                ((char *) rbuf + recv_offset), rcount, rtype, rank, tag, comm,
57                MPI_STATUS_IGNORE);
58
59   //gather to root of each SMP
60
61   for (i = 1; i < num_core_in_current_smp; i++) {
62
63     dst =
64         (inter_rank * num_core) + (intra_rank + i) % (num_core_in_current_smp);
65     src =
66         (inter_rank * num_core) + (intra_rank - i +
67                                    num_core_in_current_smp) %
68         (num_core_in_current_smp);
69     recv_offset = src * rextent * rcount;
70
71     smpi_mpi_sendrecv(sbuf, scount, stype, dst, tag,
72                  ((char *) rbuf + recv_offset), rcount, rtype, src, tag, comm,
73                  MPI_STATUS_IGNORE);
74
75   }
76
77   // INTER-SMP-ALLGATHER 
78   // Every root of each SMP node post INTER-Sendrecv, then do INTRA-Bcast for each receiving message
79   // Use logical ring algorithm
80
81   // root of each SMP
82   if (intra_rank == 0) {
83     MPI_Request *rrequest_array = xbt_new(MPI_Request, inter_comm_size - 1);
84     MPI_Request *srequest_array = xbt_new(MPI_Request, inter_comm_size - 1);
85
86     src = ((inter_rank - 1 + inter_comm_size) % inter_comm_size) * num_core;
87     dst = ((inter_rank + 1) % inter_comm_size) * num_core;
88
89     // post all inter Irecv
90     for (i = 0; i < inter_comm_size - 1; i++) {
91       recv_offset =
92           ((inter_rank - i - 1 +
93             inter_comm_size) % inter_comm_size) * num_core * sextent * scount;
94       rrequest_array[i] = smpi_mpi_irecv((char *)rbuf + recv_offset, rcount * num_core,
95                                          rtype, src, tag + i, comm);
96     }
97
98     // send first message
99     send_offset =
100         ((inter_rank +
101           inter_comm_size) % inter_comm_size) * num_core * sextent * scount;
102     srequest_array[0] = smpi_mpi_isend((char *)rbuf + send_offset, scount * num_core,
103                                        stype, dst, tag, comm);
104
105     // loop : recv-inter , send-inter, send-intra (linear-bcast)
106     for (i = 0; i < inter_comm_size - 2; i++) {
107       recv_offset =
108           ((inter_rank - i - 1 +
109             inter_comm_size) % inter_comm_size) * num_core * sextent * scount;
110       smpi_mpi_wait(&rrequest_array[i], MPI_STATUS_IGNORE);
111       srequest_array[i + 1] = smpi_mpi_isend((char *)rbuf + recv_offset, scount * num_core,
112                                              stype, dst, tag + i + 1, comm);
113       if (num_core_in_current_smp > 1) {
114         smpi_mpi_send((char *)rbuf + recv_offset, scount * num_core,
115                       stype, (rank + 1), tag + i + 1, comm);
116       }
117     }
118
119     // recv last message and send_intra
120     recv_offset =
121         ((inter_rank - i - 1 +
122           inter_comm_size) % inter_comm_size) * num_core * sextent * scount;
123     //recv_offset = ((inter_rank + 1) % inter_comm_size) * num_core * sextent * scount;
124     //i=inter_comm_size-2;
125     smpi_mpi_wait(&rrequest_array[i], MPI_STATUS_IGNORE);
126     if (num_core_in_current_smp > 1) {
127       smpi_mpi_send((char *)rbuf + recv_offset, scount * num_core,
128                                   stype, (rank + 1), tag + i + 1, comm);
129     }
130
131     smpi_mpi_waitall(inter_comm_size - 1, srequest_array, MPI_STATUSES_IGNORE);
132     xbt_free(rrequest_array);
133     xbt_free(srequest_array);
134   }
135   // last rank of each SMP
136   else if (intra_rank == (num_core_in_current_smp - 1)) {
137     for (i = 0; i < inter_comm_size - 1; i++) {
138       recv_offset =
139           ((inter_rank - i - 1 +
140             inter_comm_size) % inter_comm_size) * num_core * sextent * scount;
141       smpi_mpi_recv((char *) rbuf + recv_offset, (rcount * num_core), rtype,
142                     rank - 1, tag + i + 1, comm, MPI_STATUS_IGNORE);
143     }
144   }
145   // intermediate rank of each SMP
146   else {
147     for (i = 0; i < inter_comm_size - 1; i++) {
148       recv_offset =
149           ((inter_rank - i - 1 +
150             inter_comm_size) % inter_comm_size) * num_core * sextent * scount;
151       smpi_mpi_recv((char *) rbuf + recv_offset, (rcount * num_core), rtype,
152                     rank - 1, tag + i + 1, comm, MPI_STATUS_IGNORE);
153       smpi_mpi_send((char *) rbuf + recv_offset, (scount * num_core), stype,
154                     (rank + 1), tag + i + 1, comm);
155     }
156   }
157
158   return MPI_SUCCESS;
159 }