Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'killgraskill'
[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 three diffrent sizes
48     this is the code that handles these cases
49     if (sendsize < 200 && size > 12) {
50       retval =
51           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
52                                          recvbuf, recvcount, recvtype,
53                                          comm);
54     } else if (sendsize < 3000) {
55       retval =
56           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
57                                                 sendtype, recvbuf,
58                                                 recvcount, recvtype, comm);
59     } else {
60       retval =
61           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
62                                             recvbuf, recvcount, recvtype,
63                                             comm);
64     }
65     
66     
67     */
68     
69     
70     int sizes [3] ={ 4096, 64, 32};
71     for ( j=0 ; j < 3 ; j++ ) {
72       sb = (int *)malloc(size*sizes[j]*sizeof(int));
73       if ( !sb ) {
74         perror( "can't allocate send buffer" );
75         MPI_Abort(MPI_COMM_WORLD,EXIT_FAILURE);
76       }
77       rb = (int *)malloc(size*sizes[j]*sizeof(int));
78       if ( !rb ) {
79         perror( "can't allocate recv buffer");
80         free(sb);
81         MPI_Abort(MPI_COMM_WORLD,EXIT_FAILURE);
82       }
83       for ( i=0 ; i < size*sizes[j] ; ++i ) {
84         sb[i] = rank + 1;
85         rb[i] = 0;
86       }
87
88       /* fputs("Before MPI_Alltoall\n",stdout); */
89       MPI_Barrier(MPI_COMM_WORLD );
90       /* This should really send MPI_CHAR, but since sb and rb were allocated
91        as chunk*size*sizeof(int), the buffers are large enough */
92       status = MPI_Alltoall(sb,sizes[j],MPI_INT,rb,sizes[j],MPI_INT,
93                             MPI_COMM_WORLD);
94
95       /* fputs("Before MPI_Allreduce\n",stdout); */
96       MPI_Allreduce( &status, &gstatus, 1, MPI_INT, MPI_SUM, 
97                     MPI_COMM_WORLD );
98
99       MPI_Barrier(MPI_COMM_WORLD );
100     /* fputs("After MPI_Allreduce\n",stdout); */
101       if (rank == 0 && gstatus != 0) endstatus ++;
102
103       free(sb);
104       free(rb);
105     }
106     
107     if (rank == 0) {
108       if (endstatus == 0) printf( " No Errors\n" );
109       else 
110         printf("all_to_all returned %d erros\n",endstatus);
111     }
112     MPI_Finalize();
113
114     return(EXIT_SUCCESS);
115 }
116