Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[Energy] Fix broken commit c6fc463d75e4daa67ddf9c460cbc473b1aa9d6ba
[simgrid.git] / src / smpi / colls / reduce-NTSL.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 //#include <star-reduction.c>
9
10 int reduce_NTSL_segment_size_in_byte = 8192;
11
12 /* Non-topology-specific pipelined linear-bcast function 
13    0->1, 1->2 ,2->3, ....., ->last node : in a pipeline fashion
14 */
15 int smpi_coll_tuned_reduce_NTSL(void *buf, void *rbuf, int count,
16                                 MPI_Datatype datatype, MPI_Op op, int root,
17                                 MPI_Comm comm)
18 {
19   int tag = COLL_TAG_REDUCE;
20   MPI_Status status;
21   MPI_Request *send_request_array;
22   MPI_Request *recv_request_array;
23   MPI_Status *send_status_array;
24   MPI_Status *recv_status_array;
25   int rank, size;
26   int i;
27   MPI_Aint extent;
28   extent = smpi_datatype_get_extent(datatype);
29
30   rank = smpi_comm_rank(comm);
31   size = smpi_comm_size(comm);
32
33   /* source node and destination nodes (same through out the functions) */
34   int to = (rank - 1 + size) % size;
35   int from = (rank + 1) % size;
36
37   /* segment is segment size in number of elements (not bytes) */
38   int segment = reduce_NTSL_segment_size_in_byte / extent;
39
40   /* pipeline length */
41   int pipe_length = count / segment;
42
43   /* use for buffer offset for sending and receiving data = segment size in byte */
44   int increment = segment * extent;
45
46   /* if the input size is not divisible by segment size => 
47      the small remainder will be done with native implementation */
48   int remainder = count % segment;
49
50   /* if root is not zero send to rank zero first
51      this can be modified to make it faster by using logical src, dst.
52    */
53
54   /*
55      if (root != 0) {
56      if (rank == root){
57      smpi_mpi_send(buf,count,datatype,0,tag,comm);
58      }
59      else if (rank == 0) {
60      smpi_mpi_recv(buf,count,datatype,root,tag,comm,&status);        
61      }
62      }
63    */
64
65   char *tmp_buf;
66   tmp_buf = (char *) smpi_get_tmp_sendbuffer(count * extent);
67
68   smpi_mpi_sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank,
69                tag, comm, &status);
70
71   /* when a message is smaller than a block size => no pipeline */
72   if (count <= segment) {
73     if (rank == root) {
74       smpi_mpi_recv(tmp_buf, count, datatype, from, tag, comm, &status);
75       smpi_op_apply(op, tmp_buf, rbuf, &count, &datatype);
76     } else if (rank == ((root - 1 + size) % size)) {
77       smpi_mpi_send(rbuf, count, datatype, to, tag, comm);
78     } else {
79       smpi_mpi_recv(tmp_buf, count, datatype, from, tag, comm, &status);
80       smpi_op_apply(op, tmp_buf, rbuf, &count, &datatype);
81       smpi_mpi_send(rbuf, count, datatype, to, tag, comm);
82     }
83     smpi_free_tmp_buffer(tmp_buf);
84     return MPI_SUCCESS;
85   }
86
87   /* pipeline */
88   else {
89     send_request_array =
90         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
91     recv_request_array =
92         (MPI_Request *) xbt_malloc((size + pipe_length) * sizeof(MPI_Request));
93     send_status_array =
94         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
95     recv_status_array =
96         (MPI_Status *) xbt_malloc((size + pipe_length) * sizeof(MPI_Status));
97
98     /* root recv data */
99     if (rank == root) {
100       for (i = 0; i < pipe_length; i++) {
101         recv_request_array[i] = smpi_mpi_irecv((char *) tmp_buf + (i * increment), segment, datatype, from,
102                   (tag + i), comm);
103       }
104       for (i = 0; i < pipe_length; i++) {
105         smpi_mpi_wait(&recv_request_array[i], &status);
106         smpi_op_apply(op, tmp_buf + (i * increment), (char *)rbuf + (i * increment),
107                        &segment, &datatype);
108       }
109     }
110
111     /* last node only sends data */
112     else if (rank == ((root - 1 + size) % size)) {
113       for (i = 0; i < pipe_length; i++) {
114         send_request_array[i] = smpi_mpi_isend((char *)rbuf + (i * increment), segment, datatype, to, (tag + i),
115                   comm);
116       }
117       smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
118     }
119
120     /* intermediate nodes relay (receive, reduce, then send) data */
121     else {
122       for (i = 0; i < pipe_length; i++) {
123         recv_request_array[i] = smpi_mpi_irecv((char *) tmp_buf + (i * increment), segment, datatype, from,
124                   (tag + i), comm);
125       }
126       for (i = 0; i < pipe_length; i++) {
127         smpi_mpi_wait(&recv_request_array[i], &status);
128         smpi_op_apply(op, tmp_buf + (i * increment), (char *)rbuf + (i * increment),
129                        &segment, &datatype);
130         send_request_array[i] = smpi_mpi_isend((char *) rbuf + (i * increment), segment, datatype, to,
131                   (tag + i), comm);
132       }
133       smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
134     }
135
136     free(send_request_array);
137     free(recv_request_array);
138     free(send_status_array);
139     free(recv_status_array);
140   }                             /* end pipeline */
141
142   /* when count is not divisible by block size, use default BCAST for the remainder */
143   if ((remainder != 0) && (count > segment)) {
144     XBT_WARN("MPI_reduce_NTSL use default MPI_reduce.");          
145     smpi_mpi_reduce((char *)buf + (pipe_length * increment),
146                (char *)rbuf + (pipe_length * increment), remainder, datatype, op, root,
147                comm);
148   }
149
150   free(tmp_buf);
151   return MPI_SUCCESS;
152 }