Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
correct a few bcast algorithms and change the test to actually use them all
authorAugustin Degomme <degomme@idpann.imag.fr>
Wed, 19 Jun 2013 15:08:45 +0000 (17:08 +0200)
committerAugustin Degomme <degomme@idpann.imag.fr>
Wed, 19 Jun 2013 15:08:45 +0000 (17:08 +0200)
src/smpi/colls/bcast-arrival-nb.c
src/smpi/colls/bcast-arrival-pattern-aware.c
teshsuite/smpi/bcast_coll.c
teshsuite/smpi/bcast_coll.tesh

index 9ff27b4..be088ed 100644 (file)
@@ -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 to_clean[MAX_NODE];
   int header_buf[HEADER_SIZE];
   char temp_buf[MAX_NODE];
 
@@ -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;
+    to_clean[i]=0;
   }
   //  printf("YYY\n");
 
@@ -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++;
+              to_clean[i]=0;
               break;
             }
           }
@@ -130,6 +132,9 @@ int smpi_coll_tuned_bcast_arrival_nb(void *buf, int count,
 
 
       }                         /* while loop */
+      
+      for(i=0; i<size; i++)
+        if(to_clean[i]!=0)smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
     }
 
     /* non-root */
@@ -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);
+      }else{
+        for (i = 0; i < pipe_length; i++) 
+          smpi_mpi_wait(&recv_request_array[i], &recv_status_array[i]);
       }
 
     }
index f4a482c..a9df449 100644 (file)
@@ -27,7 +27,7 @@ int smpi_coll_tuned_bcast_arrival_pattern_aware(void *buf, int count,
   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];
 
@@ -70,6 +70,7 @@ int smpi_coll_tuned_bcast_arrival_pattern_aware(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;
+    to_clean[i] = 0;
   }
 
   /* when a message is smaller than a block size => no pipeline */
@@ -274,6 +275,7 @@ int smpi_coll_tuned_bcast_arrival_pattern_aware(void *buf, int count,
 
 
               already_sent[i] = 1;
+              to_clean[i]=1;
               sent_count++;
               break;
             }
@@ -282,6 +284,9 @@ int smpi_coll_tuned_bcast_arrival_pattern_aware(void *buf, int count,
 
       }                         /* while loop */
 
+      for(i=0; i<size; i++)
+        if(to_clean[i]!=0)smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD,
+                     &status);
       //total = MPI_Wtime() - start2;
       //total *= 1000;
       //printf("Node zero iter = %d time = %.2f\n",iteration,total);
@@ -331,8 +336,10 @@ int smpi_coll_tuned_bcast_arrival_pattern_aware(void *buf, int count,
                     header_buf[myordering + 1], tag, comm);
         }
         smpi_mpi_waitall((pipe_length), send_request_array, send_status_array);
-      }
-
+      }else{
+          smpi_mpi_waitall(pipe_length, recv_request_array, recv_status_array);
+          }
+    
     }
 
     free(send_request_array);
index 9bb6a04..67c0c13 100644 (file)
@@ -31,7 +31,24 @@ int main(int argc, char **argv)
   printf("[%d] number of values equals to 17: %d\n", rank, good);
 
   MPI_Barrier(MPI_COMM_WORLD);
+  xbt_free(values);
 
+  count = 4096;
+  values = (int *) xbt_malloc(count * sizeof(int));  
+
+  for (i = 0; i < count; i++)
+    values[i] = (size -1 == rank) ? 17 : 3;
+
+  status = MPI_Bcast(values, count, MPI_INT, size-1, MPI_COMM_WORLD);
+
+  good = 0;
+  for (i = 0; i < count; i++)
+    if (values[i]==17) good++;
+  printf("[%d] number of values equals to 17: %d\n", rank, good);
+
+
+  
+  
   if (rank == 0) {
     if (status != MPI_SUCCESS) {
       printf("bcast returned %d\n", status);
index 7dc08a5..5fb6b97 100644 (file)
@@ -22,18 +22,35 @@ $ ../../bin/smpirun -map -hostfile ${srcdir:=.}/hostfile -platform ${srcdir:=.}/
 > [rank 15] -> Tremblay
 > [0.000000] [surf_config/INFO] Switching workstation model to compound since you changed the network and/or cpu model(s)
 > [0] number of values equals to 17: 2048
+> [0] number of values equals to 17: 4096
+> [10] number of values equals to 17: 2048
+> [10] number of values equals to 17: 4096
+> [11] number of values equals to 17: 2048
+> [11] number of values equals to 17: 4096
+> [12] number of values equals to 17: 2048
+> [12] number of values equals to 17: 4096
+> [13] number of values equals to 17: 2048
+> [13] number of values equals to 17: 4096
+> [14] number of values equals to 17: 2048
+> [14] number of values equals to 17: 4096
+> [15] number of values equals to 17: 2048
+> [15] number of values equals to 17: 4096
 > [1] number of values equals to 17: 2048
+> [1] number of values equals to 17: 4096
 > [2] number of values equals to 17: 2048
+> [2] number of values equals to 17: 4096
 > [3] number of values equals to 17: 2048
+> [3] number of values equals to 17: 4096
 > [4] number of values equals to 17: 2048
+> [4] number of values equals to 17: 4096
 > [5] number of values equals to 17: 2048
+> [5] number of values equals to 17: 4096
 > [6] number of values equals to 17: 2048
+> [6] number of values equals to 17: 4096
 > [7] number of values equals to 17: 2048
+> [7] number of values equals to 17: 4096
 > [8] number of values equals to 17: 2048
+> [8] number of values equals to 17: 4096
 > [9] number of values equals to 17: 2048
-> [10] number of values equals to 17: 2048
-> [11] number of values equals to 17: 2048
-> [12] number of values equals to 17: 2048
-> [13] number of values equals to 17: 2048
-> [14] number of values equals to 17: 2048
-> [15] number of values equals to 17: 2048
+> [9] number of values equals to 17: 4096
+