Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
de425868eb13b250453f9f166d9ac0f5ac9fba3d
[simgrid.git] / teshsuite / smpi / mpich3-test / topo / distgraph1.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *
4  *  (C) 2003 by Argonne National Laboratory.
5  *      See COPYRIGHT in top-level directory.
6  */
7
8 #include "mpi.h"
9 #include <stdio.h>
10 #include <stdlib.h>
11 #include <string.h>
12 #include <assert.h>
13 #include "mpitest.h"
14
15 #define NUM_GRAPHS 10
16 #define MAX_WEIGHT 100
17
18 /* convenience globals */
19 int size, rank;
20
21 /* We need MPI 2.2 to be able to compile the following routines. */
22 #if MTEST_HAVE_MIN_MPI_VERSION(2,2)
23
24 /* Maybe use a bit vector instead? */
25 int **layout;
26
27 #define MAX_LAYOUT_NAME_LEN 256
28 char graph_layout_name[MAX_LAYOUT_NAME_LEN] = {'\0'};
29
30 static void create_graph_layout(int graph_num)
31 {
32     int i, j;
33
34     if (rank == 0) {
35         switch (graph_num) {
36             case 0:
37                 strncpy(graph_layout_name, "deterministic complete graph", MAX_LAYOUT_NAME_LEN);
38                 for (i = 0; i < size; i++)
39                     for (j = 0; j < size; j++)
40                         layout[i][j] = (i + 2) * (j + 1);
41                 break;
42             case 1:
43                 strncpy(graph_layout_name, "every other edge deleted", MAX_LAYOUT_NAME_LEN);
44                 for (i = 0; i < size; i++)
45                     for (j = 0; j < size; j++)
46                         layout[i][j] = (j % 2 ? (i + 2) * (j + 1) : 0);
47                 break;
48             case 2:
49                 strncpy(graph_layout_name, "only self-edges", MAX_LAYOUT_NAME_LEN);
50                 for (i = 0; i < size; i++) {
51                     for (j = 0; j < size; j++) {
52                         if (i == rank && j == rank)
53                             layout[i][j] = 10 * (i + 1);
54                         else
55                             layout[i][j] = 0;
56                     }
57                 }
58                 break;
59             case 3:
60                 strncpy(graph_layout_name, "no edges", MAX_LAYOUT_NAME_LEN);
61                 for (i = 0; i < size; i++)
62                     for (j = 0; j < size; j++)
63                         layout[i][j] = 0;
64                 break;
65             default:
66                 strncpy(graph_layout_name, "a random incomplete graph", MAX_LAYOUT_NAME_LEN);
67                 srand(graph_num);
68
69                 /* Create a connectivity graph; layout[i,j]==w represents an outward
70                  * connectivity from i to j with weight w, w==0 is no edge. */
71                 for (i = 0; i < size; i++) {
72                     for (j = 0; j < size; j++) {
73                         /* disable about a third of the edges */
74                         if (((rand() * 1.0) / RAND_MAX) < 0.33)
75                             layout[i][j] = 0;
76                         else
77                             layout[i][j] = rand() % MAX_WEIGHT;
78                     }
79                 }
80                 break;
81         }
82     }
83
84     /* because of the randomization we must determine the graph on rank 0 and
85      * send the layout to all other processes */
86     MPI_Bcast(graph_layout_name, MAX_LAYOUT_NAME_LEN, MPI_CHAR, 0, MPI_COMM_WORLD);
87     for (i = 0; i < size; ++i) {
88         MPI_Bcast(layout[i], size, MPI_INT, 0, MPI_COMM_WORLD);
89     }
90 }
91
92 static int verify_comm(MPI_Comm comm)
93 {
94     int local_errs = 0;
95     int i, j;
96     int indegree, outdegree, weighted;
97     int *sources, *sweights, *destinations, *dweights;
98     int use_dup;
99     int topo_type = MPI_UNDEFINED;
100     MPI_Comm dupcomm = MPI_COMM_NULL;
101
102     sources = (int *) malloc(size * sizeof(int));
103     sweights = (int *) malloc(size * sizeof(int));
104     destinations = (int *) malloc(size * sizeof(int));
105     dweights = (int *) malloc(size * sizeof(int));
106
107     for (use_dup = 0; use_dup <= 1; ++use_dup) {
108         if (!use_dup) {
109             MPI_Dist_graph_neighbors_count(comm, &indegree, &outdegree, &weighted);
110         }
111         else {
112             MPI_Comm_dup(comm, &dupcomm);
113             comm = dupcomm; /* caller retains original comm value */
114         }
115
116         MPI_Topo_test(comm, &topo_type);
117         if (topo_type != MPI_DIST_GRAPH) {
118             fprintf(stderr, "topo_type != MPI_DIST_GRAPH\n");
119             ++local_errs;
120         }
121
122         j = 0;
123         for (i = 0; i < size; i++)
124             if (layout[i][rank])
125                 j++;
126         if (j != indegree) {
127             fprintf(stderr, "indegree does not match, expected=%d got=%d, layout='%s'\n", indegree, j, graph_layout_name);
128             ++local_errs;
129         }
130
131         j = 0;
132         for (i = 0; i < size; i++)
133             if (layout[rank][i])
134                 j++;
135         if (j != outdegree) {
136             fprintf(stderr, "outdegree does not match, expected=%d got=%d, layout='%s'\n", outdegree, j, graph_layout_name);
137             ++local_errs;
138         }
139
140         if ((indegree || outdegree) && (weighted == 0)) {
141             fprintf(stderr, "MPI_Dist_graph_neighbors_count thinks the graph is not weighted\n");
142             ++local_errs;
143         }
144
145
146         MPI_Dist_graph_neighbors(comm, indegree, sources, sweights, outdegree, destinations, dweights);
147
148         /* For each incoming and outgoing edge in the matrix, search if
149          * the query function listed it in the sources. */
150         for (i = 0; i < size; i++) {
151             if (layout[i][rank]) {
152                 for (j = 0; j < indegree; j++) {
153                     assert(sources[j] >= 0);
154                     assert(sources[j] < size);
155                     if (sources[j] == i)
156                         break;
157                 }
158                 if (j == indegree) {
159                     fprintf(stderr, "no edge from %d to %d specified\n", i, rank);
160                     ++local_errs;
161                 }
162                 else {
163                     if (sweights[j] != layout[i][rank]) {
164                         fprintf(stderr, "incorrect weight for edge (%d,%d): %d instead of %d\n",
165                                 i, rank, sweights[j], layout[i][rank]);
166                         ++local_errs;
167                     }
168                 }
169             }
170             if (layout[rank][i]) {
171                 for (j = 0; j < outdegree; j++) {
172                     assert(destinations[j] >= 0);
173                     assert(destinations[j] < size);
174                     if (destinations[j] == i)
175                         break;
176                 }
177                 if (j == outdegree) {
178                     fprintf(stderr, "no edge from %d to %d specified\n", rank, i);
179                     ++local_errs;
180                 }
181                 else {
182                     if (dweights[j] != layout[rank][i]) {
183                         fprintf(stderr, "incorrect weight for edge (%d,%d): %d instead of %d\n",
184                                 rank, i, dweights[j], layout[rank][i]);
185                         ++local_errs;
186                     }
187                 }
188             }
189         }
190
191         /* For each incoming and outgoing edge in the sources, we should
192          * have an entry in the matrix */
193         for (i = 0; i < indegree; i++) {
194             if (layout[sources[i]][rank] != sweights[i]) {
195                 fprintf(stderr, "edge (%d,%d) has a weight %d instead of %d\n", i, rank,
196                         sweights[i], layout[sources[i]][rank]);
197                 ++local_errs;
198             }
199         }
200         for (i = 0; i < outdegree; i++) {
201             if (layout[rank][destinations[i]] != dweights[i]) {
202                 fprintf(stderr, "edge (%d,%d) has a weight %d instead of %d\n", rank, i,
203                         dweights[i], layout[rank][destinations[i]]);
204                 ++local_errs;
205             }
206         }
207
208     }
209
210     if (dupcomm != MPI_COMM_NULL)
211         MPI_Comm_free(&dupcomm);
212
213     return local_errs;
214 }
215
216 #endif /* At least MPI 2.2 */
217
218 int main(int argc, char *argv[])
219 {
220     int errs = 0;
221     int i, j, k, p;
222     int indegree, outdegree, reorder;
223     int check_indegree, check_outdegree, check_weighted;
224     int *sources, *sweights, *destinations, *dweights, *degrees;
225     MPI_Comm comm;
226
227     MTest_Init(&argc, &argv);
228
229     MPI_Comm_size(MPI_COMM_WORLD, &size);
230     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
231
232 #if MTEST_HAVE_MIN_MPI_VERSION(2,2)
233     layout = (int **) malloc(size * sizeof(int *));
234     assert(layout);
235     for (i = 0; i < size; i++) {
236         layout[i] = (int *) malloc(size * sizeof(int));
237         assert(layout[i]);
238     }
239     /* alloc size*size ints to handle the all-on-one-process case */
240     sources = (int *) malloc(size * size * sizeof(int));
241     sweights = (int *) malloc(size * size * sizeof(int));
242     destinations = (int *) malloc(size * size * sizeof(int));
243     dweights = (int *) malloc(size * size * sizeof(int));
244     degrees = (int *) malloc(size * size * sizeof(int));
245
246     for (i = 0; i < NUM_GRAPHS; i++) {
247         create_graph_layout(i);
248         if (rank == 0) {
249             MTestPrintfMsg( 1, "using graph layout '%s'\n", graph_layout_name );
250         }
251
252         /* MPI_Dist_graph_create_adjacent */
253         if (rank == 0) {
254             MTestPrintfMsg( 1, "testing MPI_Dist_graph_create_adjacent\n" );
255         }
256         indegree = 0;
257         k = 0;
258         for (j = 0; j < size; j++) {
259             if (layout[j][rank]) {
260                 indegree++;
261                 sources[k] = j;
262                 sweights[k++] = layout[j][rank];
263             }
264         }
265
266         outdegree = 0;
267         k = 0;
268         for (j = 0; j < size; j++) {
269             if (layout[rank][j]) {
270                 outdegree++;
271                 destinations[k] = j;
272                 dweights[k++] = layout[rank][j];
273             }
274         }
275
276         for (reorder = 0; reorder <= 1; reorder++) {
277             MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, indegree, sources, sweights,
278                                            outdegree, destinations, dweights, MPI_INFO_NULL,
279                                            reorder, &comm);
280             MPI_Barrier(comm);
281             errs += verify_comm(comm);
282             MPI_Comm_free(&comm);
283         }
284
285         /* a weak check that passing MPI_UNWEIGHTED doesn't cause
286          * create_adjacent to explode */
287         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, indegree, sources, MPI_UNWEIGHTED,
288                                        outdegree, destinations, MPI_UNWEIGHTED, MPI_INFO_NULL,
289                                        reorder, &comm);
290         MPI_Barrier(comm);
291         /* intentionally no verify here, weights won't match */
292         MPI_Comm_free(&comm);
293
294
295         /* MPI_Dist_graph_create() where each process specifies its
296          * outgoing edges */
297         if (rank == 0) {
298             MTestPrintfMsg( 1, 
299                           "testing MPI_Dist_graph_create w/ outgoing only\n" );
300         }
301         sources[0] = rank;
302         k = 0;
303         for (j = 0; j < size; j++) {
304             if (layout[rank][j]) {
305                 destinations[k] = j;
306                 dweights[k++] = layout[rank][j];
307             }
308         }
309         degrees[0] = k;
310         for (reorder = 0; reorder <= 1; reorder++) {
311             MPI_Dist_graph_create(MPI_COMM_WORLD, 1, sources, degrees, destinations, dweights,
312                                   MPI_INFO_NULL, reorder, &comm);
313             MPI_Barrier(comm);
314             errs += verify_comm(comm);
315             MPI_Comm_free(&comm);
316         }
317
318
319         /* MPI_Dist_graph_create() where each process specifies its
320          * incoming edges */
321         if (rank == 0) {
322             MTestPrintfMsg( 1, 
323                          "testing MPI_Dist_graph_create w/ incoming only\n" );
324         }
325         k = 0;
326         for (j = 0; j < size; j++) {
327             if (layout[j][rank]) {
328                 sources[k] = j;
329                 sweights[k] = layout[j][rank];
330                 degrees[k] = 1;
331                 destinations[k++] = rank;
332             }
333         }
334         for (reorder = 0; reorder <= 1; reorder++) {
335             MPI_Dist_graph_create(MPI_COMM_WORLD, k, sources, degrees, destinations, sweights,
336                                   MPI_INFO_NULL, reorder, &comm);
337             MPI_Barrier(comm);
338             errs += verify_comm(comm);
339             MPI_Comm_free(&comm);
340         }
341
342
343         /* MPI_Dist_graph_create() where rank 0 specifies the entire
344          * graph */
345         if (rank == 0) {
346             MTestPrintfMsg( 1, 
347                "testing MPI_Dist_graph_create w/ rank 0 specifies only\n" );
348         }
349         p = 0;
350         for (j = 0; j < size; j++) {
351             for (k = 0; k < size; k++) {
352                 if (layout[j][k]) {
353                     sources[p] = j;
354                     sweights[p] = layout[j][k];
355                     degrees[p] = 1;
356                     destinations[p++] = k;
357                 }
358             }
359         }
360         for (reorder = 0; reorder <= 1; reorder++) {
361             MPI_Dist_graph_create(MPI_COMM_WORLD, (rank == 0) ? p : 0, sources, degrees,
362                                   destinations, sweights, MPI_INFO_NULL, reorder, &comm);
363             MPI_Barrier(comm);
364             errs += verify_comm(comm);
365             MPI_Comm_free(&comm);
366         }
367
368         /* MPI_Dist_graph_create() where rank 0 specifies the entire
369          * graph and all other ranks pass NULL.  Can catch implementation
370          * problems when MPI_UNWEIGHTED==NULL. */
371         if (rank == 0) {
372             MTestPrintfMsg( 1, 
373            "testing MPI_Dist_graph_create w/ rank 0 specifies only -- NULLs\n");
374         }
375         p = 0;
376         for (j = 0; j < size; j++) {
377             for (k = 0; k < size; k++) {
378                 if (layout[j][k]) {
379                     sources[p] = j;
380                     sweights[p] = layout[j][k];
381                     degrees[p] = 1;
382                     destinations[p++] = k;
383                 }
384             }
385         }
386         for (reorder = 0; reorder <= 1; reorder++) {
387             if (rank == 0) {
388                 MPI_Dist_graph_create(MPI_COMM_WORLD, p, sources, degrees,
389                                       destinations, sweights, MPI_INFO_NULL, reorder, &comm);
390             }
391             else {
392                 MPI_Dist_graph_create(MPI_COMM_WORLD, 0, NULL, NULL,
393                                       NULL, NULL, MPI_INFO_NULL, reorder, &comm);
394             }
395             MPI_Barrier(comm);
396             errs += verify_comm(comm);
397             MPI_Comm_free(&comm);
398         }
399
400     }
401
402     /* now tests that don't depend on the layout[][] array */
403
404     /* The MPI-2.2 standard recommends implementations set
405      * MPI_UNWEIGHTED==NULL, but this leads to an ambiguity.  The draft
406      * MPI-3.0 standard specifically recommends _not_ setting it equal
407      * to NULL. */
408     if (MPI_UNWEIGHTED == NULL) {
409         fprintf(stderr, "MPI_UNWEIGHTED should not be NULL\n");
410         ++errs;
411     }
412
413     /* MPI_Dist_graph_create() with no graph */
414     if (rank == 0) {
415         MTestPrintfMsg( 1, "testing MPI_Dist_graph_create w/ no graph\n" );
416     }
417     for (reorder = 0; reorder <= 1; reorder++) {
418         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, sources, degrees,
419                               destinations, sweights, MPI_INFO_NULL, reorder, &comm);
420         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
421         if (!check_weighted) {
422             fprintf(stderr, "expected weighted == TRUE for the \"no graph\" case\n");
423             ++errs;
424         }
425         MPI_Comm_free(&comm);
426     }
427
428     /* MPI_Dist_graph_create() with no graph -- passing MPI_WEIGHTS_EMPTY 
429        instead */
430     /* NOTE that MPI_WEIGHTS_EMPTY was added in MPI-3 and does not 
431        appear before then.  This part of the test thus requires a check
432        on the MPI major version */
433 #if MPI_VERSION >= 3
434     if (rank == 0) {
435         MTestPrintfMsg( 1, "testing MPI_Dist_graph_create w/ no graph\n" );
436     }
437     for (reorder = 0; reorder <= 1; reorder++) {
438         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, sources, degrees,
439                               destinations, MPI_WEIGHTS_EMPTY, MPI_INFO_NULL, reorder, &comm);
440         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
441         if (!check_weighted) {
442             fprintf(stderr, "expected weighted == TRUE for the \"no graph -- MPI_WEIGHTS_EMPTY\" case\n");
443             ++errs;
444         }
445         MPI_Comm_free(&comm);
446     }
447 #endif
448
449     /* MPI_Dist_graph_create() with no graph -- passing NULLs instead */
450     if (rank == 0) {
451         MTestPrintfMsg( 1, 
452                       "testing MPI_Dist_graph_create w/ no graph -- NULLs\n" );
453     }
454     for (reorder = 0; reorder <= 1; reorder++) {
455         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, NULL, NULL,
456                               NULL, NULL, MPI_INFO_NULL, reorder, &comm);
457         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
458         /* ambiguous if they are equal, only check when they are distinct values. */
459         if (MPI_UNWEIGHTED != NULL) {
460             if (!check_weighted) {
461                 fprintf(stderr, "expected weighted == TRUE for the \"no graph -- NULLs\" case\n");
462                 ++errs;
463             }
464         }
465         MPI_Comm_free(&comm);
466     }
467
468     /* MPI_Dist_graph_create() with no graph -- passing NULLs+MPI_UNWEIGHTED instead */
469     if (rank == 0) {
470         MTestPrintfMsg( 1, 
471         "testing MPI_Dist_graph_create w/ no graph -- NULLs+MPI_UNWEIGHTED\n" );
472     }
473     for (reorder = 0; reorder <= 1; reorder++) {
474         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, NULL, NULL,
475                               NULL, MPI_UNWEIGHTED, MPI_INFO_NULL, reorder, &comm);
476         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
477         /* ambiguous if they are equal, only check when they are distinct values. */
478         if (MPI_UNWEIGHTED != NULL) {
479             if (check_weighted) {
480                 fprintf(stderr, "expected weighted == FALSE for the \"no graph -- NULLs+MPI_UNWEIGHTED\" case\n");
481                 ++errs;
482             }
483         }
484         MPI_Comm_free(&comm);
485     }
486
487     /* MPI_Dist_graph_create_adjacent() with no graph */
488     if (rank == 0) {
489         MTestPrintfMsg( 1, 
490                      "testing MPI_Dist_graph_create_adjacent w/ no graph\n" );
491     }
492     for (reorder = 0; reorder <= 1; reorder++) {
493         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, sources, sweights,
494                               0, destinations, dweights, MPI_INFO_NULL, reorder, &comm);
495         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
496         if (!check_weighted) {
497             fprintf(stderr, "expected weighted == TRUE for the \"no graph\" case\n");
498             ++errs;
499         }
500         MPI_Comm_free(&comm);
501     }
502
503     /* MPI_Dist_graph_create_adjacent() with no graph -- passing MPI_WEIGHTS_EMPTY instead */
504     /* NOTE that MPI_WEIGHTS_EMPTY was added in MPI-3 and does not 
505        appear before then.  This part of the test thus requires a check
506        on the MPI major version */
507 #if MPI_VERSION >= 3
508     if (rank == 0) {
509         MTestPrintfMsg( 1, 
510   "testing MPI_Dist_graph_create_adjacent w/ no graph -- MPI_WEIGHTS_EMPTY\n" );
511     }
512     for (reorder = 0; reorder <= 1; reorder++) {
513         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, sources, MPI_WEIGHTS_EMPTY,
514                               0, destinations, MPI_WEIGHTS_EMPTY, MPI_INFO_NULL, reorder, &comm);
515         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
516         if (!check_weighted) {
517             fprintf(stderr, "expected weighted == TRUE for the \"no graph -- MPI_WEIGHTS_EMPTY\" case\n");
518             ++errs;
519         }
520         MPI_Comm_free(&comm);
521     }
522 #endif
523
524     /* MPI_Dist_graph_create_adjacent() with no graph -- passing NULLs instead */
525     if (rank == 0) {
526         MTestPrintfMsg( 1, 
527               "testing MPI_Dist_graph_create_adjacent w/ no graph -- NULLs\n" );
528     }
529     for (reorder = 0; reorder <= 1; reorder++) {
530         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, NULL, NULL,
531                               0, NULL, NULL, MPI_INFO_NULL, reorder, &comm);
532         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
533         /* ambiguous if they are equal, only check when they are distinct values. */
534         if (MPI_UNWEIGHTED != NULL) {
535             if (!check_weighted) {
536                 fprintf(stderr, "expected weighted == TRUE for the \"no graph -- NULLs\" case\n");
537                 ++errs;
538             }
539         }
540         MPI_Comm_free(&comm);
541     }
542
543     /* MPI_Dist_graph_create_adjacent() with no graph -- passing NULLs+MPI_UNWEIGHTED instead */
544     if (rank == 0) {
545         MTestPrintfMsg( 1, 
546 "testing MPI_Dist_graph_create_adjacent w/ no graph -- NULLs+MPI_UNWEIGHTED\n");
547     }
548     for (reorder = 0; reorder <= 1; reorder++) {
549         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, NULL, MPI_UNWEIGHTED,
550                               0, NULL, MPI_UNWEIGHTED, MPI_INFO_NULL, reorder, &comm);
551         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
552         /* ambiguous if they are equal, only check when they are distinct values. */
553         if (MPI_UNWEIGHTED != NULL) {
554             if (check_weighted) {
555                 fprintf(stderr, "expected weighted == FALSE for the \"no graph -- NULLs+MPI_UNWEIGHTED\" case\n");
556                 ++errs;
557             }
558         }
559         MPI_Comm_free(&comm);
560     }
561
562
563     for (i = 0; i < size; i++)
564         free(layout[i]);
565     free(layout);
566 #endif
567
568     MTest_Finalize(errs);
569     MPI_Finalize();
570
571     return 0;
572 }