Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Change include order for smpi tests/examples to avoid conflicts
[simgrid.git] / teshsuite / smpi / mpich3-test / coll / allred4.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2004 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 #include "mpi.h"
7 #include <stdio.h>
8 #include <stdlib.h>
9 #include "mpitest.h"
10 #include <assert.h>
11
12 /*
13 static char MTEST_Descrip[] = "Test MPI_Allreduce with non-commutative user-defined operations using matrix rotations";
14 */
15
16 /* This example is similar to allred3.c, but uses only 3x3 matrics with 
17    integer-valued entries.  This is an associative but not commutative
18    operation.
19    The number of matrices is the count argument. The matrix is stored
20    in C order, so that
21      c(i,j) is cin[j+i*3]
22
23    Three different matrices are used:
24    I = identity matrix
25    A = (1 0 0    B = (0 1 0
26         0 0 1         1 0 0
27         0 1 0)        0 0 1)
28
29    The product 
30
31          I^k A I^(p-2-k-j) B I^j
32
33    is 
34
35    ( 0 1 0 
36      0 0 1
37      1 0 0 )
38
39    for all values of k, p, and j.  
40  */
41
42 void matmult( void *cinPtr, void *coutPtr, int *count, MPI_Datatype *dtype );
43
44 void matmult( void *cinPtr, void *coutPtr, int *count, MPI_Datatype *dtype )
45 {
46     const int *cin = (const int *)cinPtr;
47     int *cout = (int *)coutPtr;
48     int i, j, k, nmat;
49     int tempcol[3];
50     int offset1, offset2;
51
52     for (nmat = 0; nmat < *count; nmat++) {
53         for (j=0; j<3; j++) {
54             for (i=0; i<3; i++) {
55                 tempcol[i] = 0;
56                 for (k=0; k<3; k++) {
57                     /* col[i] += cin(i,k) * cout(k,j) */
58                     offset1 = k+i*3;
59                     offset2 = j+k*3;
60                     tempcol[i] += cin[offset1] * cout[offset2];
61                 }
62             }
63             for (i=0; i<3; i++) {
64                 offset1 = j+i*3;
65                 cout[offset1] = tempcol[i];
66             }
67         }
68         /* Advance to the next matrix */
69         cin += 9;
70         cout += 9;
71     }
72 }
73
74 /* Initialize the integer matrix as one of the 
75    above matrix entries, as a function of count.
76    We guarantee that both the A and B matrices are included.
77 */   
78 static void initMat( int rank, int size, int nmat, int mat[] )
79 {
80     int i, kind;
81
82     /* Zero the matrix */
83     for (i=0; i<9; i++) {
84         mat[i] = 0;
85     }
86
87     /* Decide which matrix to create (I, A, or B) */
88     if ( size == 2) {
89         /* rank 0 is A, 1 is B */
90         kind = 1 + rank;
91     }
92     else {
93         int tmpA, tmpB;
94         /* Most ranks are identity matrices */
95         kind = 0;
96         /* Make sure exactly one rank gets the A matrix
97            and one the B matrix */
98         tmpA = size / 4;
99         tmpB = (3 * size) / 4;
100         
101         if (rank == tmpA) kind = 1;
102         if (rank == tmpB) kind = 2;
103     }
104     
105     switch (kind) {
106     case 0: /* Identity */
107         mat[0] = 1;
108         mat[4] = 1;
109         mat[8] = 1;
110         break;
111     case 1: /* A */
112         mat[0] = 1;
113         mat[5] = 1;
114         mat[7] = 1;
115         break;
116     case 2: /* B */
117         mat[1] = 1;
118         mat[3] = 1;
119         mat[8] = 1;
120         break;
121     }
122 }
123
124 /* Compare a matrix with the known result */
125 static int checkResult( int nmat, int mat[], const char *msg )
126 {
127     int n, k, errs = 0, wrank;
128     static int solution[9] = { 0, 1, 0, 
129                                0, 0, 1, 
130                                1, 0, 0 };
131
132     MPI_Comm_rank( MPI_COMM_WORLD, &wrank );
133
134     for (n=0; n<nmat; n++) {
135         for (k=0; k<9; k++) {
136             if (mat[k] != solution[k]) {
137                 errs ++;
138                 if (errs == 1) {
139                     printf( "Errors for communicators %s\n", 
140                             MTestGetIntracommName() ); fflush(stdout);
141
142                 }
143                 if (errs < 10) {
144                     printf( "[%d]matrix #%d(%s): Expected mat[%d,%d] = %d, got %d\n",
145                             wrank, n, msg, k / 3, k % 3, solution[k], mat[k] );
146                     fflush(stdout);
147                 }
148             }
149         }
150         /* Advance to the next matrix */
151         mat += 9;
152     }
153     return errs;
154 }
155
156 int main( int argc, char *argv[] )
157 {
158     int errs = 0;
159     int size, rank;
160     int minsize = 2, count; 
161     MPI_Comm      comm;
162     int *buf, *bufout;
163     MPI_Op op;
164     MPI_Datatype mattype;
165     int i;
166
167     MTest_Init( &argc, &argv );
168
169     MPI_Op_create( matmult, 0, &op );
170     
171     /* A single rotation matrix (3x3, stored as 9 consequetive elements) */
172     MPI_Type_contiguous( 9, MPI_INT, &mattype );
173     MPI_Type_commit( &mattype );
174
175     /* Sanity check: test that our routines work properly */
176     { int one = 1;
177     buf = (int *)malloc( 4*9 * sizeof(int) );
178     initMat( 0, 4, 0, &buf[0] );
179     initMat( 1, 4, 0, &buf[9] );
180     initMat( 2, 4, 0, &buf[18] );
181     initMat( 3, 4, 0, &buf[27] );
182     matmult( &buf[0], &buf[9], &one, &mattype );
183     matmult( &buf[9], &buf[18], &one, &mattype );
184     matmult( &buf[18], &buf[27], &one, &mattype );
185     checkResult( 1, &buf[27], "Sanity Check" );
186     free(buf);
187     }
188     
189     while (MTestGetIntracommGeneral( &comm, minsize, 1 )) {
190         if (comm == MPI_COMM_NULL) continue;
191
192         MPI_Comm_size( comm, &size );
193         MPI_Comm_rank( comm, &rank );
194
195         for (count = 1; count < size; count ++ ) {
196             
197             /* Allocate the matrices */
198             buf = (int *)malloc( count * 9 * sizeof(int) );
199             if (!buf) {
200                 MPI_Abort( MPI_COMM_WORLD, 1 );
201                 exit(1);
202             }
203
204             bufout = (int *)malloc( count * 9 * sizeof(int) );
205             if (!bufout) {
206                 MPI_Abort( MPI_COMM_WORLD, 1 );
207                 exit(1);
208             }
209
210             for (i=0; i < count; i++) {
211                 initMat( rank, size, i, &buf[i*9] );
212             }
213             
214             MPI_Allreduce( buf, bufout, count, mattype, op, comm );
215             errs += checkResult( count, bufout, "" );
216
217             /* Try the same test, but using MPI_IN_PLACE */
218             for (i=0; i < count; i++) {
219                 initMat( rank, size, i, &bufout[i*9] );
220             }
221             MPI_Allreduce( MPI_IN_PLACE, bufout, count, mattype, op, comm );
222             errs += checkResult( count, bufout, "IN_PLACE" );
223
224             free( buf );
225             free( bufout );
226         }
227         MTestFreeComm( &comm );
228     }
229         
230     MPI_Op_free( &op );
231     MPI_Type_free( &mattype );
232
233     MTest_Finalize( errs );
234     MPI_Finalize();
235     return 0;
236 }