Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add tesh files to test all new collectives
[simgrid.git] / src / smpi / colls / reduce-flat-tree.c
1 #include "colls.h"
2 //#include <star-reduction.c>
3
4 int
5 smpi_coll_tuned_reduce_flat_tree(void *sbuf, void *rbuf, int count,
6                                  MPI_Datatype dtype, MPI_Op op,
7                                  int root, MPI_Comm comm)
8 {
9   int i, tag = 4321;
10   int size;
11   int rank;
12   MPI_Aint extent;
13   char *origin = 0;
14   char *inbuf;
15   MPI_Status status;
16
17   MPI_Comm_rank(comm, &rank);
18   MPI_Comm_size(comm, &size);
19
20   /* If not root, send data to the root. */
21   MPI_Type_extent(dtype, &extent);
22
23   if (rank != root) {
24     MPI_Send(sbuf, count, dtype, root, tag, comm);
25     return 0;
26   }
27
28   /* Root receives and reduces messages.  Allocate buffer to receive
29      messages. */
30
31   if (size > 1)
32     origin = (char *) malloc(count * extent);
33
34
35   /* Initialize the receive buffer. */
36   if (rank == (size - 1))
37     MPI_Sendrecv(sbuf, count, dtype, rank, tag,
38                  rbuf, count, dtype, rank, tag, comm, &status);
39   else
40     MPI_Recv(rbuf, count, dtype, size - 1, tag, comm, &status);
41
42   /* Loop receiving and calling reduction function (C or Fortran). */
43
44   for (i = size - 2; i >= 0; --i) {
45     if (rank == i)
46       inbuf = sbuf;
47     else {
48       MPI_Recv(origin, count, dtype, i, tag, comm, &status);
49       inbuf = origin;
50     }
51
52     /* Call reduction function. */
53     star_reduction(op, inbuf, rbuf, &count, &dtype);
54
55   }
56
57   if (origin)
58     free(origin);
59
60   /* All done */
61   return 0;
62 }