Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
MPI_Scatter() ok
authorgenaud <genaud@48e7efb5-ca39-0410-a469-dd3cf9ba447f>
Fri, 10 Jul 2009 14:05:11 +0000 (14:05 +0000)
committergenaud <genaud@48e7efb5-ca39-0410-a469-dd3cf9ba447f>
Fri, 10 Jul 2009 14:05:11 +0000 (14:05 +0000)
git-svn-id: svn+ssh://scm.gforge.inria.fr/svn/simgrid/simgrid/trunk@6468 48e7efb5-ca39-0410-a469-dd3cf9ba447f

examples/smpi/scatter.c
src/smpi/smpi_coll_private.h
src/smpi/smpi_mpi.c

index 9647aee..2f68c91 100644 (file)
@@ -29,14 +29,20 @@ int ibm_test(int rank, int size)  {
 #define MAXLEN  10000
 
          int success = 1;
-         int root,i,j,k;
-         int out[MAXLEN*64];
-         int in[MAXLEN];
-
-
-         for(j=1,root=0;j<=MAXLEN;j*=10,root=(root+1)%size)  {
-                   if(rank == root)
-                               for(i=0;i<j*size;i++)  out[i] = i;
+         int root=0;
+        int i,j,k;
+         int *out;
+         int *in;
+        
+        
+        out = malloc(MAXLEN*64*sizeof(int));
+        in  = malloc(MAXLEN*sizeof(int));
+
+         for(j=1;j<=MAXLEN;j*=10){
+                root=(root+1)%size;
+                   if (rank == root)
+                               for(i=0;i<j*size;i++)  
+                                out[i] = i;
 
                    MPI_Scatter(out,j,MPI_INT,in,j,MPI_INT,root,MPI_COMM_WORLD);
 
@@ -47,60 +53,77 @@ int ibm_test(int rank, int size)  {
                                }
                    }
          }
-        return( success );
+        free(out);
+        free(in);
          MPI_Barrier( MPI_COMM_WORLD );
+        return( success );
 }
 
 /**
  * small test: the root sends a single distinct double to other processes
  **/
 int small_test( int rank, int size ) {
-  int success=1;
-  int sendcount=size;
-  int recvcount=1;
-  int i;
-  double *sndbuf =NULL;
-  double rcvd;
-  int root=0; // arbitrary choice 
-
-  // on root, initialize sendbuf
-  if (root == rank )  {
-           sndbuf = malloc( sendcount*sizeof(double));
-           for (i=0; i<sendcount;i++) {
-                       sndbuf[i] = (double) rank;
-           }
-  }
-
-   MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd,recvcount,MPI_DOUBLE,root,MPI_COMM_WORLD);
-
-  return(success);
+        int success=1;
+        int retval;
+        int sendcount=1; // one double to each process
+        int recvcount=1;
+        int i;
+        double *sndbuf =NULL;
+        double rcvd;
+        int root=0; // arbitrary choice 
+
+        // on root, initialize sendbuf
+        if (root == rank )  {
+                sndbuf = malloc( size * sizeof(double));
+                for (i=0; i< size;i++) {
+                        sndbuf[i] = (double) i;
+                }
+        }
+
+        retval=MPI_Scatter(sndbuf, sendcount, MPI_DOUBLE, &rcvd,recvcount,MPI_DOUBLE,root,MPI_COMM_WORLD);
+        if (retval != MPI_SUCCESS) {
+                fprintf(stderr,"(%s:%d) MPI_Scatter() returned retval=%d\n",__FILE__,__LINE__,retval);
+                return 0;
+        }
+
+        // verification
+        if ((double)rank != rcvd) {
+                fprintf(stderr,"[%d] has %lf instead of %d\n",rank,rcvd,rank);
+                success=0;
+        }
+        return(success);
 }
 
 
 
 int main(int argc, char **argv)
 {
-  int size, rank;
+        int size, rank;
 
 
-  MPI_Init(&argc, &argv);
-  MPI_Comm_size(MPI_COMM_WORLD, &size);
-  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
+        MPI_Init(&argc, &argv);
+        MPI_Comm_size(MPI_COMM_WORLD, &size);
+        MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
+        /* test 1 */
+        if (0 == rank)
+                printf("** Small Test Result: ... \n");
+        if (! small_test( rank, size ))
+                printf("\t[%d] failed.\n", rank);
+        else
+                printf("\t[%d] ok.\n", rank);
 
 
-  if (0 == rank)
-    printf("** Small Test Result: ... \n");
-  small_test( rank, size );
+        MPI_Barrier( MPI_COMM_WORLD );
 
+        /* test 2 */
+        if (0 == rank)
+                printf("** IBM Test Result: ... \n");
+        if (!ibm_test(rank, size))
+                printf("\t[%d] failed.\n", rank);
+        else
+                printf("\t[%d] ok.\n", rank);
 
-  if (0 == rank)
-    printf("** IBM Test Result: ... \n");
-/*  if (!ibm_test(rank, size))
-    printf("\t[%d] failed.\n", rank);
-  else
-    printf("\t[%d] ok.\n", rank);
-*/
-  MPI_Finalize();
-  return 0;
+        MPI_Finalize();
+        return 0;
 }
