Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines.
[simgrid.git] / src / smpi / colls / reduce / reduce-NTSL.cpp
1 /* Copyright (c) 2013-2021. 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 reduce__NTSL(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   unsigned char* tmp_buf = smpi_get_tmp_sendbuffer(count * extent);
64
65   Request::sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank,
66                tag, comm, &status);
67
68   /* when a message is smaller than a block size => no pipeline */
69   if (count <= segment) {
70     if (rank == root) {
71       Request::recv(tmp_buf, count, datatype, from, tag, comm, &status);
72       if(op!=MPI_OP_NULL) op->apply( tmp_buf, rbuf, &count, datatype);
73     } else if (rank == ((root - 1 + size) % size)) {
74       Request::send(rbuf, count, datatype, to, tag, comm);
75     } else {
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       Request::send(rbuf, count, datatype, to, tag, comm);
79     }
80     smpi_free_tmp_buffer(tmp_buf);
81     return MPI_SUCCESS;
82   }
83
84   /* pipeline */
85   else {
86     auto* send_request_array = new MPI_Request[size + pipe_length];
87     auto* recv_request_array = new MPI_Request[size + pipe_length];
88     auto* send_status_array  = new MPI_Status[size + pipe_length];
89     auto* recv_status_array  = new MPI_Status[size + pipe_length];
90
91     /* root recv data */
92     if (rank == root) {
93       for (i = 0; i < pipe_length; i++) {
94         recv_request_array[i] = Request::irecv(tmp_buf + (i * increment), segment, datatype, from, (tag + i), comm);
95       }
96       for (i = 0; i < pipe_length; i++) {
97         Request::wait(&recv_request_array[i], &status);
98         if(op!=MPI_OP_NULL) op->apply( tmp_buf + (i * increment), (char *)rbuf + (i * increment),
99                        &segment, datatype);
100       }
101     }
102
103     /* last node only sends data */
104     else if (rank == ((root - 1 + size) % size)) {
105       for (i = 0; i < pipe_length; i++) {
106         send_request_array[i] = Request::isend((char *)rbuf + (i * increment), segment, datatype, to, (tag + i),
107                   comm);
108       }
109       Request::waitall((pipe_length), send_request_array, send_status_array);
110     }
111
112     /* intermediate nodes relay (receive, reduce, then send) data */
113     else {
114       for (i = 0; i < pipe_length; i++) {
115         recv_request_array[i] = Request::irecv(tmp_buf + (i * increment), segment, datatype, from, (tag + i), comm);
116       }
117       for (i = 0; i < pipe_length; i++) {
118         Request::wait(&recv_request_array[i], &status);
119         if(op!=MPI_OP_NULL) op->apply( tmp_buf + (i * increment), (char *)rbuf + (i * increment),
120                        &segment, datatype);
121         send_request_array[i] = Request::isend((char *) rbuf + (i * increment), segment, datatype, to,
122                   (tag + i), comm);
123       }
124       Request::waitall((pipe_length), send_request_array, send_status_array);
125     }
126
127     delete[] send_request_array;
128     delete[] recv_request_array;
129     delete[] send_status_array;
130     delete[] recv_status_array;
131   }                             /* end pipeline */
132
133   /* when count is not divisible by block size, use default BCAST for the remainder */
134   if ((remainder != 0) && (count > segment)) {
135     XBT_WARN("MPI_reduce_NTSL use default MPI_reduce.");
136     reduce__default((char *)buf + (pipe_length * increment),
137                     (char *)rbuf + (pipe_length * increment), remainder, datatype, op, root,
138                     comm);
139   }
140
141   smpi_free_tmp_buffer(tmp_buf);
142   return MPI_SUCCESS;
143 }
144 }
145 }