Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
In smpi_mpi.c:
[simgrid.git] / src / smpi / smpi_mpi.c
index 2b9caeb..2586730 100644 (file)
@@ -1,3 +1,14 @@
+/* $Id$tag */
+
+/* smpi_mpi.c -- 
+ *
+ * Eventually will contain the user level MPI primitives and its corresponding 
+ * internal wrapper. The implementations of these primitives should go to specific
+ * files. For example, a SMPI_MPI_Bcast() in this file, should call the wrapper 
+ * smpi_mpi_bcast(), which decides which implementation to call. Currently, it
+ * calls nary_tree_bcast() in smpi_coll.c.  (Stéphane Genaud).
+ * */
+
 
 
 #include "private.h"
@@ -102,10 +113,11 @@ int SMPI_MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
                    int tag, MPI_Comm comm, MPI_Request * request)
 {
   int retval = MPI_SUCCESS;
+  int rank;
 
   smpi_bench_end();
-
-  retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
+  rank = smpi_mpi_comm_rank(comm);
+  retval = smpi_create_request(buf, count, datatype, src, rank, tag, comm,
                                request);
   if (NULL != *request && MPI_SUCCESS == retval) {
     retval = smpi_mpi_irecv(*request);
@@ -120,11 +132,13 @@ int SMPI_MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src,
                   int tag, MPI_Comm comm, MPI_Status * status)
 {
   int retval = MPI_SUCCESS;
+  int rank;
   smpi_mpi_request_t request;
 
   smpi_bench_end();
 
-  retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
+  rank = smpi_mpi_comm_rank(comm);
+  retval = smpi_create_request(buf, count, datatype, src, rank, tag, comm,
                                &request);
   if (NULL != request && MPI_SUCCESS == retval) {
     retval = smpi_mpi_irecv(request);
@@ -143,10 +157,12 @@ int SMPI_MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
                    int tag, MPI_Comm comm, MPI_Request * request)
 {
   int retval = MPI_SUCCESS;
+  int rank;
 
   smpi_bench_end();
 
-  retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
+  rank = smpi_mpi_comm_rank(comm);
+  retval = smpi_create_request(buf, count, datatype, rank, dst, tag, comm,
                                request);
   if (NULL != *request && MPI_SUCCESS == retval) {
     retval = smpi_mpi_isend(*request);
@@ -164,11 +180,13 @@ int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
                   int tag, MPI_Comm comm)
 {
   int retval = MPI_SUCCESS;
+  int rank;
   smpi_mpi_request_t request;
 
   smpi_bench_end();
 
-  retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
+  rank = smpi_mpi_comm_rank(comm);
+  retval = smpi_create_request(buf, count, datatype, rank, dst, tag, comm,
                                &request);
   if (NULL != request && MPI_SUCCESS == retval) {
     retval = smpi_mpi_isend(request);
@@ -242,18 +260,33 @@ int retval = MPI_SUCCESS;
  **/
 int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
 {
-  return smpi_mpi_wait(*request, status);
+  int retval;
+
+  smpi_bench_end();
+  retval = smpi_mpi_wait(*request, status);
+  smpi_bench_begin();
+  return retval;
 }
 
 int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
 {
-  return smpi_mpi_waitall(count, requests, status);
+  int retval;
+
+  smpi_bench_end();
+  retval = smpi_mpi_waitall(count, requests, status);
+  smpi_bench_begin();
+  return retval;
 }
 
 int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index,
                      MPI_Status status[])
 {
-  return smpi_mpi_waitany(count, requests, index, status);
+  int retval;
+
+  smpi_bench_end();
+  retval = smpi_mpi_waitany(count, requests, index, status);
+  smpi_bench_begin();
+  return retval;
 }
 
 /**
@@ -296,6 +329,9 @@ int smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
                    MPI_Comm comm)
 {
   int retval = MPI_SUCCESS;
+  int rank = smpi_mpi_comm_rank(comm);
+
+  DEBUG1("<%d> entered smpi_mpi_bcast(). Calls nary_tree_bcast()",rank);
   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
   return retval;
@@ -358,7 +394,7 @@ int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
         int rank;
         int size;
         int i;
-        int tag = 0;
+        int system_tag = 666;
         smpi_mpi_request_t *requests;
         smpi_mpi_request_t request;
 
@@ -366,13 +402,14 @@ int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
 
         rank = smpi_mpi_comm_rank(comm);
         size = comm->size;
+        DEBUG1("<%d> entered smpi_mpi_reduce()",rank);
 
         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
 
 #ifdef DEBUG_REDUCE
                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
 #endif
-                retval = smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
+                retval = smpi_create_request(sendbuf, count, datatype, rank, root, system_tag, comm,
                                         &request);
                 smpi_mpi_isend(request);
                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
@@ -397,7 +434,7 @@ int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
                         // reminder: for smpi_create_request() the src is always the process sending.
                         src = i < root ? i : i + 1;
                         retval = smpi_create_request(tmpbufs[i], count, datatype,
-                                        src, root, tag, comm, &(requests[i]));
+                                        src, root, system_tag, comm, &(requests[i]));
                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
                                 if (MPI_SUCCESS == retval) {
                                         smpi_mpi_irecv(requests[i]);
@@ -408,8 +445,9 @@ int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
                 for (i = 0; i < size-1; i++) {
                         int index = MPI_UNDEFINED;
                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
+                        DEBUG3("<%d> waitany() unblocked by reception (completes request[%d]) (%d reqs remaining)",
+                                        rank,index,size-i-2);
 #ifdef DEBUG_REDUCE
-                        printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
                                         rank);
 #endif
@@ -452,7 +490,7 @@ int retval = MPI_SUCCESS;
 /**
  * MPI_Allreduce
  *
- * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
+ * Same as MPI_Reduce except that the result appears in the receive buffer of all the group members.
  **/
 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
                    MPI_Op op, MPI_Comm comm )
@@ -550,8 +588,8 @@ int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype,
  * ompi/mca/coll/tuned/coll_tuned_module.c
  **/
 int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
-                        void *recvbuf, int recvcount, MPI_Datatype recvtype,
-                          MPI_Comm comm)
+                         void *recvbuf, int recvcount, MPI_Datatype recvtype,
+                           MPI_Comm comm)
 {
   int retval = MPI_SUCCESS;
   int block_dsize;
@@ -561,7 +599,7 @@ int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype,
 
   rank = smpi_mpi_comm_rank(comm);
   block_dsize = datatype->size * sendcount;
-  INFO2("[%d] optimized alltoall() called. Block size sent to each rank=%d.\n",rank,block_dsize);
+  DEBUG2("<%d> optimized alltoall() called. Block size sent to each rank: %d bytes.",rank,block_dsize);
 
   if ((block_dsize < 200) && (comm->size > 12)) {
            retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
@@ -581,6 +619,30 @@ int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype,
   return retval;
 }
 
+/**
+ * MPI_Alltoallv user entry point
+ * 
+ * As in OpenMPI, alltoallv is not optimized
+ * ompi/mca/coll/basic/coll_basic_alltoallv.c 
+ **/
+int SMPI_MPI_Alltoallv(void *sendbuf, int *scounts, int *sdisps, MPI_Datatype datatype, 
+                          void *recvbuf, int *rcounts, int *rdisps, MPI_Datatype recvtype,
+                            MPI_Comm comm)
+{
+  int retval = MPI_SUCCESS;
+  int rank;
+
+  smpi_bench_end();
+  rank = smpi_mpi_comm_rank(comm);
+  DEBUG1("<%d> basic alltoallv() called.",rank);
+
+  retval = smpi_coll_basic_alltoallv(sendbuf, scounts, sdisps, datatype, 
+                                     recvbuf, rcounts, rdisps, recvtype,
+                                     comm); 
+  smpi_bench_begin();
+  return retval;
+}
+
 
 
 
@@ -718,7 +780,12 @@ int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
 
 double SMPI_MPI_Wtime(void)
 {
-  return (SIMIX_get_clock());
+  double time;
+
+  smpi_bench_end();
+  time = SIMIX_get_clock();
+  smpi_bench_begin();
+  return time;
 }