Logo AND Algorithmique Numérique Distribuée

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