Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
ff300cb51adaf969ca3291efe1a230e4f4b16b1f
[simgrid.git] / teshsuite / smpi / replay-ti-colls / replay-ti-colls.c
1 #include "mpi.h"
2 #include <assert.h>
3 #include <stdio.h>
4 #include <string.h>
5
6 #define BUFSIZE (1024 * 1024)
7 #define BOUNDED(sz) ((sz) < BUFSIZE ? (sz) : BUFSIZE)
8
9 static void setup_recvbuf(int nprocs, int** recvbuf, int** displs, int** counts, int** rcounts)
10 {
11   *recvbuf = malloc(BUFSIZE * nprocs * sizeof(int));
12   assert(*recvbuf != NULL);
13   for (int i = 0; i < BUFSIZE * nprocs; i++)
14     (*recvbuf)[i] = i;
15
16   *displs  = malloc(nprocs * sizeof(int));
17   *counts  = malloc(nprocs * sizeof(int));
18   *rcounts = malloc(nprocs * sizeof(int));
19   for (int i = 0; i < nprocs; i++) {
20     (*displs)[i]  = i * BUFSIZE;
21     (*counts)[i]  = BOUNDED(i);
22     (*rcounts)[i] = (*counts)[i];
23   }
24 }
25
26 int main(int argc, char** argv)
27 {
28   int nprocs = -1;
29   int rank   = -1;
30
31   /* init */
32   MPI_Init(&argc, &argv);
33   MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
34   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
35
36   int* sendbuf = malloc(BUFSIZE * nprocs * sizeof(int));
37   for (int i = 0; i < BUFSIZE * nprocs; i++)
38     sendbuf[i] = rank;
39
40   int* alltoallvcounts = malloc(nprocs * sizeof(int));
41   for (int i = 0; i < nprocs; i++)
42     alltoallvcounts[i] = BOUNDED(i + rank);
43
44   int* dummy_buffer = malloc(sizeof(int));
45   // initialize buffers with an invalid value (we want to trigger a valgrind error if they are used)
46   int* recvbuf      = dummy_buffer + 1;
47   int* displs       = dummy_buffer + 1;
48   int* counts       = dummy_buffer + 1;
49   int* rcounts      = dummy_buffer + 1;
50   if (rank == 0)
51     setup_recvbuf(nprocs, &recvbuf, &displs, &counts, &rcounts);
52
53   // first test, with unallocated non significative buffers
54   MPI_Barrier(MPI_COMM_WORLD);
55   MPI_Bcast(sendbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
56   MPI_Gather(&sendbuf[rank * BUFSIZE], BUFSIZE, MPI_INT, recvbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
57   MPI_Scatter(recvbuf, BUFSIZE, MPI_INT, sendbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
58   MPI_Gatherv(&sendbuf[rank * BUFSIZE], BOUNDED(rank), MPI_INT, recvbuf, rcounts, displs, MPI_INT, 0, MPI_COMM_WORLD);
59   MPI_Scatterv(recvbuf, counts, displs, MPI_INT, sendbuf, BOUNDED(rank), MPI_INT, 0, MPI_COMM_WORLD);
60   MPI_Reduce(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, 0, MPI_COMM_WORLD);
61
62   free(dummy_buffer);
63   if (rank != 0)
64     setup_recvbuf(nprocs, &recvbuf, &displs, &counts, &rcounts);
65
66   MPI_Barrier(MPI_COMM_WORLD);
67   MPI_Bcast(sendbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
68   MPI_Gather(&sendbuf[rank * BUFSIZE], BUFSIZE, MPI_INT, recvbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
69   MPI_Scatter(recvbuf, BUFSIZE, MPI_INT, sendbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
70   MPI_Gatherv(&sendbuf[rank * BUFSIZE], BOUNDED(rank), MPI_INT, recvbuf, rcounts, displs, MPI_INT, 0, MPI_COMM_WORLD);
71   MPI_Scatterv(recvbuf, counts, displs, MPI_INT, sendbuf, BOUNDED(rank), MPI_INT, 0, MPI_COMM_WORLD);
72   MPI_Reduce(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, 0, MPI_COMM_WORLD);
73   MPI_Allgather(sendbuf, BUFSIZE, MPI_INT, recvbuf, BUFSIZE, MPI_INT, MPI_COMM_WORLD);
74   MPI_Alltoall(recvbuf, BUFSIZE, MPI_INT, sendbuf, BUFSIZE, MPI_INT, MPI_COMM_WORLD);
75   MPI_Allgatherv(sendbuf, BOUNDED(rank), MPI_INT, recvbuf, rcounts, displs, MPI_INT, MPI_COMM_WORLD);
76   MPI_Alltoallv(recvbuf, alltoallvcounts, displs, MPI_INT, sendbuf, alltoallvcounts, displs, MPI_INT, MPI_COMM_WORLD);
77   MPI_Allreduce(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
78   MPI_Reduce_scatter(sendbuf, recvbuf, rcounts, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
79   MPI_Scan(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
80   MPI_Exscan(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
81   MPI_Barrier(MPI_COMM_WORLD);
82
83   free(alltoallvcounts);
84   free(sendbuf);
85   free(recvbuf);
86   free(displs);
87   free(counts);
88   free(rcounts);
89
90   MPI_Finalize();
91   return 0;
92 }