Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
4e0cb638d6e138ae7fa8dc4ff01ea875ade2697c
[simgrid.git] / teshsuite / smpi / replay-ti-colls / replay-ti-colls.c
1 //exclude from clang static analysis, as there is an intentional uninitialized value passed to MPI calls.
2 #ifndef __clang_analyzer__
3
4 #include "mpi.h"
5 #include <stdio.h>
6 #include <string.h>
7
8 #define BUFSIZE (1024 * 1024)
9 #define BOUNDED(sz) ((sz) < BUFSIZE ? (sz) : BUFSIZE)
10
11 static void setup_recvbuf(int nprocs, int** recvbuf, int** displs, int** counts, int** rcounts)
12 {
13   *recvbuf = malloc(BUFSIZE * nprocs * sizeof(int));
14   for (int i = 0; i < BUFSIZE * nprocs; i++)
15     (*recvbuf)[i] = i;
16
17   *displs  = malloc(nprocs * sizeof(int));
18   *counts  = malloc(nprocs * sizeof(int));
19   *rcounts = malloc(nprocs * sizeof(int));
20   for (int i = 0; i < nprocs; i++) {
21     (*displs)[i]  = i * BUFSIZE;
22     (*counts)[i]  = BOUNDED(i);
23     (*rcounts)[i] = (*counts)[i];
24   }
25 }
26
27 int main(int argc, char** argv)
28 {
29   int nprocs = -1;
30   int rank   = -1;
31
32   /* init */
33   MPI_Init(&argc, &argv);
34   MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
35   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
36
37   int* sendbuf = malloc(BUFSIZE * nprocs * sizeof(int));
38   for (int i = 0; i < BUFSIZE * nprocs; i++)
39     sendbuf[i] = rank;
40
41   int* alltoallvcounts = malloc(nprocs * sizeof(int));
42   for (int i = 0; i < nprocs; i++)
43     alltoallvcounts[i] = BOUNDED(i + rank);
44
45   int* recvbuf;
46   int* displs;
47   int* counts;
48   int* rcounts;
49   if (rank == 0)
50     setup_recvbuf(nprocs, &recvbuf, &displs, &counts, &rcounts);
51
52   // first test, with unallocated non significative buffers
53   MPI_Barrier(MPI_COMM_WORLD);
54   MPI_Bcast(sendbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
55   MPI_Gather(&sendbuf[rank * BUFSIZE], BUFSIZE, MPI_INT, recvbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
56   MPI_Scatter(recvbuf, BUFSIZE, MPI_INT, sendbuf, BUFSIZE, MPI_INT, 0, MPI_COMM_WORLD);
57   MPI_Gatherv(&sendbuf[rank * BUFSIZE], BOUNDED(rank), MPI_INT, recvbuf, rcounts, displs, MPI_INT, 0, MPI_COMM_WORLD);
58   MPI_Scatterv(recvbuf, counts, displs, MPI_INT, sendbuf, BOUNDED(rank), MPI_INT, 0, MPI_COMM_WORLD);
59   MPI_Reduce(sendbuf, recvbuf, BUFSIZE, MPI_INT, MPI_MAX, 0, MPI_COMM_WORLD);
60
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   free(alltoallvcounts);
82   free(sendbuf);
83   free(recvbuf);
84   free(displs);
85   free(counts);
86   free(rcounts);
87
88   MPI_Finalize();
89   return 0;
90 }
91
92 #endif