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
leaks --
[simgrid.git]
/
src
/
smpi
/
colls
/
allreduce-smp-rsag-lr.c
diff --git
a/src/smpi/colls/allreduce-smp-rsag-lr.c
b/src/smpi/colls/allreduce-smp-rsag-lr.c
index
713ae31
..
6a928bc
100644
(file)
--- a/
src/smpi/colls/allreduce-smp-rsag-lr.c
+++ b/
src/smpi/colls/allreduce-smp-rsag-lr.c
@@
-1,11
+1,11
@@
-#include "colls.h"
-//#include <star-reduction.c>
+/* Copyright (c) 2013-2014. The SimGrid Team.
+ * All rights reserved. */
+
+/* This program is free software; you can redistribute it and/or modify it
+ * under the terms of the license (GNU LGPL) which comes with this package. */
-/* change number of core per smp-node
- we assume that number of core per process will be the same for all implementations */
-#ifndef NUM_CORE
-#define NUM_CORE 8
-#endif
+#include "colls_private.h"
+//#include <star-reduction.c>
/*
This fucntion performs all-reduce operation as follow.
/*
This fucntion performs all-reduce operation as follow.
@@
-20,10
+20,16
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_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 mask, src, dst;
MPI_Status status;
- int num_core = NUM_CORE;
+ if(smpi_comm_get_leaders_comm(comm)==MPI_COMM_NULL){
+ smpi_comm_init_smp(comm);
+ }
+ int num_core=1;
+ if (smpi_comm_is_uniform(comm)){
+ num_core = smpi_comm_size(smpi_comm_get_intra_comm(comm));
+ }
/*
#ifdef MPICH2_REDUCTION
MPI_User_function * uop = MPIR_Op_table[op % 16 - 1];
/*
#ifdef MPICH2_REDUCTION
MPI_User_function * uop = MPIR_Op_table[op % 16 - 1];
@@
-34,11
+40,11
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
uop = op_ptr->op;
#endif
*/
uop = op_ptr->op;
#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 *)
smpi_get_tmp_sendbuffer
(count * extent);
int intra_rank, inter_rank;
intra_rank = rank % num_core;
int intra_rank, inter_rank;
intra_rank = rank % num_core;
@@
-54,7
+60,7
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
}
}
-
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);
recv_buf, count, dtype, rank, tag, comm, &status);
@@
-65,14
+71,14
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
src = (inter_rank * num_core) + (intra_rank | mask);
// if (src < ((inter_rank + 1) * num_core)) {
if (src < comm_size) {
src = (inter_rank * num_core) + (intra_rank | mask);
// if (src < ((inter_rank + 1) * num_core)) {
if (src < comm_size) {
-
MPI_R
ecv(tmp_buf, count, dtype, src, tag, comm, &status);
- s
tar_reduction
(op, tmp_buf, recv_buf, &count, &dtype);
+
smpi_mpi_r
ecv(tmp_buf, count, dtype, src, tag, comm, &status);
+ s
mpi_op_apply
(op, tmp_buf, recv_buf, &count, &dtype);
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
}
} else {
dst = (inter_rank * num_core) + (intra_rank & (~mask));
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
}
} else {
dst = (inter_rank * num_core) + (intra_rank & (~mask));
-
MPI_S
end(recv_buf, count, dtype, dst, tag, comm);
+
smpi_mpi_s
end(recv_buf, count, dtype, dst, tag, comm);
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
break;
}
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
break;
}
@@
-119,12
+125,12
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
else
recv_count = curr_size + curr_remainder;
else
recv_count = curr_size + curr_remainder;
-
MPI_S
endrecv((char *) recv_buf + send_offset, send_count, dtype, to,
+
smpi_mpi_s
endrecv((char *) recv_buf + send_offset, send_count, dtype, to,
tag + i, tmp_buf, recv_count, dtype, from, tag + i, comm,
&status);
// result is in rbuf
tag + i, tmp_buf, recv_count, dtype, from, tag + i, comm,
&status);
// result is in rbuf
- s
tar_reduction
(op, tmp_buf, (char *) recv_buf + recv_offset, &recv_count,
+ s
mpi_op_apply
(op, tmp_buf, (char *) recv_buf + recv_offset, &recv_count,
&dtype);
}
&dtype);
}
@@
-149,7
+155,7
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
else
recv_count = curr_size + curr_remainder;
else
recv_count = curr_size + curr_remainder;
-
MPI_S
endrecv((char *) recv_buf + send_offset, send_count, dtype, to,
+
smpi_mpi_s
endrecv((char *) recv_buf + send_offset, send_count, dtype, to,
tag + i, (char *) recv_buf + recv_offset, recv_count, dtype,
from, tag + i, comm, &status);
tag + i, (char *) recv_buf + recv_offset, recv_count, dtype,
from, tag + i, comm, &status);
@@
-169,14
+175,14
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
if ((mask & inter_rank) == 0) {
src = (inter_rank | mask) * num_core;
if (src < comm_size) {
if ((mask & inter_rank) == 0) {
src = (inter_rank | mask) * num_core;
if (src < comm_size) {
-
MPI_R
ecv(tmp_buf, count, dtype, src, tag, comm, &status);
+
smpi_mpi_r
ecv(tmp_buf, count, dtype, src, tag, comm, &status);
(* uop) (tmp_buf, recv_buf, &count, &dtype);
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
}
}
else {
dst = (inter_rank & (~mask)) * num_core;
(* uop) (tmp_buf, recv_buf, &count, &dtype);
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
}
}
else {
dst = (inter_rank & (~mask)) * num_core;
-
MPI_S
end(recv_buf, count, dtype, dst, tag, comm);
+
smpi_mpi_s
end(recv_buf, count, dtype, dst, tag, comm);
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
break;
}
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
break;
}
@@
-195,7
+201,7
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
if (inter_rank & mask) {
src = (inter_rank - mask) * num_core;
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
if (inter_rank & mask) {
src = (inter_rank - mask) * num_core;
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
-
MPI_R
ecv(recv_buf, count, dtype, src, tag, comm, &status);
+
smpi_mpi_r
ecv(recv_buf, count, dtype, src, tag, comm, &status);
break;
}
mask <<= 1;
break;
}
mask <<= 1;
@@
-209,7
+215,7
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
dst = (inter_rank + mask) * num_core;
if (dst < comm_size) {
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
dst = (inter_rank + mask) * num_core;
if (dst < comm_size) {
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
-
MPI_S
end(recv_buf, count, dtype, dst, tag, comm);
+
smpi_mpi_s
end(recv_buf, count, dtype, dst, tag, comm);
}
}
mask >>= 1;
}
}
mask >>= 1;
@@
-230,7
+236,7
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
if (intra_rank & mask) {
src = (inter_rank * num_core) + (intra_rank - mask);
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
if (intra_rank & mask) {
src = (inter_rank * num_core) + (intra_rank - mask);
//printf("Node %d recv from node %d when mask is %d\n", rank, src, mask);
-
MPI_R
ecv(recv_buf, count, dtype, src, tag, comm, &status);
+
smpi_mpi_r
ecv(recv_buf, count, dtype, src, tag, comm, &status);
break;
}
mask <<= 1;
break;
}
mask <<= 1;
@@
-243,12
+249,12
@@
int smpi_coll_tuned_allreduce_smp_rsag_lr(void *send_buf, void *recv_buf,
dst = (inter_rank * num_core) + (intra_rank + mask);
if (dst < comm_size) {
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
dst = (inter_rank * num_core) + (intra_rank + mask);
if (dst < comm_size) {
//printf("Node %d send to node %d when mask is %d\n", rank, dst, mask);
-
MPI_S
end(recv_buf, count, dtype, dst, tag, comm);
+
smpi_mpi_s
end(recv_buf, count, dtype, dst, tag, comm);
}
mask >>= 1;
}
}
mask >>= 1;
}
-
free
(tmp_buf);
+
smpi_free_tmp_buffer
(tmp_buf);
return MPI_SUCCESS;
}
return MPI_SUCCESS;
}