Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
more leaks
[simgrid.git] / src / smpi / colls / reduce-binomial.c
1 #include "colls_private.h"
2
3 //#include <star-reduction.c>
4
5 int smpi_coll_tuned_reduce_binomial(void *sendbuf, void *recvbuf, int count,
6                                     MPI_Datatype datatype, MPI_Op op, int root,
7                                     MPI_Comm comm)
8 {
9   MPI_Status status;
10   int comm_size, rank;
11   int mask, relrank, source;
12   int dst;
13   int tag = COLL_TAG_REDUCE;
14   MPI_Aint extent;
15   void *tmp_buf;
16   MPI_Aint true_lb, true_extent;
17   if (count == 0)
18     return 0;
19   rank = smpi_comm_rank(comm);
20   comm_size = smpi_comm_size(comm);
21
22   extent = smpi_datatype_get_extent(datatype);
23
24   tmp_buf = (void *) xbt_malloc(count * extent);
25   int is_commutative = smpi_op_is_commute(op);
26   mask = 1;
27   
28   int lroot;
29   if (is_commutative) 
30         lroot   = root;
31   else
32         lroot   = 0;
33   relrank = (rank - lroot + comm_size) % comm_size;
34
35   smpi_datatype_extent(datatype, &true_lb, &true_extent);
36
37   /* adjust for potential negative lower bound in datatype */
38   tmp_buf = (void *)((char*)tmp_buf - true_lb);
39     
40   /* If I'm not the root, then my recvbuf may not be valid, therefore
41      I have to allocate a temporary one */
42   if (rank != root) {
43       recvbuf = (void *) malloc(count*(max(extent,true_extent)));
44       recvbuf = (void *)((char*)recvbuf - true_lb);
45   }
46    if ((rank != root) || (sendbuf != MPI_IN_PLACE)) {
47       smpi_datatype_copy(sendbuf, count, datatype, recvbuf,count, datatype);
48   }
49
50   while (mask < comm_size) {
51     /* Receive */
52     if ((mask & relrank) == 0) {
53       source = (relrank | mask);
54       if (source < comm_size) {
55         source = (source + lroot) % comm_size;
56         smpi_mpi_recv(tmp_buf, count, datatype, source, tag, comm, &status);
57         
58         if (is_commutative) {
59           smpi_op_apply(op, tmp_buf, recvbuf, &count, &datatype);
60         } else {
61           smpi_op_apply(op, recvbuf, tmp_buf, &count, &datatype);
62           smpi_datatype_copy(tmp_buf, count, datatype,recvbuf, count, datatype);
63         }
64       }
65     } else {
66       dst = ((relrank & (~mask)) + lroot) % comm_size;
67       smpi_mpi_send(recvbuf, count, datatype, dst, tag, comm);
68       break;
69     }
70     mask <<= 1;
71   }
72
73   if (!is_commutative && (root != 0)){
74     if (rank == 0){
75       smpi_mpi_send(recvbuf, count, datatype, root,tag, comm);
76     }else if (rank == root){
77       smpi_mpi_recv(recvbuf, count, datatype, 0, tag, comm, &status);
78     }
79   }
80
81   if (rank != root) {
82         xbt_free(recvbuf);
83   }
84   free(tmp_buf);
85
86   return 0;
87 }