Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Reduce the size of partial shared malloc tests.
[simgrid.git] / teshsuite / smpi / mpich3-test / perf / transp-datatype.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2006 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 /*    modified 01/23/2011 by Jim Hoekstra - ISU
7  *      changed test to follow mtest_init/mtest_finalize convention
8  *      The following changes are based on suggestions from Chris Sadlo:
9  *        variable row changed to col.
10  *        manual transpose - code added to perform 'swap'.
11  *        MPI_Send/MPI_Recv involving xpose changed.
12  */
13
14 /* This is based on an example in the MPI standard and a bug report submitted
15    by Alexandr Konovalov of Intel */
16
17 #include "mpi.h"
18 #include <stdio.h>
19 #include "mpitest.h"
20
21 #define SIZE 100
22 #define ITER 100
23
24 int main(int argc, char *argv[])
25 {
26     int i, j, k;
27     static double a[SIZE][SIZE], b[SIZE][SIZE];
28     double t1, t2, t, ts, tst;
29     double temp;
30     int myrank, mysize, errs = 0;
31     MPI_Status status;
32     MPI_Aint sizeofreal;
33
34     MPI_Datatype col, xpose;
35
36     MTest_Init(&argc, &argv);
37     MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
38     MPI_Comm_size(MPI_COMM_WORLD, &mysize);
39     if (mysize != 2) {
40         fprintf(stderr, "This test must be run with 2 processes\n");
41         MPI_Abort(MPI_COMM_WORLD, 1);
42     }
43
44     MPI_Type_extent(MPI_DOUBLE, &sizeofreal);
45
46     MPI_Type_vector(SIZE, 1, SIZE, MPI_DOUBLE, &col);
47     MPI_Type_hvector(SIZE, 1, sizeofreal, col, &xpose);
48     MPI_Type_commit(&xpose);
49
50     /* Preset the arrays so that they're in memory */
51     for (i = 0; i < SIZE; i++)
52         for (j = 0; j < SIZE; j++) {
53             a[i][j] = 0;
54             b[i][j] = 0;
55         }
56     a[SIZE - 1][0] = 1;
57
58     /* Time the transpose example */
59     MPI_Barrier(MPI_COMM_WORLD);
60     t1 = MPI_Wtime();
61     for (i = 0; i < ITER; i++) {
62         if (myrank == 0)
63             MPI_Send(&a[0][0], SIZE * SIZE, MPI_DOUBLE, 1, 0, MPI_COMM_WORLD);
64         else
65             MPI_Recv(&b[0][0], 1, xpose, 0, 0, MPI_COMM_WORLD, &status);
66     }
67     t2 = MPI_Wtime();
68     t = (t2 - t1) / ITER;
69
70     /* Time sending the same amount of data, but without the transpose */
71     MPI_Barrier(MPI_COMM_WORLD);
72     t1 = MPI_Wtime();
73     for (i = 0; i < ITER; i++) {
74         if (myrank == 0) {
75             MPI_Send(&a[0][0], sizeof(a), MPI_BYTE, 1, 0, MPI_COMM_WORLD);
76         }
77         else {
78             MPI_Recv(&b[0][0], sizeof(b), MPI_BYTE, 0, 0, MPI_COMM_WORLD, &status);
79         }
80     }
81     t2 = MPI_Wtime();
82     ts = (t2 - t1) / ITER;
83
84     /* Time sending the same amount of data, with the transpose done
85      * as a separate step */
86     MPI_Barrier(MPI_COMM_WORLD);
87     t1 = MPI_Wtime();
88     for (k = 0; k < ITER; k++) {
89         if (myrank == 0) {
90             MPI_Send(&a[0][0], sizeof(a), MPI_BYTE, 1, 0, MPI_COMM_WORLD);
91         }
92         else {
93             MPI_Recv(&b[0][0], sizeof(b), MPI_BYTE, 0, 0, MPI_COMM_WORLD, &status);
94             for (i = 0; i < SIZE; i++)
95                 for (j = i; j < SIZE; j++) {
96                     temp = b[j][i];
97                     b[j][i] = b[i][j];
98                     b[i][j] = temp;
99                 }
100         }
101     }
102     t2 = MPI_Wtime();
103     tst = (t2 - t1) / ITER;
104
105     /* Print out the results */
106     if (myrank == 1) {
107         /* if t and tst are too different, then there is a performance
108          * problem in the handling of the datatypes */
109
110         if (t > 2 * tst) {
111             errs++;
112             fprintf(stderr,
113                     "Transpose time with datatypes is more than twice time without datatypes\n");
114             fprintf(stderr, "%f\t%f\t%f\n", t, ts, tst);
115         }
116     }
117
118     MPI_Type_free(&col);
119     MPI_Type_free(&xpose);
120
121     MTest_Finalize(errs);
122     MPI_Finalize();
123     return 0;
124 }