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
Add support for various MPI_Type functions, to handle creation of new MPI types
[simgrid.git]
/
src
/
smpi
/
smpi_pmpi.c
diff --git
a/src/smpi/smpi_pmpi.c
b/src/smpi/smpi_pmpi.c
index
ba07a45
..
16fccc0
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;
@@
-39,6
+42,8
@@
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();
@@
-87,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;
}
@@
-101,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)
{
@@
-124,8
+138,8
@@
int PMPI_Type_free(MPI_Datatype * datatype)
if (!datatype) {
retval = MPI_ERR_ARG;
} else {
- // FIXME: always fail for now
- retval = MPI_
ERR_TYPE
;
+ smpi_datatype_free(datatype);
+ retval = MPI_
SUCCESS
;
}
smpi_bench_begin();
return retval;
@@
-908,6
+922,7
@@
int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
return retval;
}
+
int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
int tag, MPI_Comm comm, MPI_Request * request)
{
@@
-916,6
+931,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);
@@
-931,11
+947,14
@@
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;
}
+
+
int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
MPI_Comm comm, MPI_Status * status)
{
@@
-945,6
+964,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) {
@@
-954,8
+975,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;
@@
-969,6
+993,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);
@@
-981,6
+1006,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;
@@
-996,6
+1022,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__);
@@
-1016,6
+1043,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;
@@
-1072,6
+1101,56
@@
int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
return retval;
}
+int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
+{
+ int retval;
+
+ smpi_bench_end();
+ if (flag == NULL) {
+ retval = MPI_ERR_ARG;
+ } else {
+ *flag = smpi_mpi_testall(count, requests, statuses);
+ retval = MPI_SUCCESS;
+ }
+ smpi_bench_begin();
+ 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;
@@
-1081,6
+1160,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);
@@
-1100,6
+1181,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;
@@
-1113,9
+1196,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) {
@@
-1128,15
+1211,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;
@@
-1154,9
+1244,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;
@@
-1169,9
+1261,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);
@@
-1183,8
+1275,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);
@@
-1200,9
+1297,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;
@@
-1231,6
+1329,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
@@
-1242,6
+1341,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;
@@
-1254,6
+1354,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) {
@@
-1264,6
+1365,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;
@@
-1278,6
+1380,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
@@
-1293,6
+1396,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;
@@
-1307,6
+1411,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
@@
-1324,6
+1429,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;
@@
-1338,6
+1444,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) {
@@
-1366,6
+1473,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) {
@@
-1382,6
+1490,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;
@@
-1396,6
+1505,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
@@
-1411,6
+1521,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;
@@
-1425,6
+1536,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
@@
-1442,6
+1554,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;
@@
-1455,6
+1568,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
@@
-1468,6
+1582,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;
@@
-1481,6
+1596,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) {
@@
-1495,6
+1611,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;
@@
-1508,6
+1625,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) {
@@
-1522,6
+1640,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;
@@
-1536,6
+1655,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) {
@@
-1564,6
+1684,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;
@@
-1578,6
+1699,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) {
@@
-1607,6
+1729,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;
@@
-1621,6
+1744,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) {
@@
-1639,6
+1763,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;
@@
-1684,12
+1809,117
@@
int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
return retval;
}
+int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
+ int retval;
+
+ smpi_bench_end();
+ if (old_type == MPI_DATATYPE_NULL) {
+ retval = MPI_ERR_TYPE;
+ } else if (count<=0){
+ retval = MPI_ERR_COUNT;
+ } else {
+ retval = smpi_datatype_contiguous(count, old_type, new_type);
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+int PMPI_Type_commit(MPI_Datatype* datatype) {
+ int retval;
+
+ smpi_bench_end();
+ if (datatype == MPI_DATATYPE_NULL) {
+ retval = MPI_ERR_TYPE;
+ } else {
+ smpi_datatype_commit(datatype);
+ retval = MPI_SUCCESS;
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+
+int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
+ int retval;
+
+ smpi_bench_end();
+ if (old_type == MPI_DATATYPE_NULL) {
+ retval = MPI_ERR_TYPE;
+ } else if (count<=0 || blocklen<=0){
+ retval = MPI_ERR_COUNT;
+ } else {
+ retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
+ int retval;
+
+ smpi_bench_end();
+ if (old_type == MPI_DATATYPE_NULL) {
+ retval = MPI_ERR_TYPE;
+ } else if (count<=0 || blocklen<=0){
+ retval = MPI_ERR_COUNT;
+ } else {
+ retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+
+int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
+ int retval;
+
+ smpi_bench_end();
+ if (old_type == MPI_DATATYPE_NULL) {
+ retval = MPI_ERR_TYPE;
+ } else if (count<=0){
+ retval = MPI_ERR_COUNT;
+ } else {
+ retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
+ int retval;
+
+ smpi_bench_end();
+ if (old_type == MPI_DATATYPE_NULL) {
+ retval = MPI_ERR_TYPE;
+ } else if (count<=0){
+ retval = MPI_ERR_COUNT;
+ } else {
+ retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
+ }
+ smpi_bench_begin();
+ return retval;
+}
+
+
+int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
+ int retval;
+
+ smpi_bench_end();
+ if (count<=0){
+ retval = MPI_ERR_COUNT;
+ } else {
+ retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
+ }
+ smpi_bench_begin();
+ return retval;}
+
+
/* The following calls are not yet implemented and will fail at runtime. */
/* Once implemented, please move them above this notice. */
static int not_yet_implemented(void) {
-
xbt_die
("Not yet implemented");
- return MPI_
ERR_OTHER
;
+
XBT_WARN
("Not yet implemented");
+ return MPI_
SUCCESS
;
}
int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
@@
-1780,9
+2010,6
@@
int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
return not_yet_implemented();
}
-int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
- return not_yet_implemented();
-}
int PMPI_Cancel(MPI_Request* request) {
return not_yet_implemented();
@@
-1808,30
+2035,6
@@
int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcou
return not_yet_implemented();
}
-int PMPI_Type_commit(MPI_Datatype* datatype) {
- return not_yet_implemented();
-}
-
-int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
- return not_yet_implemented();
-}
-
-int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
- return not_yet_implemented();
-}
-
-int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
- return not_yet_implemented();
-}
-
-int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
- return not_yet_implemented();
-}
-
-int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
- return not_yet_implemented();
-}
-
int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
return not_yet_implemented();
}
@@
-1872,9
+2075,6
@@
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();
@@
-1916,10
+2116,6
@@
int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int out
return not_yet_implemented();
}
-int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses) {
- return not_yet_implemented();
-}
-
int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
return not_yet_implemented();
}
@@
-1928,10
+2124,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();
}