index 95a0676..1df83cc 100644 (file)
@@ -1,5 +1,5 @@
 /**
- * $Id$tag 
+ * $Id$tag 
  *
  * smpi_coll_private.h -- functions of smpi_coll.c that are exported to other SMPI modules.
  *
index 6aae3c6..6b767a9 100644 (file)
@@ -264,11 +264,11 @@ int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
 
 
 
-#ifdef DEBUG_REDUCE
+//#ifdef DEBUG_REDUCE
 /**
  * debugging helper function
  **/
-static void print_buffer_int(void *buf, int len, const char *msg, int rank)
+static void print_buffer_int(void *buf, int len, char *msg, int rank)
 {
   int tmp, *v;
   printf("**[%d] %s: ", rank, msg);
@@ -279,8 +279,21 @@ static void print_buffer_int(void *buf, int len, const char *msg, int rank)
   printf("\n");
   free(msg);
 }
-#endif
+static void print_buffer_double(void *buf, int len, char *msg, int rank)
+{
+  int tmp;
+  double *v;
+  printf("**[%d] %s: ", rank, msg);
+  for (tmp = 0; tmp < len; tmp++) {
+    v = buf;
+    printf("[%lf]", v[tmp]);
+  }
+  printf("\n");
+  free(msg);
+}
+
 
+//#endif
 /**
  * MPI_Reduce
  **/
@@ -406,7 +419,7 @@ int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype,
   int cnt=0;  
   int rank;
   int tag=0;
-  char *cbuf;  // to manipulate the void * buffers
+  char *cptr;  // to manipulate the void * buffers
   smpi_mpi_request_t *requests;
   smpi_mpi_request_t request;
   smpi_mpi_status_t status;
@@ -418,24 +431,25 @@ int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype,
 
   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
   if (rank == root) {
-           // i am the root: distribute my sendbuf
-           for (i=0; i < comm->size; i++) {
-                                 cbuf = sendbuf;
-                                 cbuf += i*sendcount*datatype->size;
-                       if ( i!=root ) { // send to processes ...
-
-                                 retval = smpi_create_request((void *)cbuf, sendcount, 
-                                                       datatype, root, i, tag, comm, &(requests[cnt++]));
-                                 if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
-                                           if (MPI_SUCCESS == retval) {
-                                                       smpi_mpi_isend(requests[cnt]);
-                                           }
+          // i am the root: distribute my sendbuf
+          //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
+          cptr = sendbuf;
+          for (i=0; i < comm->size; i++) {
+                  if ( i!=root ) { // send to processes ...
+
+                          retval = smpi_create_request((void *)cptr, sendcount, 
+                                          datatype, root, i, tag, comm, &(requests[cnt]));
+                          if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
+                                  if (MPI_SUCCESS == retval) {
+                                          smpi_mpi_isend(requests[cnt]);
+                                  }
                                  }
                                  cnt++;
                        } 
                        else { // ... except if it's me.
-                                 memcpy(recvbuf, (void *)cbuf, recvcount*recvtype->size*sizeof(char));
+                                 memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
                        }
+                  cptr += sendcount*datatype->size;
            }
            for(i=0; i<cnt; i++) { // wait for send to complete
                            /* FIXME: waitall() should be slightly better */