Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
added tesh tests for DVFS
[simgrid.git] / src / smpi / colls / bcast-SMP-binomial.c
1 #include "colls_private.h"
2 #ifndef NUM_CORE
3 #define NUM_CORE 8
4 #endif
5
6 int smpi_coll_tuned_bcast_SMP_binomial(void *buf, int count,
7                                        MPI_Datatype datatype, int root,
8                                        MPI_Comm comm)
9 {
10   int mask = 1;
11   int size;
12   int rank;
13   MPI_Status status;
14   int tag = 50;
15
16   size = smpi_comm_size(comm);
17   rank = smpi_comm_rank(comm);
18
19   int to_intra, to_inter;
20   int from_intra, from_inter;
21   int inter_rank = rank / NUM_CORE;
22   int inter_size = (size - 1) / NUM_CORE + 1;
23   int intra_rank = rank % NUM_CORE;
24   int intra_size = NUM_CORE;
25   if (((rank / NUM_CORE) * NUM_CORE) == ((size / NUM_CORE) * NUM_CORE))
26     intra_size = size - (rank / NUM_CORE) * NUM_CORE;
27
28   // if root is not zero send to rank zero first
29   if (root != 0) {
30     if (rank == root)
31       smpi_mpi_send(buf, count, datatype, 0, tag, comm);
32     else if (rank == 0)
33       smpi_mpi_recv(buf, count, datatype, root, tag, comm, &status);
34   }
35   //FIRST STEP node 0 send to every root-of-each-SMP with binomial tree
36
37   //printf("node %d inter_rank = %d, inter_size = %d\n",rank,inter_rank, inter_size);
38
39   if (intra_rank == 0) {
40     mask = 1;
41     while (mask < inter_size) {
42       if (inter_rank & mask) {
43         from_inter = (inter_rank - mask) * NUM_CORE;
44         //printf("Node %d recv from node %d when mask is %d\n", rank, from_inter, mask);
45         smpi_mpi_recv(buf, count, datatype, from_inter, tag, comm, &status);
46         break;
47       }
48       mask <<= 1;
49     }
50
51     mask >>= 1;
52     //printf("My rank = %d my mask = %d\n", rank,mask);
53
54     while (mask > 0) {
55       if (inter_rank < inter_size) {
56         to_inter = (inter_rank + mask) * NUM_CORE;
57         if (to_inter < size) {
58           //printf("Node %d send to node %d when mask is %d\n", rank, to_inter, mask);
59           smpi_mpi_send(buf, count, datatype, to_inter, tag, comm);
60         }
61       }
62       mask >>= 1;
63     }
64   }
65   // SECOND STEP every root-of-each-SMP send to all children with binomial tree
66   // base is a rank of root-of-each-SMP
67   int base = (rank / NUM_CORE) * NUM_CORE;
68   mask = 1;
69   while (mask < intra_size) {
70     if (intra_rank & mask) {
71       from_intra = base + (intra_rank - mask);
72       //printf("Node %d recv from node %d when mask is %d\n", rank, from_inter, mask);
73       smpi_mpi_recv(buf, count, datatype, from_intra, tag, comm, &status);
74       break;
75     }
76     mask <<= 1;
77   }
78
79   mask >>= 1;
80
81   //printf("My rank = %d my mask = %d\n", rank,mask);
82
83   while (mask > 0) {
84     if (intra_rank < intra_size) {
85       to_intra = base + (intra_rank + mask);
86       if (to_intra < size) {
87         //printf("Node %d send to node %d when mask is %d\n", rank, to_inter, mask);
88         smpi_mpi_send(buf, count, datatype, to_intra, tag, comm);
89       }
90     }
91     mask >>= 1;
92   }
93
94   return MPI_SUCCESS;
95 }