Logo AND Algorithmique Numérique Distribuée

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