Logo AND Algorithmique Numérique Distribuée

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