Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
484ebf736d488d7c35a1521fba21cd2e89b2942b
[simgrid.git] / teshsuite / smpi / mpich3-test / coll / alltoallv.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2001 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 #include "mpi.h"
7 #include "mpitest.h"
8 #include <stdlib.h>
9 #include <stdio.h>
10 #include <string.h>
11
12 /*
13   This program tests MPI_Alltoallv by having processor i send different
14   amounts of data to each processor.
15
16   Because there are separate send and receive types to alltoallv,
17   there need to be tests to rearrange data on the fly.  Not done yet.
18   
19   The first test sends i items to processor i from all processors.
20
21   Currently, the test uses only MPI_INT; this is adequate for testing systems
22   that use point-to-point operations
23  */
24
25 int main( int argc, char **argv )
26 {
27
28     MPI_Comm comm;
29     int      *sbuf, *rbuf;
30     int      rank, size;
31     int      *sendcounts, *recvcounts, *rdispls, *sdispls;
32     int      i, j, *p, err;
33     
34     MTest_Init( &argc, &argv );
35     err = 0;
36     
37     while (MTestGetIntracommGeneral( &comm, 2, 1 )) {
38       if (comm == MPI_COMM_NULL) continue;
39
40       /* Create the buffer */
41       MPI_Comm_size( comm, &size );
42       MPI_Comm_rank( comm, &rank );
43       sbuf = (int *)malloc( size * size * sizeof(int) );
44       rbuf = (int *)malloc( size * size * sizeof(int) );
45       if (!sbuf || !rbuf) {
46         fprintf( stderr, "Could not allocated buffers!\n" );
47         MPI_Abort( comm, 1 );
48         exit(1);
49       }
50       
51       /* Load up the buffers */
52       for (i=0; i<size*size; i++) {
53         sbuf[i] = i + 100*rank;
54         rbuf[i] = -i;
55       }
56       
57       /* Create and load the arguments to alltoallv */
58       sendcounts = (int *)malloc( size * sizeof(int) );
59       recvcounts = (int *)malloc( size * sizeof(int) );
60       rdispls    = (int *)malloc( size * sizeof(int) );
61       sdispls    = (int *)malloc( size * sizeof(int) );
62       if (!sendcounts || !recvcounts || !rdispls || !sdispls) {
63         fprintf( stderr, "Could not allocate arg items!\n" );
64         MPI_Abort( comm, 1 );
65         exit(1);
66       }
67       for (i=0; i<size; i++) {
68         sendcounts[i] = i;
69         recvcounts[i] = rank;
70         rdispls[i]    = i * rank;
71         sdispls[i]    = (i * (i+1))/2;
72       }
73       MPI_Alltoallv( sbuf, sendcounts, sdispls, MPI_INT,
74                      rbuf, recvcounts, rdispls, MPI_INT, comm );
75       
76       /* Check rbuf */
77       for (i=0; i<size; i++) {
78         p = rbuf + rdispls[i];
79         for (j=0; j<rank; j++) {
80           if (p[j] != i * 100 + (rank*(rank+1))/2 + j) {
81             fprintf( stderr, "[%d] got %d expected %d for %dth\n",
82                      rank, p[j],(i*(i+1))/2 + j, j );
83             err++;
84           }
85         }
86       }
87
88       free( sdispls );
89       free( sendcounts );
90       free( sbuf );
91
92 #if MTEST_HAVE_MIN_MPI_VERSION(2,2)
93       /* check MPI_IN_PLACE, added in MPI-2.2 */
94       free( rbuf );
95       rbuf = (int *)malloc( size * (2 * size) * sizeof(int) );
96       if (!rbuf) {
97         fprintf( stderr, "Could not reallocate rbuf!\n" );
98         MPI_Abort( comm, 1 );
99         exit(1);
100       }
101
102       /* Load up the buffers */
103       for (i = 0; i < size; i++) {
104         recvcounts[i] = i + rank;
105         rdispls[i]    = i * (2 * size);
106       }
107       memset(rbuf, -1, size * (2 * size) * sizeof(int));
108       for (i=0; i < size; i++) {
109         p = rbuf + rdispls[i];
110         for (j = 0; j < recvcounts[i]; ++j) {
111           p[j] = 100 * rank + 10 * i + j;
112         }
113       }
114       MPI_Alltoallv( MPI_IN_PLACE, NULL, NULL, MPI_INT,
115                      rbuf, recvcounts, rdispls, MPI_INT, comm );
116       /* Check rbuf */
117       for (i=0; i<size; i++) {
118         p = rbuf + rdispls[i];
119         for (j=0; j<recvcounts[i]; j++) {
120           int expected = 100 * i + 10 * rank + j;
121           if (p[j] != expected) {
122             fprintf(stderr, "[%d] got %d expected %d for block=%d, element=%dth\n",
123                     rank, p[j], expected, i, j);
124             ++err;
125           }
126         }
127       }
128 #endif
129
130       free( rdispls );
131       free( recvcounts );
132       free( rbuf );
133       MTestFreeComm( &comm );
134     }
135
136     MTest_Finalize( err );
137     MPI_Finalize();
138     return 0;
139 }