Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
02c92a90929181b7015d56cb7e98845abcc2ea13
[simgrid.git] / teshsuite / smpi / mpich-test / coll / coll13.c
1 #include "mpi.h"
2 #include "test.h"
3
4 /* 
5 From: hook@nas.nasa.gov (Edward C. Hook)
6  */
7
8 #include <stdlib.h>
9 #include <stdio.h>
10
11 #include <string.h>
12 #include <errno.h>
13 #ifndef EXIT_SUCCESS
14 #define EXIT_SUCCESS 0
15 #define EXIT_FAILURE 1
16 #endif
17
18 int main( int argc, char *argv[] )
19 {
20     int rank, size;
21     int chunk = 4096;
22     int i,j;
23     int *sb;
24     int *rb;
25     int status, gstatus, endstatus;
26     endstatus=0;
27     MPI_Init(&argc,&argv);
28     MPI_Comm_rank(MPI_COMM_WORLD,&rank);
29     MPI_Comm_size(MPI_COMM_WORLD,&size);
30
31     for ( i=1 ; i < argc ; ++i ) {
32         if ( argv[i][0] != '-' )
33             continue;
34         switch(argv[i][1]) {
35         case 'm':
36             chunk = atoi(argv[++i]);
37             break;
38         default:
39             fprintf(stderr,"Unrecognized argument %s\n",
40                     argv[i]);
41             MPI_Abort(MPI_COMM_WORLD,EXIT_FAILURE);
42         }
43     }
44      
45     
46     /*
47     SMPI addition: we want to test all three alltoall algorithms, thus we use
48     three different sizes.  This is the code that handles these cases:
49
50     if (sendsize < 200 && size > 12) {
51       retval =
52           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
53                                          recvbuf, recvcount, recvtype,
54                                          comm);
55     } else if (sendsize < 3000) {
56       retval =
57           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
58                                                 sendtype, recvbuf,
59                                                 recvcount, recvtype, comm);
60     } else {
61       retval =
62           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
63                                             recvbuf, recvcount, recvtype,
64                                             comm);
65     }
66     */
67     
68     
69     int sizes [3] ={ 4096, 64, 32};
70     for ( j=0 ; j < 3 ; j++ ) {
71       sb = (int *)malloc(size*sizes[j]*sizeof(int));
72       if ( !sb ) {
73         perror( "can't allocate send buffer" );
74         MPI_Abort(MPI_COMM_WORLD,EXIT_FAILURE);
75       }
76       rb = (int *)malloc(size*sizes[j]*sizeof(int));
77       if ( !rb ) {
78         perror( "can't allocate recv buffer");
79         free(sb);
80         MPI_Abort(MPI_COMM_WORLD,EXIT_FAILURE);
81       }
82       for ( i=0 ; i < size*sizes[j] ; ++i ) {
83         sb[i] = rank + 1;
84         rb[i] = 0;
85       }
86
87       /* fputs("Before MPI_Alltoall\n",stdout); */
88       MPI_Barrier(MPI_COMM_WORLD );
89       /* This should really send MPI_CHAR, but since sb and rb were allocated
90        as chunk*size*sizeof(int), the buffers are large enough */
91       status = MPI_Alltoall(sb,sizes[j],MPI_INT,rb,sizes[j],MPI_INT,
92                             MPI_COMM_WORLD);
93
94       /* fputs("Before MPI_Allreduce\n",stdout); */
95       MPI_Allreduce( &status, &gstatus, 1, MPI_INT, MPI_SUM, 
96                     MPI_COMM_WORLD );
97
98       MPI_Barrier(MPI_COMM_WORLD );
99     /* fputs("After MPI_Allreduce\n",stdout); */
100       if (rank == 0 && gstatus != 0) endstatus ++;
101
102       free(sb);
103       free(rb);
104     }
105     
106     if (rank == 0) {
107       if (endstatus == 0) printf( " No Errors\n" );
108       else 
109         printf("all_to_all returned %d erros\n",endstatus);
110     }
111     MPI_Finalize();
112
113     return(EXIT_SUCCESS);
114 }
115