Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
223996096b7450041c1d5f6b28ecf2e35c27147e
[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 = COLL_TAG_BCAST;
15
16   size = smpi_comm_size(comm);
17   rank = smpi_comm_rank(comm);
18
19   if(size%NUM_CORE)
20     THROWF(arg_error,0, "bcast SMP binomial can't be used with non multiple of NUM_CORE=%d number of processes ! ",NUM_CORE);
21
22   int to_intra, to_inter;
23   int from_intra, from_inter;
24   int inter_rank = rank / NUM_CORE;
25   int inter_size = (size - 1) / NUM_CORE + 1;
26   int intra_rank = rank % NUM_CORE;
27   int intra_size = NUM_CORE;
28   if (((rank / NUM_CORE) * NUM_CORE) == ((size / NUM_CORE) * NUM_CORE))
29     intra_size = size - (rank / NUM_CORE) * NUM_CORE;
30
31   // if root is not zero send to rank zero first
32   if (root != 0) {
33     if (rank == root)
34       smpi_mpi_send(buf, count, datatype, 0, tag, comm);
35     else if (rank == 0)
36       smpi_mpi_recv(buf, count, datatype, root, tag, comm, &status);
37   }
38   //FIRST STEP node 0 send to every root-of-each-SMP with binomial tree
39
40   //printf("node %d inter_rank = %d, inter_size = %d\n",rank,inter_rank, inter_size);
41
42   if (intra_rank == 0) {
43     mask = 1;
44     while (mask < inter_size) {
45       if (inter_rank & mask) {
46         from_inter = (inter_rank - mask) * NUM_CORE;
47         //printf("Node %d recv from node %d when mask is %d\n", rank, from_inter, mask);
48         smpi_mpi_recv(buf, count, datatype, from_inter, tag, comm, &status);
49         break;
50       }
51       mask <<= 1;
52     }
53
54     mask >>= 1;
55     //printf("My rank = %d my mask = %d\n", rank,mask);
56
57     while (mask > 0) {
58       if (inter_rank < inter_size) {
59         to_inter = (inter_rank + mask) * NUM_CORE;
60         if (to_inter < size) {
61           //printf("Node %d send to node %d when mask is %d\n", rank, to_inter, mask);
62           smpi_mpi_send(buf, count, datatype, to_inter, tag, comm);
63         }
64       }
65       mask >>= 1;
66     }
67   }
68   // SECOND STEP every root-of-each-SMP send to all children with binomial tree
69   // base is a rank of root-of-each-SMP
70   int base = (rank / NUM_CORE) * NUM_CORE;
71   mask = 1;
72   while (mask < intra_size) {
73     if (intra_rank & mask) {
74       from_intra = base + (intra_rank - mask);
75       //printf("Node %d recv from node %d when mask is %d\n", rank, from_inter, mask);
76       smpi_mpi_recv(buf, count, datatype, from_intra, tag, comm, &status);
77       break;
78     }
79     mask <<= 1;
80   }
81
82   mask >>= 1;
83
84   //printf("My rank = %d my mask = %d\n", rank,mask);
85
86   while (mask > 0) {
87     if (intra_rank < intra_size) {
88       to_intra = base + (intra_rank + mask);
89       if (to_intra < size) {
90         //printf("Node %d send to node %d when mask is %d\n", rank, to_inter, mask);
91         smpi_mpi_send(buf, count, datatype, to_intra, tag, comm);
92       }
93     }
94     mask >>= 1;
95   }
96
97   return MPI_SUCCESS;
98 }