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
get the right sender for tracing, because it might not be known if MPI_ANY_SOURCE...
[simgrid.git]
/
src
/
smpi
/
smpi_pmpi.c
diff --git
a/src/smpi/smpi_pmpi.c
b/src/smpi/smpi_pmpi.c
index 0422b98b767936fa3bfb19a8ca026f992e7832c6..561ad0d507894fe88c9f47724cba9690d3643a5f 100644
(file)
--- a/
src/smpi/smpi_pmpi.c
+++ b/
src/smpi/smpi_pmpi.c
@@
-28,7
+28,10
@@
int PMPI_Init(int *argc, char ***argv)
{
smpi_process_init(argc, argv);
#ifdef HAVE_TRACING
- TRACE_smpi_init(smpi_process_index());
+ int rank = smpi_process_index();
+ TRACE_smpi_init(rank);
+
+ TRACE_smpi_computing_init(rank);
#endif
smpi_bench_begin();
return MPI_SUCCESS;
@@
-36,8
+39,11
@@
int PMPI_Init(int *argc, char ***argv)
int PMPI_Finalize(void)
{
+ smpi_process_finalize();
smpi_bench_end();
#ifdef HAVE_TRACING
+ int rank = smpi_process_index();
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_finalize(smpi_process_index());
#endif
smpi_process_destroy();
@@
-86,8
+92,12
@@
int PMPI_Abort(MPI_Comm comm, int errorcode)
{
smpi_bench_end();
smpi_process_destroy();
+#ifdef HAVE_TRACING
+ int rank = smpi_process_index();
+ TRACE_smpi_computing_out(rank);
+#endif
// FIXME: should kill all processes in comm instead
-
SIMIX_req
_process_kill(SIMIX_process_self());
+
simcall
_process_kill(SIMIX_process_self());
return MPI_SUCCESS;
}
@@
-100,6
+110,11
@@
double PMPI_Wtime(void)
smpi_bench_begin();
return time;
}
+extern double sg_maxmin_precision;
+double PMPI_Wtick(void)
+{
+ return sg_maxmin_precision;
+}
int PMPI_Address(void *location, MPI_Aint * address)
{
@@
-763,6
+778,25
@@
int PMPI_Comm_free(MPI_Comm * comm)
return retval;
}
+int PMPI_Comm_disconnect(MPI_Comm * comm)
+{
+ /* TODO: wait until all communication in comm are done */
+ int retval;
+
+ smpi_bench_end();
+ if (comm == NULL) {
+ retval = MPI_ERR_ARG;
+ } else if (*comm == MPI_COMM_NULL) {
+ retval = MPI_ERR_COMM;
+ } else {
+ smpi_comm_destroy(*comm);
+ *comm = MPI_COMM_NULL;
+ retval = MPI_SUCCESS;
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
{
int retval;
@@
-896,6
+930,7
@@
int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
TRACE_smpi_send(rank, rank, dst_traced);
@@
-911,6
+946,7
@@
int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
#ifdef HAVE_TRACING
TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
(*request)->send = 1;
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-925,6
+961,8
@@
int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
+ TRACE_smpi_computing_out(rank);
+
TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-934,8
+972,11
@@
int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
retval = MPI_SUCCESS;
}
#ifdef HAVE_TRACING
+ //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
+ if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE);
TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
TRACE_smpi_recv(rank, src_traced, rank);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-949,6
+990,7
@@
int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
TRACE_smpi_send(rank, rank, dst_traced);
@@
-961,6
+1003,7
@@
int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
}
#ifdef HAVE_TRACING
TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-976,6
+1019,7
@@
int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
@@
-996,6
+1040,8
@@
int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
TRACE_smpi_recv(rank, rank, dst_traced);
TRACE_smpi_recv(rank, src_traced, rank);
+ TRACE_smpi_computing_in(rank);
+
#endif
smpi_bench_begin();
return retval;
@@
-1052,6
+1098,43
@@
int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
return retval;
}
+
+
+int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
+ int retval;
+ smpi_bench_end();
+
+ if (status == NULL) {
+ retval = MPI_ERR_ARG;
+ }else if (comm == MPI_COMM_NULL) {
+ retval = MPI_ERR_COMM;
+ } else {
+ smpi_mpi_probe(source, tag, comm, status);
+ retval = MPI_SUCCESS;
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+
+int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
+ int retval;
+ smpi_bench_end();
+
+ if (flag == NULL) {
+ retval = MPI_ERR_ARG;
+ }else if (status == NULL) {
+ retval = MPI_ERR_ARG;
+ }else if (comm == MPI_COMM_NULL) {
+ retval = MPI_ERR_COMM;
+ } else {
+ smpi_mpi_iprobe(source, tag, comm, flag, status);
+ retval = MPI_SUCCESS;
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
int PMPI_Wait(MPI_Request * request, MPI_Status * status)
{
int retval;
@@
-1061,6
+1144,8
@@
int PMPI_Wait(MPI_Request * request, MPI_Status * status)
int rank = request && (*request)->comm != MPI_COMM_NULL
? smpi_comm_rank((*request)->comm)
: -1;
+ TRACE_smpi_computing_out(rank);
+
MPI_Group group = smpi_comm_group((*request)->comm);
int src_traced = smpi_group_rank(group, (*request)->src);
int dst_traced = smpi_group_rank(group, (*request)->dst);
@@
-1080,6
+1165,8
@@
int PMPI_Wait(MPI_Request * request, MPI_Status * status)
if (is_wait_for_receive) {
TRACE_smpi_recv(rank, src_traced, dst_traced);
}
+ TRACE_smpi_computing_in(rank);
+
#endif
smpi_bench_begin();
return retval;
@@
-1093,9
+1180,9
@@
int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
#ifdef HAVE_TRACING
//save requests information for tracing
int i;
- xbt_dynar_t srcs = xbt_dynar_new(sizeof(int),
xbt_free
);
- xbt_dynar_t dsts = xbt_dynar_new(sizeof(int),
xbt_free
);
- xbt_dynar_t recvs = xbt_dynar_new(sizeof(int),
xbt_free
);
+ xbt_dynar_t srcs = xbt_dynar_new(sizeof(int),
NULL
);
+ xbt_dynar_t dsts = xbt_dynar_new(sizeof(int),
NULL
);
+ xbt_dynar_t recvs = xbt_dynar_new(sizeof(int),
NULL
);
for (i = 0; i < count; i++) {
MPI_Request req = requests[i]; //already received requests are no longer valid
if (req) {
@@
-1108,15
+1195,22
@@
int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
xbt_dynar_insert_at(srcs, i, asrc);
xbt_dynar_insert_at(dsts, i, adst);
xbt_dynar_insert_at(recvs, i, arecv);
+ xbt_free(asrc);
+ xbt_free(adst);
+ xbt_free(arecv);
} else {
int *t = xbt_new(int, 1);
xbt_dynar_insert_at(srcs, i, t);
xbt_dynar_insert_at(dsts, i, t);
xbt_dynar_insert_at(recvs, i, t);
+ xbt_free(t);
}
}
int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
+ TRACE_smpi_computing_out(rank_traced);
+
TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
+
#endif
if (index == NULL) {
retval = MPI_ERR_ARG;
@@
-1134,9
+1228,11
@@
int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * sta
}
TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
//clean-up of dynars
- xbt_free(srcs);
- xbt_free(dsts);
- xbt_free(recvs);
+ xbt_dynar_free(&srcs);
+ xbt_dynar_free(&dsts);
+ xbt_dynar_free(&recvs);
+ TRACE_smpi_computing_in(rank_traced);
+
#endif
smpi_bench_begin();
return retval;
@@
-1149,9
+1245,9
@@
int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
#ifdef HAVE_TRACING
//save information from requests
int i;
- xbt_dynar_t srcs = xbt_dynar_new(sizeof(int),
xbt_free
);
- xbt_dynar_t dsts = xbt_dynar_new(sizeof(int),
xbt_free
);
- xbt_dynar_t recvs = xbt_dynar_new(sizeof(int),
xbt_free
);
+ xbt_dynar_t srcs = xbt_dynar_new(sizeof(int),
NULL
);
+ xbt_dynar_t dsts = xbt_dynar_new(sizeof(int),
NULL
);
+ xbt_dynar_t recvs = xbt_dynar_new(sizeof(int),
NULL
);
for (i = 0; i < count; i++) {
MPI_Request req = requests[i]; //all req should be valid in Waitall
int *asrc = xbt_new(int, 1);
@@
-1163,8
+1259,13
@@
int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
xbt_dynar_insert_at(srcs, i, asrc);
xbt_dynar_insert_at(dsts, i, adst);
xbt_dynar_insert_at(recvs, i, arecv);
+ xbt_free(asrc);
+ xbt_free(adst);
+ xbt_free(arecv);
}
int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
+ TRACE_smpi_computing_out(rank_traced);
+
TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
#endif
smpi_mpi_waitall(count, requests, status);
@@
-1180,9
+1281,10
@@
int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
}
TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
//clean-up of dynars
- xbt_free(srcs);
- xbt_free(dsts);
- xbt_free(recvs);
+ xbt_dynar_free(&srcs);
+ xbt_dynar_free(&dsts);
+ xbt_dynar_free(&recvs);
+ TRACE_smpi_computing_in(rank_traced);
#endif
smpi_bench_begin();
return MPI_SUCCESS;
@@
-1211,6
+1313,7
@@
int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm c
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
#endif
@@
-1222,6
+1325,7
@@
int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm c
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1234,6
+1338,7
@@
int PMPI_Barrier(MPI_Comm comm)
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1244,6
+1349,7
@@
int PMPI_Barrier(MPI_Comm comm)
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1258,6
+1364,7
@@
int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
#endif
@@
-1273,6
+1380,7
@@
int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1287,6
+1395,7
@@
int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
#endif
@@
-1304,6
+1413,7
@@
int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1318,6
+1428,7
@@
int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1346,6
+1457,7
@@
int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1362,6
+1474,7
@@
int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1376,6
+1489,7
@@
int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
#endif
@@
-1391,6
+1505,7
@@
int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1405,6
+1520,7
@@
int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
#endif
@@
-1422,6
+1538,7
@@
int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1435,6
+1552,7
@@
int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
#endif
@@
-1448,6
+1566,7
@@
int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1461,6
+1580,7
@@
int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1475,6
+1595,7
@@
int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1488,6
+1609,7
@@
int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1502,6
+1624,7
@@
int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1516,6
+1639,7
@@
int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
smpi_bench_end();
#ifdef HAVE_TRACING
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1544,6
+1668,7
@@
int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1558,6
+1683,7
@@
int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1587,6
+1713,7
@@
int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1601,6
+1728,7
@@
int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
smpi_bench_end();
#ifdef HAVE_TRACING
int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
+ TRACE_smpi_computing_out(rank);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
#endif
if (comm == MPI_COMM_NULL) {
@@
-1619,6
+1747,7
@@
int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
}
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
+ TRACE_smpi_computing_in(rank);
#endif
smpi_bench_begin();
return retval;
@@
-1852,9
+1981,7
@@
int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag,
return not_yet_implemented();
}
-int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
- return not_yet_implemented();
-}
+
int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
return not_yet_implemented();
@@
-1908,10
+2035,6
@@
int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
return not_yet_implemented();
}
-int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
- return not_yet_implemented();
-}
-
int PMPI_Initialized(int* flag) {
return not_yet_implemented();
}