Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[examples,smpi,MM] get some positive value for the time
authorjean-noel quintin <jnquintin@dhcp-892b9b4c.ucd.ie>
Wed, 10 Oct 2012 10:27:57 +0000 (11:27 +0100)
committerjean-noel quintin <jnquintin@dhcp-892b9b4c.ucd.ie>
Wed, 10 Oct 2012 10:27:57 +0000 (11:27 +0100)
examples/smpi/MM/2.5D_MM.c
examples/smpi/MM/Summa.c

index eb61096..054cdf1 100644 (file)
@@ -195,7 +195,7 @@ double two_dot_five(
     MPI_Barrier(my_world);
   }
   end_time_intern = MPI_Wtime();
-  communication_time += start_time - end_time_intern;
+  communication_time += end_time_intern - start_time;
 
   XBT_INFO( "group %zu NB_block: %zu, NB_groups %zu\n"
               ,group,NB_Block, NB_groups);
@@ -237,10 +237,10 @@ double two_dot_five(
 
   MPI_Barrier(my_world);
   end_time = MPI_Wtime();
-  time = start_time - end_time;
-  double reduce_time = start_time_reduce - end_time_reduce;
-  printf("communication time: %le reduce time: %le nanoseconds, "
-         "total time: %le nanoseconds\n",communication_time,reduce_time,time);
+  time = end_time - start_time;
+  double reduce_time = end_time_reduce - start_time_reduce;
+  printf("communication time: %le reduce time: %le seconds, "
+         "total time: %le seconds\n",communication_time,reduce_time,time);
   MPI_Barrier(my_world);
 
 #if CHECK_25D
index 9edf3c8..12c2370 100644 (file)
@@ -125,7 +125,7 @@ inline double Summa(
       XBT_DEBUG("position of B_b: %zu \n", pos_b);
     }
     end_time_intern = MPI_Wtime();
-    communication_time += start_time_intern - end_time_intern;
+    communication_time += end_time_intern - start_time_intern;
 
     MPI_Barrier(row_comm);
     MPI_Barrier(col_comm);
@@ -142,16 +142,16 @@ inline double Summa(
           c[i*ldc+j] += B_a[i*lda_local+k]*B_b[k*ldb_local+j];
 
     end_time_intern = MPI_Wtime();
-    computation_time += start_time_intern - end_time_intern;
+    computation_time += end_time_intern - start_time_intern;
 
   }
   MPI_Barrier(row_comm);
   MPI_Barrier(col_comm);
 
   end_time = MPI_Wtime();
-  time = start_time - end_time;
-  printf("communication time: %le nanoseconds, "
-         "computation time: %le nanoseconds\n",
+  time = end_time - start_time ;
+  printf("communication time: %le seconds, "
+         "computation time: %le seconds\n",
          communication_time, computation_time);