Logo AND Algorithmique Numérique Distribuée

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