Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update topo
[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,
128                     j, graph_layout_name);
129             ++local_errs;
130         }
131
132         j = 0;
133         for (i = 0; i < size; i++)
134             if (layout[rank][i])
135                 j++;
136         if (j != outdegree) {
137             fprintf(stderr, "outdegree does not match, expected=%d got=%d, layout='%s'\n",
138                     outdegree, j, graph_layout_name);
139             ++local_errs;
140         }
141
142         if ((indegree || outdegree) && (weighted == 0)) {
143             fprintf(stderr, "MPI_Dist_graph_neighbors_count thinks the graph is not weighted\n");
144             ++local_errs;
145         }
146
147
148         MPI_Dist_graph_neighbors(comm, indegree, sources, sweights, outdegree, destinations,
149                                  dweights);
150
151         /* For each incoming and outgoing edge in the matrix, search if
152          * the query function listed it in the sources. */
153         for (i = 0; i < size; i++) {
154             if (layout[i][rank]) {
155                 for (j = 0; j < indegree; j++) {
156                     assert(sources[j] >= 0);
157                     assert(sources[j] < size);
158                     if (sources[j] == i)
159                         break;
160                 }
161                 if (j == indegree) {
162                     fprintf(stderr, "no edge from %d to %d specified\n", i, rank);
163                     ++local_errs;
164                 }
165                 else {
166                     if (sweights[j] != layout[i][rank]) {
167                         fprintf(stderr, "incorrect weight for edge (%d,%d): %d instead of %d\n",
168                                 i, rank, sweights[j], layout[i][rank]);
169                         ++local_errs;
170                     }
171                 }
172             }
173             if (layout[rank][i]) {
174                 for (j = 0; j < outdegree; j++) {
175                     assert(destinations[j] >= 0);
176                     assert(destinations[j] < size);
177                     if (destinations[j] == i)
178                         break;
179                 }
180                 if (j == outdegree) {
181                     fprintf(stderr, "no edge from %d to %d specified\n", rank, i);
182                     ++local_errs;
183                 }
184                 else {
185                     if (dweights[j] != layout[rank][i]) {
186                         fprintf(stderr, "incorrect weight for edge (%d,%d): %d instead of %d\n",
187                                 rank, i, dweights[j], layout[rank][i]);
188                         ++local_errs;
189                     }
190                 }
191             }
192         }
193
194         /* For each incoming and outgoing edge in the sources, we should
195          * have an entry in the matrix */
196         for (i = 0; i < indegree; i++) {
197             if (layout[sources[i]][rank] != sweights[i]) {
198                 fprintf(stderr, "edge (%d,%d) has a weight %d instead of %d\n", i, rank,
199                         sweights[i], layout[sources[i]][rank]);
200                 ++local_errs;
201             }
202         }
203         for (i = 0; i < outdegree; i++) {
204             if (layout[rank][destinations[i]] != dweights[i]) {
205                 fprintf(stderr, "edge (%d,%d) has a weight %d instead of %d\n", rank, i,
206                         dweights[i], layout[rank][destinations[i]]);
207                 ++local_errs;
208             }
209         }
210
211     }
212
213     if (dupcomm != MPI_COMM_NULL)
214         MPI_Comm_free(&dupcomm);
215
216     free(sources);
217     free(sweights);
218     free(destinations);
219     free(dweights);
220     return local_errs;
221 }
222
223 #endif /* At least MPI 2.2 */
224
225 int main(int argc, char *argv[])
226 {
227     int errs = 0;
228     int i, j, k, p;
229     int indegree, outdegree, reorder;
230     int check_indegree, check_outdegree, check_weighted;
231     int *sources, *sweights, *destinations, *dweights, *degrees;
232     MPI_Comm comm;
233
234     MTest_Init(&argc, &argv);
235
236     MPI_Comm_size(MPI_COMM_WORLD, &size);
237     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
238
239 #if MTEST_HAVE_MIN_MPI_VERSION(2,2)
240     layout = (int **) malloc(size * sizeof(int *));
241     assert(layout);
242     for (i = 0; i < size; i++) {
243         layout[i] = (int *) malloc(size * sizeof(int));
244         assert(layout[i]);
245     }
246     /* alloc size*size ints to handle the all-on-one-process case */
247     sources = (int *) malloc(size * size * sizeof(int));
248     sweights = (int *) malloc(size * size * sizeof(int));
249     destinations = (int *) malloc(size * size * sizeof(int));
250     dweights = (int *) malloc(size * size * sizeof(int));
251     degrees = (int *) malloc(size * size * sizeof(int));
252
253     for (i = 0; i < NUM_GRAPHS; i++) {
254         create_graph_layout(i);
255         if (rank == 0) {
256             MTestPrintfMsg(1, "using graph layout '%s'\n", graph_layout_name);
257         }
258
259         /* MPI_Dist_graph_create_adjacent */
260         if (rank == 0) {
261             MTestPrintfMsg(1, "testing MPI_Dist_graph_create_adjacent\n");
262         }
263         indegree = 0;
264         k = 0;
265         for (j = 0; j < size; j++) {
266             if (layout[j][rank]) {
267                 indegree++;
268                 sources[k] = j;
269                 sweights[k++] = layout[j][rank];
270             }
271         }
272
273         outdegree = 0;
274         k = 0;
275         for (j = 0; j < size; j++) {
276             if (layout[rank][j]) {
277                 outdegree++;
278                 destinations[k] = j;
279                 dweights[k++] = layout[rank][j];
280             }
281         }
282
283         for (reorder = 0; reorder <= 1; reorder++) {
284             MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, indegree, sources, sweights,
285                                            outdegree, destinations, dweights, MPI_INFO_NULL,
286                                            reorder, &comm);
287             MPI_Barrier(comm);
288             errs += verify_comm(comm);
289             MPI_Comm_free(&comm);
290         }
291
292         /* a weak check that passing MPI_UNWEIGHTED doesn't cause
293          * create_adjacent to explode */
294         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, indegree, sources, MPI_UNWEIGHTED,
295                                        outdegree, destinations, MPI_UNWEIGHTED, MPI_INFO_NULL,
296                                        reorder, &comm);
297         MPI_Barrier(comm);
298         /* intentionally no verify here, weights won't match */
299         MPI_Comm_free(&comm);
300
301
302         /* MPI_Dist_graph_create() where each process specifies its
303          * outgoing edges */
304         if (rank == 0) {
305             MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ outgoing only\n");
306         }
307         sources[0] = rank;
308         k = 0;
309         for (j = 0; j < size; j++) {
310             if (layout[rank][j]) {
311                 destinations[k] = j;
312                 dweights[k++] = layout[rank][j];
313             }
314         }
315         degrees[0] = k;
316         for (reorder = 0; reorder <= 1; reorder++) {
317             MPI_Dist_graph_create(MPI_COMM_WORLD, 1, sources, degrees, destinations, dweights,
318                                   MPI_INFO_NULL, reorder, &comm);
319             MPI_Barrier(comm);
320             errs += verify_comm(comm);
321             MPI_Comm_free(&comm);
322         }
323
324
325         /* MPI_Dist_graph_create() where each process specifies its
326          * incoming edges */
327         if (rank == 0) {
328             MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ incoming only\n");
329         }
330         k = 0;
331         for (j = 0; j < size; j++) {
332             if (layout[j][rank]) {
333                 sources[k] = j;
334                 sweights[k] = layout[j][rank];
335                 degrees[k] = 1;
336                 destinations[k++] = rank;
337             }
338         }
339         for (reorder = 0; reorder <= 1; reorder++) {
340             MPI_Dist_graph_create(MPI_COMM_WORLD, k, sources, degrees, destinations, sweights,
341                                   MPI_INFO_NULL, reorder, &comm);
342             MPI_Barrier(comm);
343             errs += verify_comm(comm);
344             MPI_Comm_free(&comm);
345         }
346
347
348         /* MPI_Dist_graph_create() where rank 0 specifies the entire
349          * graph */
350         if (rank == 0) {
351             MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ rank 0 specifies only\n");
352         }
353         p = 0;
354         for (j = 0; j < size; j++) {
355             for (k = 0; k < size; k++) {
356                 if (layout[j][k]) {
357                     sources[p] = j;
358                     sweights[p] = layout[j][k];
359                     degrees[p] = 1;
360                     destinations[p++] = k;
361                 }
362             }
363         }
364         for (reorder = 0; reorder <= 1; reorder++) {
365             MPI_Dist_graph_create(MPI_COMM_WORLD, (rank == 0) ? p : 0, sources, degrees,
366                                   destinations, sweights, MPI_INFO_NULL, reorder, &comm);
367             MPI_Barrier(comm);
368             errs += verify_comm(comm);
369             MPI_Comm_free(&comm);
370         }
371
372         /* MPI_Dist_graph_create() where rank 0 specifies the entire
373          * graph and all other ranks pass NULL.  Can catch implementation
374          * problems when MPI_UNWEIGHTED==NULL. */
375         if (rank == 0) {
376             MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ rank 0 specifies only -- NULLs\n");
377         }
378         p = 0;
379         for (j = 0; j < size; j++) {
380             for (k = 0; k < size; k++) {
381                 if (layout[j][k]) {
382                     sources[p] = j;
383                     sweights[p] = layout[j][k];
384                     degrees[p] = 1;
385                     destinations[p++] = k;
386                 }
387             }
388         }
389         for (reorder = 0; reorder <= 1; reorder++) {
390             if (rank == 0) {
391                 MPI_Dist_graph_create(MPI_COMM_WORLD, p, sources, degrees,
392                                       destinations, sweights, MPI_INFO_NULL, reorder, &comm);
393             }
394             else {
395                 MPI_Dist_graph_create(MPI_COMM_WORLD, 0, NULL, NULL,
396                                       NULL, NULL, MPI_INFO_NULL, reorder, &comm);
397             }
398             MPI_Barrier(comm);
399             errs += verify_comm(comm);
400             MPI_Comm_free(&comm);
401         }
402
403     }
404
405     /* now tests that don't depend on the layout[][] array */
406
407     /* The MPI-2.2 standard recommends implementations set
408      * MPI_UNWEIGHTED==NULL, but this leads to an ambiguity.  The draft
409      * MPI-3.0 standard specifically recommends _not_ setting it equal
410      * to NULL. */
411     if (MPI_UNWEIGHTED == NULL) {
412         fprintf(stderr, "MPI_UNWEIGHTED should not be NULL\n");
413         ++errs;
414     }
415
416     /* MPI_Dist_graph_create() with no graph */
417     if (rank == 0) {
418         MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ no graph\n");
419     }
420     for (reorder = 0; reorder <= 1; reorder++) {
421         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, sources, degrees,
422                               destinations, sweights, MPI_INFO_NULL, reorder, &comm);
423         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
424         if (!check_weighted) {
425             fprintf(stderr, "expected weighted == TRUE for the \"no graph\" case\n");
426             ++errs;
427         }
428         MPI_Comm_free(&comm);
429     }
430
431     /* MPI_Dist_graph_create() with no graph -- passing MPI_WEIGHTS_EMPTY
432      * instead */
433     /* NOTE that MPI_WEIGHTS_EMPTY was added in MPI-3 and does not
434      * appear before then.  This part of the test thus requires a check
435      * on the MPI major version */
436 #if MPI_VERSION >= 3
437     if (rank == 0) {
438         MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ no graph\n");
439     }
440     for (reorder = 0; reorder <= 1; reorder++) {
441         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, sources, degrees,
442                               destinations, MPI_WEIGHTS_EMPTY, MPI_INFO_NULL, reorder, &comm);
443         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
444         if (!check_weighted) {
445             fprintf(stderr,
446                     "expected weighted == TRUE for the \"no graph -- MPI_WEIGHTS_EMPTY\" case\n");
447             ++errs;
448         }
449         MPI_Comm_free(&comm);
450     }
451 #endif
452
453     /* MPI_Dist_graph_create() with no graph -- passing NULLs instead */
454     if (rank == 0) {
455         MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ no graph -- NULLs\n");
456     }
457     for (reorder = 0; reorder <= 1; reorder++) {
458         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, NULL, NULL,
459                               NULL, NULL, MPI_INFO_NULL, reorder, &comm);
460         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
461         /* ambiguous if they are equal, only check when they are distinct values. */
462         if (MPI_UNWEIGHTED != NULL) {
463             if (!check_weighted) {
464                 fprintf(stderr, "expected weighted == TRUE for the \"no graph -- NULLs\" case\n");
465                 ++errs;
466             }
467         }
468         MPI_Comm_free(&comm);
469     }
470
471     /* MPI_Dist_graph_create() with no graph -- passing NULLs+MPI_UNWEIGHTED instead */
472     if (rank == 0) {
473         MTestPrintfMsg(1, "testing MPI_Dist_graph_create w/ no graph -- NULLs+MPI_UNWEIGHTED\n");
474     }
475     for (reorder = 0; reorder <= 1; reorder++) {
476         MPI_Dist_graph_create(MPI_COMM_WORLD, 0, NULL, NULL,
477                               NULL, MPI_UNWEIGHTED, MPI_INFO_NULL, reorder, &comm);
478         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
479         /* ambiguous if they are equal, only check when they are distinct values. */
480         if (MPI_UNWEIGHTED != NULL) {
481             if (check_weighted) {
482                 fprintf(stderr,
483                         "expected weighted == FALSE for the \"no graph -- NULLs+MPI_UNWEIGHTED\" case\n");
484                 ++errs;
485             }
486         }
487         MPI_Comm_free(&comm);
488     }
489
490     /* MPI_Dist_graph_create_adjacent() with no graph */
491     if (rank == 0) {
492         MTestPrintfMsg(1, "testing MPI_Dist_graph_create_adjacent w/ no graph\n");
493     }
494     for (reorder = 0; reorder <= 1; reorder++) {
495         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, sources, sweights,
496                                        0, destinations, dweights, MPI_INFO_NULL, reorder, &comm);
497         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
498         if (!check_weighted) {
499             fprintf(stderr, "expected weighted == TRUE for the \"no graph\" case\n");
500             ++errs;
501         }
502         MPI_Comm_free(&comm);
503     }
504
505     /* MPI_Dist_graph_create_adjacent() with no graph -- passing MPI_WEIGHTS_EMPTY instead */
506     /* NOTE that MPI_WEIGHTS_EMPTY was added in MPI-3 and does not
507      * appear before then.  This part of the test thus requires a check
508      * on the MPI major version */
509 #if MPI_VERSION >= 3
510     if (rank == 0) {
511         MTestPrintfMsg(1,
512                        "testing MPI_Dist_graph_create_adjacent w/ no graph -- MPI_WEIGHTS_EMPTY\n");
513     }
514     for (reorder = 0; reorder <= 1; reorder++) {
515         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, sources, MPI_WEIGHTS_EMPTY,
516                                        0, destinations, MPI_WEIGHTS_EMPTY, MPI_INFO_NULL, reorder,
517                                        &comm);
518         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
519         if (!check_weighted) {
520             fprintf(stderr,
521                     "expected weighted == TRUE for the \"no graph -- MPI_WEIGHTS_EMPTY\" case\n");
522             ++errs;
523         }
524         MPI_Comm_free(&comm);
525     }
526 #endif
527
528     /* MPI_Dist_graph_create_adjacent() with no graph -- passing NULLs instead */
529     if (rank == 0) {
530         MTestPrintfMsg(1, "testing MPI_Dist_graph_create_adjacent w/ no graph -- NULLs\n");
531     }
532     for (reorder = 0; reorder <= 1; reorder++) {
533         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, NULL, NULL,
534                                        0, NULL, NULL, MPI_INFO_NULL, reorder, &comm);
535         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
536         /* ambiguous if they are equal, only check when they are distinct values. */
537         if (MPI_UNWEIGHTED != NULL) {
538             if (!check_weighted) {
539                 fprintf(stderr, "expected weighted == TRUE for the \"no graph -- NULLs\" case\n");
540                 ++errs;
541             }
542         }
543         MPI_Comm_free(&comm);
544     }
545
546     /* MPI_Dist_graph_create_adjacent() with no graph -- passing NULLs+MPI_UNWEIGHTED instead */
547     if (rank == 0) {
548         MTestPrintfMsg(1,
549                        "testing MPI_Dist_graph_create_adjacent w/ no graph -- NULLs+MPI_UNWEIGHTED\n");
550     }
551     for (reorder = 0; reorder <= 1; reorder++) {
552         MPI_Dist_graph_create_adjacent(MPI_COMM_WORLD, 0, NULL, MPI_UNWEIGHTED,
553                                        0, NULL, MPI_UNWEIGHTED, MPI_INFO_NULL, reorder, &comm);
554         MPI_Dist_graph_neighbors_count(comm, &check_indegree, &check_outdegree, &check_weighted);
555         /* ambiguous if they are equal, only check when they are distinct values. */
556         if (MPI_UNWEIGHTED != NULL) {
557             if (check_weighted) {
558                 fprintf(stderr,
559                         "expected weighted == FALSE for the \"no graph -- NULLs+MPI_UNWEIGHTED\" case\n");
560                 ++errs;
561             }
562         }
563         MPI_Comm_free(&comm);
564     }
565
566
567     for (i = 0; i < size; i++)
568         free(layout[i]);
569     free(layout);
570     free(sources);
571     free(sweights);
572     free(destinations);
573     free(dweights);
574     free(degrees);
575 #endif
576
577     MTest_Finalize(errs);
578     MPI_Finalize();
579
580     return 0;
581 }