Logo AND Algorithmique Numérique Distribuée

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