Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Use xbt_malloc that never return null.
[simgrid.git] / teshsuite / smpi / replay-ti-colls / replay-ti-colls.c
index 4e0cb63..c14b802 100644 (file)
@@ -1,6 +1,3 @@
-//exclude from clang static analysis, as there is an intentional uninitialized value passed to MPI calls.
-#ifndef __clang_analyzer__
-
 #include "mpi.h"
 #include <stdio.h>
 #include <string.h>
 
 static void setup_recvbuf(int nprocs, int** recvbuf, int** displs, int** counts, int** rcounts)
 {
-  *recvbuf = malloc(BUFSIZE * nprocs * sizeof(int));
+  *recvbuf = xbt_malloc(BUFSIZE * nprocs * sizeof(int));
   for (int i = 0; i < BUFSIZE * nprocs; i++)
     (*recvbuf)[i] = i;
 
-  *displs  = malloc(nprocs * sizeof(int));
-  *counts  = malloc(nprocs * sizeof(int));
-  *rcounts = malloc(nprocs * sizeof(int));
+  *displs  = xbt_malloc(nprocs * sizeof(int));
+  *counts  = xbt_malloc(nprocs * sizeof(int));
+  *rcounts = xbt_malloc(nprocs * sizeof(int));
   for (int i = 0; i < nprocs; i++) {
     (*displs)[i]  = i * BUFSIZE;
     (*counts)[i]  = BOUNDED(i);
@@ -34,18 +31,20 @@ int main(int argc, char** argv)
   MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
-  int* sendbuf = malloc(BUFSIZE * nprocs * sizeof(int));
+  int* sendbuf = xbt_malloc(BUFSIZE * nprocs * sizeof(int));
   for (int i = 0; i < BUFSIZE * nprocs; i++)
     sendbuf[i] = rank;
 
-  int* alltoallvcounts = malloc(nprocs * sizeof(int));
+  int* alltoallvcounts = xbt_malloc(nprocs * sizeof(int));
   for (int i = 0; i < nprocs; i++)
     alltoallvcounts[i] = BOUNDED(i + rank);
 
-  int* recvbuf;
-  int* displs;
-  int* counts;
-  int* rcounts;
+  int* dummy_buffer = xbt_malloc(sizeof(int));
+  // initialize buffers with an invalid value (we want to trigger a valgrind error if they are used)
+  int* recvbuf      = dummy_buffer + 1;
+  int* displs       = dummy_buffer + 1;
+  int* counts       = dummy_buffer + 1;
+  int* rcounts      = dummy_buffer + 1;
   if (rank == 0)
     setup_recvbuf(nprocs, &recvbuf, &displs, &counts, &rcounts);
 
@@ -58,6 +57,7 @@ int main(int argc, char** argv)
   MPI_Scatterv(recvbuf, counts, displs, MPI_INT, sendbuf, BOUNDED(rank), MPI_INT, 0, MPI_COMM_WORLD);
   MPI_Reduce(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, 0, MPI_COMM_WORLD);
 
+  xbt_free(dummy_buffer);
   if (rank != 0)
     setup_recvbuf(nprocs, &recvbuf, &displs, &counts, &rcounts);
 
@@ -78,15 +78,13 @@ int main(int argc, char** argv)
   MPI_Exscan(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
   MPI_Barrier(MPI_COMM_WORLD);
 
-  free(alltoallvcounts);
-  free(sendbuf);
-  free(recvbuf);
-  free(displs);
-  free(counts);
-  free(rcounts);
+  xbt_free(alltoallvcounts);
+  xbt_free(sendbuf);
+  xbt_free(recvbuf);
+  xbt_free(displs);
+  xbt_free(counts);
+  xbt_free(rcounts);
 
   MPI_Finalize();
   return 0;
 }
-
-#endif