Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
9fdb21d651a50ce9637ec16a31f4b898b48bc637
[simgrid.git] / src / smpi / colls / bcast / bcast-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
9 static int bcast_NTSL_segment_size_in_byte = 8192;
10
11 /* Non-topology-specific pipelined linear-bcast function
12    0->1, 1->2 ,2->3, ....., ->last node : in a pipeline fashion
13 */
14 namespace simgrid::smpi {
15 int bcast__NTSL(void *buf, int count, MPI_Datatype datatype,
16                 int root, MPI_Comm comm)
17 {
18   int tag = COLL_TAG_BCAST;
19   MPI_Status status;
20   MPI_Request request;
21   int rank, size;
22   int i;
23   MPI_Aint extent;
24   extent = datatype->get_extent();
25
26   rank = comm->rank();
27   size = comm->size();
28
29   /* source node and destination nodes (same through out the functions) */
30   int to = (rank + 1) % size;
31   int from = (rank + size - 1) % size;
32
33   /* segment is segment size in number of elements (not bytes) */
34   int segment = extent == 0 ? 1 : (bcast_NTSL_segment_size_in_byte / extent);
35   segment =  segment == 0 ? 1 :segment;
36   /* pipeline length */
37   int pipe_length = count / segment;
38
39   /* use for buffer offset for sending and receiving data = segment size in byte */
40   int increment = segment * extent;
41
42   /* if the input size is not divisible by segment size =>
43      the small remainder will be done with native implementation */
44   int remainder = count % segment;
45
46   /* if root is not zero send to rank zero first
47      this can be modified to make it faster by using logical src, dst.
48    */
49   if (root != 0) {
50     if (rank == root) {
51       Request::send(buf, count, datatype, 0, tag, comm);
52     } else if (rank == 0) {
53       Request::recv(buf, count, datatype, root, tag, comm, &status);
54     }
55   }
56
57   /* when a message is smaller than a block size => no pipeline */
58   if (count <= segment) {
59     if (rank == 0) {
60       Request::send(buf, count, datatype, to, tag, comm);
61     } else if (rank == (size - 1)) {
62       request = Request::irecv(buf, count, datatype, from, tag, comm);
63       Request::wait(&request, &status);
64     } else {
65       request = Request::irecv(buf, count, datatype, from, tag, comm);
66       Request::wait(&request, &status);
67       Request::send(buf, count, datatype, to, tag, comm);
68     }
69     return MPI_SUCCESS;
70   }
71
72   /* pipeline bcast */
73   else {
74     auto* send_request_array = new MPI_Request[size + pipe_length];
75     auto* recv_request_array = new MPI_Request[size + pipe_length];
76     auto* send_status_array  = new MPI_Status[size + pipe_length];
77     auto* recv_status_array  = new MPI_Status[size + pipe_length];
78
79     /* root send data */
80     if (rank == 0) {
81       for (i = 0; i < pipe_length; i++) {
82         send_request_array[i] = Request::isend((char *) buf + (i * increment), segment, datatype, to,
83                   (tag + i), comm);
84       }
85       Request::waitall((pipe_length), send_request_array, send_status_array);
86     }
87
88     /* last node only receive data */
89     else if (rank == (size - 1)) {
90       for (i = 0; i < pipe_length; i++) {
91         recv_request_array[i] = Request::irecv((char *) buf + (i * increment), segment, datatype, from,
92                   (tag + i), comm);
93       }
94       Request::waitall((pipe_length), recv_request_array, recv_status_array);
95     }
96
97     /* intermediate nodes relay (receive, then send) data */
98     else {
99       for (i = 0; i < pipe_length; i++) {
100         recv_request_array[i] = Request::irecv((char *) buf + (i * increment), segment, datatype, from,
101                   (tag + i), comm);
102       }
103       for (i = 0; i < pipe_length; i++) {
104         Request::wait(&recv_request_array[i], &status);
105         send_request_array[i] = Request::isend((char *) buf + (i * increment), segment, datatype, to,
106                   (tag + i), comm);
107       }
108       Request::waitall((pipe_length), send_request_array, send_status_array);
109     }
110
111     delete[] send_request_array;
112     delete[] recv_request_array;
113     delete[] send_status_array;
114     delete[] recv_status_array;
115   }                             /* end pipeline */
116
117   if ((remainder != 0) && (count > segment)) {
118     XBT_INFO("MPI_bcast_arrival_NTSL: count is not divisible by block size, use default MPI_bcast for remainder.");
119     colls::bcast((char*)buf + (pipe_length * increment), remainder, datatype, root, comm);
120   }
121
122   return MPI_SUCCESS;
123 }
124
125 } // namespace simgrid::smpi