Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
remove unwanted files
[simgrid.git] / teshsuite / smpi / mpich-test / topol / graphtest.c
1 /*
2
3   Test for the MPI Graph routines :
4
5   MPI_Graphdims_get
6   MPI_Graph_create
7   MPI_Graph_get
8   MPI_Graph_map
9   MPI_Graph_neighbors
10   MPI_Graph_neighbors_count
11 */
12
13
14 #include "mpi.h"
15 #include <stdio.h>
16 /* stdlib.h Needed for malloc declaration */
17 #include <stdlib.h>
18 #include "test.h"
19
20 void NumberEdges ( int **, int **, int, int, int );
21 void PrintGraph  ( int, int *, int * );
22
23 int main( int argc, char **argv )
24 {
25     MPI_Comm comm, new_comm;
26     int      reorder;
27     int      nbrarray[3], baseindex;
28     int      size, i, j, nnodes, nedges, q_nnodes, q_nedges, q_nnbrs, newrank;
29     int      *index, *edges, *q_index, *q_edges, *rankbuf;
30     int      worldrank, err = 0, toterr;
31
32     MPI_Init( &argc, &argv );
33
34     MPI_Comm_rank( MPI_COMM_WORLD, &worldrank );
35
36 /* Generate the graph for a binary tree.  
37    
38    Note that EVERY process must have the SAME data
39    */
40     comm = MPI_COMM_WORLD;
41     MPI_Comm_size( comm, &size );
42
43     index = (int *)malloc( (size + 1) * sizeof(int) );
44     edges = (int *)malloc( (size + 1) * 3 * sizeof(int) );
45     reorder = 0;
46     for (i=0; i < size; i++) {
47         index[i] = 0;
48     }
49     NumberEdges( &index, &edges, -1, 0, size - 1 );
50     nedges= index[0];
51     for (i=1; i<size; i++) {
52         nedges += index[i];
53         index[i] = index[i] + index[i-1];
54     }
55     nnodes = size;
56 #ifdef DEBUG
57     PrintGraph( nnodes, index, edges );
58 #endif
59     MPI_Graph_create( comm, nnodes, index, edges, reorder, &new_comm );
60
61 /* Now, try to get the information about this graph */
62     MPI_Graphdims_get( new_comm, &q_nnodes, &q_nedges );
63     if (q_nnodes != nnodes) {
64         printf( "Wrong number of nodes, expected %d got %d\n", nnodes, q_nnodes );
65         err++;
66     }
67     if (q_nedges != nedges) {
68         printf( "Wrong number of edges; expected %d got %d\n", nedges, q_nedges );
69         err++;
70     }
71     q_index = (int *)malloc( q_nnodes * sizeof(int) );
72     q_edges = (int *)malloc( q_nedges * sizeof(int) );
73
74     MPI_Graph_get( new_comm, q_nnodes, q_nedges, q_index, q_edges );
75
76 /* Check with original */
77     if (worldrank == 0) {
78         printf( "Checking graph_get\n" );
79     }
80 /* Because reorder was set to zero, we should have the same data */
81     for (i=0; i<size; i++) {
82         if (index[i] != q_index[i]) {
83             err++;
84             printf( "index[%d] is %d, should be %d\n", i, q_index[i], index[i] );
85         }
86     }
87     for (i=0; i<nedges; i++) {
88         if (edges[i] != q_edges[i]) {
89             err++;
90             printf( "edges[%d] is %d, should be %d\n", i, q_edges[i], edges[i] );
91         }
92     }
93
94 /* Now, get each neighbor set individually */
95     for (i=0; i<size; i++) {
96         MPI_Graph_neighbors_count( new_comm, i, &q_nnbrs );
97         MPI_Graph_neighbors( new_comm, i, 3, nbrarray );
98
99         /* Need to test */
100         baseindex = (i > 0) ? index[i-1] : 0;
101         for (j=0; j<q_nnbrs; j++) {
102             if (nbrarray[j] != edges[baseindex+j]) {
103                 err++;
104                 printf( "nbrarray[%d] for rank %d should be %d, is %d\n",
105                         j, i, edges[baseindex+j], nbrarray[j] );
106             }
107         }
108     }
109
110 /* Test MPI_Graph_map by seeing what ranks are generated for this graph */
111     MPI_Graph_map( comm, nnodes, index, edges, &newrank );
112
113     if (worldrank == 0) {
114         printf( "Checking graph_map\n" );
115     }
116 /* Check that the ranks are at least disjoint among all processors. */
117     rankbuf = (int *)malloc( size * sizeof(int) );
118     MPI_Allgather( &newrank, 1, MPI_INT, rankbuf, 1, MPI_INT, comm );
119     for (i=0; i<size; i++) {
120         for (j=0; j<size; j++) {
121             if (rankbuf[j] == i) break;
122         }
123         if (j >= size) {
124             err++;
125             printf( "Rank %d missing in graph_map\n", i );
126         }
127     }
128
129     MPI_Allreduce( &err, &toterr, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD );
130     if (worldrank == 0) {
131         if (toterr == 0) 
132             printf( "No errors in MPI Graph routines\n" );
133         else
134             printf( "Found %d errors in MPI Graph routines\n", toterr );
135     }
136
137     MPI_Comm_free( &new_comm );
138     free( index );
139     free( edges );
140     free( q_index );
141     free( q_edges );
142     free( rankbuf );
143     MPI_Finalize( );
144     return 0;
145 }
146
147 /* 
148  * Routine to print out a graph for debugging
149  */
150 void PrintGraph( nnodes, index, edges )
151 int nnodes, *index, *edges;
152 {
153     int i, lastidx, j;
154     lastidx=0;
155     printf( "rank\tindex\tedges\n" );
156     for (i=0; i<nnodes; i++) {
157         printf( "%d\t%d\t", i, index[i] );
158         for (j=0; j<index[i] - lastidx; j++) {
159             printf( "%d ", *edges++ );
160         }
161         printf( "\n" );
162         lastidx = index[i];
163     }
164 }
165
166 /* 
167    Number index[0] as first, add its children, and then number them.
168    Note that because of the way the index/edge list is defined, we 
169    need to do a depth-first evaluation
170
171    Each process is connected to the processes rank+1
172    and rank + 1 + floor((size)/2), where size is the size of the subtree 
173
174    Make index[i] the DEGREE of node i.  We'll make the relative to 0
175    at the end.
176  */
177 void NumberEdges( Index, Edges, parent, first, last )
178 int **Index, **Edges, parent, first, last;
179 {
180     int *index = *Index;
181     int *edges = *Edges;
182     int right;
183
184     index[0] = 0;
185     if (parent >= 0) {
186 #ifdef DEBUG    
187         printf( "Adding parent %d to %d\n", parent, first );
188 #endif
189         *index   = *index + 1;
190         *edges++ = parent;
191     }
192     if (first >= last) {
193         /* leaf */
194         index++;
195         if (parent >= 0) {
196             *Index = index;
197             *Edges = edges;
198         }
199         return;
200     }
201
202 /* Internal node.  Always at least a left child */
203 #ifdef DEBUG
204     printf( "Adding left child %d to %d\n", first + 1, first );
205 #endif
206     *index       = *index + 1;
207     *edges++ = first + 1;
208
209 /* Try to add a right child */
210     right = (last - first)/2;
211     right = first + right + 1;
212     if (right == first + 1) 
213         right++;
214     if (right <= last) {
215         /* right child */
216 #ifdef DEBUG    
217         printf( "Adding rightchild %d to %d\n", right, first );
218 #endif
219         *index   = *index + 1;
220         *edges++ = right;
221     }
222     index++;
223     if (first + 1 <= last && right - 1 > first) {
224         NumberEdges( &index, &edges, first, first + 1, 
225                      (right <= last) ? right - 1: last );
226     }
227     if (right <= last) {
228         NumberEdges( &index, &edges, first, right, last );
229     }
230     if (parent >= 0) {
231         *Index = index;
232         *Edges = edges;
233     }
234 }