Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Include directory is in source_dir, not in binary_dir.
[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       }
49       
50       /* Load up the buffers */
51       for (i=0; i<size*size; i++) {
52         sbuf[i] = i + 100*rank;
53         rbuf[i] = -i;
54       }
55       
56       /* Create and load the arguments to alltoallv */
57       sendcounts = (int *)malloc( size * sizeof(int) );
58       recvcounts = (int *)malloc( size * sizeof(int) );
59       rdispls    = (int *)malloc( size * sizeof(int) );
60       sdispls    = (int *)malloc( size * sizeof(int) );
61       if (!sendcounts || !recvcounts || !rdispls || !sdispls) {
62         fprintf( stderr, "Could not allocate arg items!\n" );
63         MPI_Abort( comm, 1 );
64       }
65       for (i=0; i<size; i++) {
66         sendcounts[i] = i;
67         recvcounts[i] = rank;
68         rdispls[i]    = i * rank;
69         sdispls[i]    = (i * (i+1))/2;
70       }
71       MPI_Alltoallv( sbuf, sendcounts, sdispls, MPI_INT,
72                      rbuf, recvcounts, rdispls, MPI_INT, comm );
73       
74       /* Check rbuf */
75       for (i=0; i<size; i++) {
76         p = rbuf + rdispls[i];
77         for (j=0; j<rank; j++) {
78           if (p[j] != i * 100 + (rank*(rank+1))/2 + j) {
79             fprintf( stderr, "[%d] got %d expected %d for %dth\n",
80                      rank, p[j],(i*(i+1))/2 + j, j );
81             err++;
82           }
83         }
84       }
85
86       free( sdispls );
87       free( sendcounts );
88       free( sbuf );
89
90 #if MTEST_HAVE_MIN_MPI_VERSION(2,2)
91       /* check MPI_IN_PLACE, added in MPI-2.2 */
92       free( rbuf );
93       rbuf = (int *)malloc( size * (2 * size) * sizeof(int) );
94       if (!rbuf) {
95         fprintf( stderr, "Could not reallocate rbuf!\n" );
96         MPI_Abort( comm, 1 );
97       }
98
99       /* Load up the buffers */
100       for (i = 0; i < size; i++) {
101         recvcounts[i] = i + rank;
102         rdispls[i]    = i * (2 * size);
103       }
104       memset(rbuf, -1, size * (2 * size) * sizeof(int));
105       for (i=0; i < size; i++) {
106         p = rbuf + rdispls[i];
107         for (j = 0; j < recvcounts[i]; ++j) {
108           p[j] = 100 * rank + 10 * i + j;
109         }
110       }
111       MPI_Alltoallv( MPI_IN_PLACE, NULL, NULL, MPI_INT,
112                      rbuf, recvcounts, rdispls, MPI_INT, comm );
113       /* Check rbuf */
114       for (i=0; i<size; i++) {
115         p = rbuf + rdispls[i];
116         for (j=0; j<recvcounts[i]; j++) {
117           int expected = 100 * i + 10 * rank + j;
118           if (p[j] != expected) {
119             fprintf(stderr, "[%d] got %d expected %d for block=%d, element=%dth\n",
120                     rank, p[j], expected, i, j);
121             ++err;
122           }
123         }
124       }
125 #endif
126
127       free( rdispls );
128       free( recvcounts );
129       free( rbuf );
130       MTestFreeComm( &comm );
131     }
132
133     MTest_Finalize( err );
134     MPI_Finalize();
135     return 0;
136 }