Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Fix a doc error about actors (Tutorial_algorithms)
[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   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;
32   int from = (rank + size - 1) % size;
33
34   /* segment is segment size in number of elements (not bytes) */
35   int segment = extent == 0 ? 1 : (bcast_NTSL_segment_size_in_byte / extent);
36   segment =  segment == 0 ? 1 :segment;
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   if (root != 0) {
51     if (rank == root) {
52       Request::send(buf, count, datatype, 0, tag, comm);
53     } else if (rank == 0) {
54       Request::recv(buf, count, datatype, root, tag, comm, &status);
55     }
56   }
57
58   /* when a message is smaller than a block size => no pipeline */
59   if (count <= segment) {
60     if (rank == 0) {
61       Request::send(buf, count, datatype, to, tag, comm);
62     } else if (rank == (size - 1)) {
63       request = Request::irecv(buf, count, datatype, from, tag, comm);
64       Request::wait(&request, &status);
65     } else {
66       request = Request::irecv(buf, count, datatype, from, tag, comm);
67       Request::wait(&request, &status);
68       Request::send(buf, count, datatype, to, tag, comm);
69     }
70     return MPI_SUCCESS;
71   }
72
73   /* pipeline bcast */
74   else {
75     MPI_Request* send_request_array = new MPI_Request[size + pipe_length];
76     MPI_Request* recv_request_array = new MPI_Request[size + pipe_length];
77     MPI_Status* send_status_array   = new MPI_Status[size + pipe_length];
78     MPI_Status* recv_status_array   = new MPI_Status[size + pipe_length];
79
80     /* root send data */
81     if (rank == 0) {
82       for (i = 0; i < pipe_length; i++) {
83         send_request_array[i] = Request::isend((char *) buf + (i * increment), segment, datatype, to,
84                   (tag + i), comm);
85       }
86       Request::waitall((pipe_length), send_request_array, send_status_array);
87     }
88
89     /* last node only receive data */
90     else if (rank == (size - 1)) {
91       for (i = 0; i < pipe_length; i++) {
92         recv_request_array[i] = Request::irecv((char *) buf + (i * increment), segment, datatype, from,
93                   (tag + i), comm);
94       }
95       Request::waitall((pipe_length), recv_request_array, recv_status_array);
96     }
97
98     /* intermediate nodes relay (receive, then send) data */
99     else {
100       for (i = 0; i < pipe_length; i++) {
101         recv_request_array[i] = Request::irecv((char *) buf + (i * increment), segment, datatype, from,
102                   (tag + i), comm);
103       }
104       for (i = 0; i < pipe_length; i++) {
105         Request::wait(&recv_request_array[i], &status);
106         send_request_array[i] = Request::isend((char *) buf + (i * increment), segment, datatype, to,
107                   (tag + i), comm);
108       }
109       Request::waitall((pipe_length), send_request_array, send_status_array);
110     }
111
112     delete[] send_request_array;
113     delete[] recv_request_array;
114     delete[] send_status_array;
115     delete[] recv_status_array;
116   }                             /* end pipeline */
117
118   /* when count is not divisible by block size, use default BCAST for the remainder */
119   if ((remainder != 0) && (count > segment)) {
120     XBT_WARN("MPI_bcast_arrival_NTSL use default MPI_bcast.");
121     Colls::bcast((char *) buf + (pipe_length * increment), remainder, datatype,
122               root, comm);
123   }
124
125   return MPI_SUCCESS;
126 }
127
128 }
129 }