A
lgorithmique
N
umérique
D
istribuée
Public GIT Repository
projects
/
simgrid.git
/ blobdiff
commit
grep
author
committer
pickaxe
?
search:
re
summary
|
shortlog
|
log
|
commit
|
commitdiff
|
tree
raw
|
inline
| side by side
Fix smpi issue
[simgrid.git]
/
src
/
smpi
/
colls
/
allreduce-smp-binomial-pipeline.c
diff --git
a/src/smpi/colls/allreduce-smp-binomial-pipeline.c
b/src/smpi/colls/allreduce-smp-binomial-pipeline.c
index
b5efc1e
..
68c99f5
100644
(file)
--- a/
src/smpi/colls/allreduce-smp-binomial-pipeline.c
+++ b/
src/smpi/colls/allreduce-smp-binomial-pipeline.c
@@
-1,4
+1,4
@@
-#include "colls.h"
+#include "colls
_private
.h"
/* IMPLEMENTED BY PITCH PATARASUK
Non-topoloty-specific (however, number of cores/node need to be changed)
all-reduce operation designed for smp clusters
/* IMPLEMENTED BY PITCH PATARASUK
Non-topoloty-specific (however, number of cores/node need to be changed)
all-reduce operation designed for smp clusters
@@
-29,18
+29,6
@@
int allreduce_smp_binomial_pipeline_segment_size = 4096;
This code assume commutative and associative reduce operator (MPI_SUM, MPI_MAX, etc).
*/
This code assume commutative and associative reduce operator (MPI_SUM, MPI_MAX, etc).
*/
-#ifndef MPICH2
-extern void *MPIR_ToPointer();
-
-struct MPIR_OP {
- MPI_User_function *op;
- int commute, permanent;
-};
-
-#else
-extern MPI_User_function *MPIR_Op_table[];
-#endif
-
/*
This fucntion performs all-reduce operation as follow. ** in a pipeline fashion **
1) binomial_tree reduce inside each SMP node
/*
This fucntion performs all-reduce operation as follow. ** in a pipeline fashion **
1) binomial_tree reduce inside each SMP node
@@
-55,24
+43,16
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
{
int comm_size, rank;
void *tmp_buf;
{
int comm_size, rank;
void *tmp_buf;
- int tag =
50
;
+ int tag =
COLL_TAG_ALLREDUCE
;
int mask, src, dst;
MPI_Status status;
int num_core = NUM_CORE;
int mask, src, dst;
MPI_Status status;
int num_core = NUM_CORE;
- MPI_User_function *uop;
-#ifndef MPICH2
- struct MPIR_OP *op_ptr = MPIR_ToPointer(op);
- uop = (MPI_User_function *) op_ptr->op;
-#else
- uop = MPIR_Op_table[op % 16 - 1];
-#endif
-
- MPI_Comm_size(comm, &comm_size);
- MPI_Comm_rank(comm, &rank);
+ comm_size = smpi_comm_size(comm);
+ rank = smpi_comm_rank(comm);
MPI_Aint extent;
MPI_Aint extent;
-
MPI_Type_extent(dtype, &extent
);
- tmp_buf = (void *) malloc(count * extent);
+
extent = smpi_datatype_get_extent(dtype
);
+ tmp_buf = (void *)
xbt_
malloc(count * extent);
int intra_rank, inter_rank;
intra_rank = rank % num_core;
int intra_rank, inter_rank;
intra_rank = rank % num_core;
@@
-91,7
+71,7
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
int inter_comm_size = (comm_size + num_core - 1) / num_core;
/* copy input buffer to output buffer */
int inter_comm_size = (comm_size + num_core - 1) / num_core;
/* copy input buffer to output buffer */
-
MPI_S
endrecv(send_buf, count, dtype, rank, tag,
+
smpi_mpi_s
endrecv(send_buf, count, dtype, rank, tag,
recv_buf, count, dtype, rank, tag, comm, &status);
/* compute pipe length */
recv_buf, count, dtype, rank, tag, comm, &status);
/* compute pipe length */
@@
-110,13
+90,13
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
src = (inter_rank * num_core) + (intra_rank | mask);
if (src < comm_size) {
recv_offset = phase * pcount * extent;
src = (inter_rank * num_core) + (intra_rank | mask);
if (src < comm_size) {
recv_offset = phase * pcount * extent;
-
MPI_R
ecv(tmp_buf, pcount, dtype, src, tag, comm, &status);
-
(*uop) (
tmp_buf, (char *)recv_buf + recv_offset, &pcount, &dtype);
+
smpi_mpi_r
ecv(tmp_buf, pcount, dtype, src, tag, comm, &status);
+
smpi_op_apply(op,
tmp_buf, (char *)recv_buf + recv_offset, &pcount, &dtype);
}
} else {
send_offset = phase * pcount * extent;
dst = (inter_rank * num_core) + (intra_rank & (~mask));
}
} else {
send_offset = phase * pcount * extent;
dst = (inter_rank * num_core) + (intra_rank & (~mask));
-
MPI_S
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
+
smpi_mpi_s
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
break;
}
mask <<= 1;
break;
}
mask <<= 1;
@@
-134,13
+114,13
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
src = (inter_rank | mask) * num_core;
if (src < comm_size) {
recv_offset = (phase - 1) * pcount * extent;
src = (inter_rank | mask) * num_core;
if (src < comm_size) {
recv_offset = (phase - 1) * pcount * extent;
-
MPI_R
ecv(tmp_buf, pcount, dtype, src, tag, comm, &status);
-
(*uop) (
tmp_buf, (char *)recv_buf + recv_offset, &pcount, &dtype);
+
smpi_mpi_r
ecv(tmp_buf, pcount, dtype, src, tag, comm, &status);
+
smpi_op_apply(op,
tmp_buf, (char *)recv_buf + recv_offset, &pcount, &dtype);
}
} else {
dst = (inter_rank & (~mask)) * num_core;
send_offset = (phase - 1) * pcount * extent;
}
} else {
dst = (inter_rank & (~mask)) * num_core;
send_offset = (phase - 1) * pcount * extent;
-
MPI_S
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
+
smpi_mpi_s
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
break;
}
mask <<= 1;
break;
}
mask <<= 1;
@@
-157,7
+137,7
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
if (inter_rank & mask) {
src = (inter_rank - mask) * num_core;
recv_offset = (phase - 2) * pcount * extent;
if (inter_rank & mask) {
src = (inter_rank - mask) * num_core;
recv_offset = (phase - 2) * pcount * extent;
-
MPI_R
ecv((char *)recv_buf + recv_offset, pcount, dtype, src, tag, comm,
+
smpi_mpi_r
ecv((char *)recv_buf + recv_offset, pcount, dtype, src, tag, comm,
&status);
break;
}
&status);
break;
}
@@
-171,7
+151,7
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
if (dst < comm_size) {
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
send_offset = (phase - 2) * pcount * extent;
if (dst < comm_size) {
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
send_offset = (phase - 2) * pcount * extent;
-
MPI_S
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
+
smpi_mpi_s
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
}
}
mask >>= 1;
}
}
mask >>= 1;
@@
-190,7
+170,7
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
if (intra_rank & mask) {
src = (inter_rank * num_core) + (intra_rank - mask);
recv_offset = (phase - 3) * pcount * extent;
if (intra_rank & mask) {
src = (inter_rank * num_core) + (intra_rank - mask);
recv_offset = (phase - 3) * pcount * extent;
-
MPI_R
ecv((char *)recv_buf + recv_offset, pcount, dtype, src, tag, comm,
+
smpi_mpi_r
ecv((char *)recv_buf + recv_offset, pcount, dtype, src, tag, comm,
&status);
break;
}
&status);
break;
}
@@
-202,7
+182,7
@@
int smpi_coll_tuned_allreduce_smp_binomial_pipeline(void *send_buf,
dst = (inter_rank * num_core) + (intra_rank + mask);
if (dst < comm_size) {
send_offset = (phase - 3) * pcount * extent;
dst = (inter_rank * num_core) + (intra_rank + mask);
if (dst < comm_size) {
send_offset = (phase - 3) * pcount * extent;
-
MPI_S
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
+
smpi_mpi_s
end((char *)recv_buf + send_offset, pcount, dtype, dst, tag, comm);
}
mask >>= 1;
}
}
mask >>= 1;
}