Logo AND Algorithmique Numérique Distribuée

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