Logo AND Algorithmique Numérique Distribuée

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