Logo AND Algorithmique Numérique Distribuée

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