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 MSG_storage_get_properties()
[simgrid.git]
/
src
/
smpi
/
colls
/
bcast-arrival-nb.c
diff --git
a/src/smpi/colls/bcast-arrival-nb.c
b/src/smpi/colls/bcast-arrival-nb.c
index
9ff27b4
..
84a06ae
100644
(file)
--- a/
src/smpi/colls/bcast-arrival-nb.c
+++ b/
src/smpi/colls/bcast-arrival-nb.c
@@
-10,7
+10,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
MPI_Datatype datatype, int root,
MPI_Comm comm)
{
MPI_Datatype datatype, int root,
MPI_Comm comm)
{
- int tag =
50
;
+ int tag =
-COLL_TAG_BCAST
;
MPI_Status status;
MPI_Request request;
MPI_Request *send_request_array;
MPI_Status status;
MPI_Request request;
MPI_Request *send_request_array;
@@
-27,7
+27,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
int header_index;
int flag_array[MAX_NODE];
int already_sent[MAX_NODE];
int header_index;
int flag_array[MAX_NODE];
int already_sent[MAX_NODE];
-
+ int to_clean[MAX_NODE];
int header_buf[HEADER_SIZE];
char temp_buf[MAX_NODE];
int header_buf[HEADER_SIZE];
char temp_buf[MAX_NODE];
@@
-39,8
+39,8
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
- rank = smpi_comm_rank(
MPI_COMM_WORLD
);
- size = smpi_comm_size(
MPI_COMM_WORLD
);
+ rank = smpi_comm_rank(
comm
);
+ size = smpi_comm_size(
comm
);
/* segment is segment size in number of elements (not bytes) */
/* segment is segment size in number of elements (not bytes) */
@@
-70,6
+70,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
/* value == 0 means root has not send data (or header) to the node yet */
for (i = 0; i < MAX_NODE; i++) {
already_sent[i] = 0;
/* value == 0 means root has not send data (or header) to the node yet */
for (i = 0; i < MAX_NODE; i++) {
already_sent[i] = 0;
+ to_clean[i]=0;
}
// printf("YYY\n");
}
// printf("YYY\n");
@@
-83,7
+84,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
// for (j=0;j<1000;j++) {
for (i = 1; i < size; i++) {
if (already_sent[i] == 0)
// for (j=0;j<1000;j++) {
for (i = 1; i < size; i++) {
if (already_sent[i] == 0)
- smpi_mpi_iprobe(i, MPI_ANY_TAG,
MPI_COMM_WORLD
, &flag_array[i],
+ smpi_mpi_iprobe(i, MPI_ANY_TAG,
comm
, &flag_array[i],
MPI_STATUSES_IGNORE);
}
//}
MPI_STATUSES_IGNORE);
}
//}
@@
-94,7
+95,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
/* message arrive */
if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
/* message arrive */
if ((flag_array[i] == 1) && (already_sent[i] == 0)) {
- smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag,
MPI_COMM_WORLD
, &status);
+ smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag,
comm
, &status);
header_buf[header_index] = i;
header_index++;
sent_count++;
header_buf[header_index] = i;
header_index++;
sent_count++;
@@
-123,6
+124,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
smpi_mpi_send(buf, count, datatype, i, tag, comm);
already_sent[i] = 1;
sent_count++;
smpi_mpi_send(buf, count, datatype, i, tag, comm);
already_sent[i] = 1;
sent_count++;
+ to_clean[i]=0;
break;
}
}
break;
}
}
@@
-130,6
+132,9
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
} /* while loop */
} /* while loop */
+
+ for(i=0; i<size; i++)
+ if(to_clean[i]!=0)smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, comm, &status);
}
/* non-root */
}
/* non-root */
@@
-183,11
+188,11
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
for (k = 0; k < 3; k++) {
for (i = 1; i < size; i++) {
if ((already_sent[i] == 0) && (will_send[i] == 0)) {
for (k = 0; k < 3; k++) {
for (i = 1; i < size; i++) {
if ((already_sent[i] == 0) && (will_send[i] == 0)) {
- smpi_mpi_iprobe(i, MPI_ANY_TAG,
MPI_COMM_WORLD
, &flag_array[i],
+ smpi_mpi_iprobe(i, MPI_ANY_TAG,
comm
, &flag_array[i],
&temp_status_array[i]);
if (flag_array[i] == 1) {
will_send[i] = 1;
&temp_status_array[i]);
if (flag_array[i] == 1) {
will_send[i] = 1;
- smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag,
MPI_COMM_WORLD
,
+ smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag,
comm
,
&status);
i = 1;
}
&status);
i = 1;
}
@@
-310,10
+315,10
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
/* probe before exit in case there are messages to recv */
for (i = 1; i < size; i++) {
/* probe before exit in case there are messages to recv */
for (i = 1; i < size; i++) {
- smpi_mpi_iprobe(i, MPI_ANY_TAG,
MPI_COMM_WORLD
, &flag_array[i],
+ smpi_mpi_iprobe(i, MPI_ANY_TAG,
comm
, &flag_array[i],
&temp_status_array[i]);
if (flag_array[i] == 1)
&temp_status_array[i]);
if (flag_array[i] == 1)
- smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag,
MPI_COMM_WORLD
, &status);
+ smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag,
comm
, &status);
}
}
}
}
@@
-322,7
+327,7
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
else {
/* if root already send a message to this node, don't send one-byte message */
else {
/* if root already send a message to this node, don't send one-byte message */
- smpi_mpi_iprobe(0, MPI_ANY_TAG,
MPI_COMM_WORLD
, &flag_array[0], &status);
+ smpi_mpi_iprobe(0, MPI_ANY_TAG,
comm
, &flag_array[0], &status);
/* send 1-byte message to root */
if (flag_array[0] == 0)
/* send 1-byte message to root */
if (flag_array[0] == 0)
@@
-366,6
+371,9
@@
int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
header_buf[myordering + 1], tag, comm);
}
smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
header_buf[myordering + 1], tag, comm);
}
smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
+ }else{
+ for (i = 0; i < pipe_length; i++)
+ smpi_mpi_wait(&recv_request_array[i], &recv_status_array[i]);
}
}
}
}