Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
3b5bf4b71ed09ed86ad0018fd760a1b834eec24a
[simgrid.git] / src / smpi / colls / bcast / bcast-flattree-pipeline.cpp
1 /* Copyright (c) 2013-2017. 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 int flattree_segment_in_byte = 8192;
10 namespace simgrid{
11 namespace smpi{
12 int
13 Coll_bcast_flattree_pipeline::bcast(void *buff, int count,
14                                         MPI_Datatype data_type, int root,
15                                         MPI_Comm comm)
16 {
17   int i, j, rank, num_procs;
18   int tag = COLL_TAG_BCAST;
19
20   MPI_Aint extent;
21   extent = data_type->get_extent();
22
23   int segment = flattree_segment_in_byte / extent;
24   segment =  segment == 0 ? 1 :segment; 
25   int pipe_length = count / segment;
26   int increment = segment * extent;
27   if (pipe_length==0) {
28     XBT_WARN("MPI_bcast_flattree_pipeline use default MPI_bcast_flattree.");
29     return Coll_bcast_flattree::bcast(buff, count, data_type, root, comm);
30   }
31   rank = comm->rank();
32   num_procs = comm->size();
33
34   MPI_Request *request_array;
35   MPI_Status *status_array;
36
37   request_array = (MPI_Request *) xbt_malloc(pipe_length * sizeof(MPI_Request));
38   status_array = (MPI_Status *) xbt_malloc(pipe_length * sizeof(MPI_Status));
39
40   if (rank != root) {
41     for (i = 0; i < pipe_length; i++) {
42       request_array[i] = Request::irecv((char *)buff + (i * increment), segment, data_type, root, tag, comm);
43     }
44     Request::waitall(pipe_length, request_array, status_array);
45   }
46
47   else {
48     // Root sends data to all others
49     for (j = 0; j < num_procs; j++) {
50       if (j == rank)
51         continue;
52       else {
53         for (i = 0; i < pipe_length; i++) {
54           Request::send((char *)buff + (i * increment), segment, data_type, j, tag, comm);
55         }
56       }
57     }
58
59   }
60
61   free(request_array);
62   free(status_array);
63   return MPI_SUCCESS;
64 }
65
66 }
67 }