1 /*
2  * Copyright (C) 2010. See COPYRIGHT in top-level directory.
3  */
4 
5 #include <stdio.h>
6 #include <stdlib.h>
7 #include <assert.h>
8 
9 #include <mpi.h>
10 
11 #define INTERCOMM_TAG 0
12 #define NITER 1024
13 
14 
15 void pgroup_create(int grp_size, int *pid_list, MPI_Comm *group_out);
16 void pgroup_free(MPI_Comm *group);
17 
18 
19 /** Free a pgroup
20   */
pgroup_free(MPI_Comm * group)21 void pgroup_free(MPI_Comm *group) {
22   /* Note: It's ok to compare predefined handles */
23   if (*group == MPI_COMM_NULL || *group == MPI_COMM_SELF)
24     return;
25 
26   MPI_Comm_free(group);
27 }
28 
29 
30 /* Create a processor group containing the processes in pid_list.
31  *
32  * NOTE: pid_list list must be identical and sorted on all processes
33  */
pgroup_create(int grp_size,int * pid_list,MPI_Comm * group_out)34 void pgroup_create(int grp_size, int *pid_list, MPI_Comm *group_out) {
35   int       i, grp_me, me, nproc, merge_size;
36   MPI_Comm  pgroup, inter_pgroup;
37 
38   MPI_Comm_rank(MPI_COMM_WORLD, &me);
39   MPI_Comm_size(MPI_COMM_WORLD, &nproc);
40 
41   /* CASE: pgroup size 0 */
42   if (grp_size == 0) {
43     *group_out = MPI_COMM_NULL;
44     return;
45   }
46 
47   /* CASE: pgroup size 1 */
48   else if (grp_size == 1 && pid_list[0] == me) {
49     *group_out = MPI_COMM_SELF;
50     return;
51   }
52 
53   /* CHECK: If I'm not a member, return COMM_NULL */
54   grp_me = -1;
55   for (i = 0; i < grp_size; i++) {
56     if (pid_list[i] == me) {
57       grp_me = i;
58       break;
59     }
60   }
61 
62   if (grp_me < 0) {
63     *group_out = MPI_COMM_NULL;
64     return;
65   }
66 
67   pgroup = MPI_COMM_SELF;
68 
69   for (merge_size = 1; merge_size < grp_size; merge_size *= 2) {
70     int      gid        = grp_me / merge_size;
71     MPI_Comm pgroup_old = pgroup;
72 
73     if (gid % 2 == 0) {
74       /* Check if right partner doesn't exist */
75       if ((gid+1)*merge_size >= grp_size)
76         continue;
77 
78       MPI_Intercomm_create(pgroup, 0, MPI_COMM_WORLD, pid_list[(gid+1)*merge_size], INTERCOMM_TAG, &inter_pgroup);
79       MPI_Intercomm_merge(inter_pgroup, 0 /* LOW */, &pgroup);
80     } else {
81       MPI_Intercomm_create(pgroup, 0, MPI_COMM_WORLD, pid_list[(gid-1)*merge_size], INTERCOMM_TAG, &inter_pgroup);
82       MPI_Intercomm_merge(inter_pgroup, 1 /* HIGH */, &pgroup);
83     }
84 
85     MPI_Comm_free(&inter_pgroup);
86     if (pgroup_old != MPI_COMM_SELF) MPI_Comm_free(&pgroup_old);
87   }
88 
89   *group_out = pgroup;
90 }
91 
92 
main(int argc,char ** argv)93 int main(int argc, char **argv) {
94   int me, nproc, i, j, *glist;
95   MPI_Comm groups[NITER];
96   MPI_Group world_group;
97 
98   MPI_Init(&argc, &argv);
99 
100   MPI_Comm_rank(MPI_COMM_WORLD, &me);
101   MPI_Comm_size(MPI_COMM_WORLD, &nproc);
102   MPI_Comm_group(MPI_COMM_WORLD, &world_group);
103 
104   glist = (int*) malloc(nproc*sizeof(int));
105 
106   for (i = 0; i < nproc; i++)
107     glist[i] = i;
108 
109   if (me == 0)
110     printf("Gsize\tPGgroup (sec)\tComm (sec)\n");
111 
112   MPI_Barrier(MPI_COMM_WORLD);
113 
114   for (i = 1; i <= nproc; i*= 2) {
115     double t_start, t_pg, t_comm;
116     MPI_Group intracomm_group;
117 
118     /** Benchmark pgroup creation cost **/
119 
120     MPI_Barrier(MPI_COMM_WORLD);
121     t_start = MPI_Wtime();
122 
123     for (j = 0; j < NITER; j++)
124       pgroup_create(i, glist, &groups[j]);
125 
126     t_pg = MPI_Wtime() - t_start;
127 
128     for (j = 0; j < NITER; j++)
129       pgroup_free(&groups[j]);
130 
131     /** Benchmark intracommunicator creation cost **/
132 
133     MPI_Group_incl(world_group, i, glist, &intracomm_group);
134     MPI_Barrier(MPI_COMM_WORLD);
135     t_start = MPI_Wtime();
136 
137     for (j = 0; j < NITER; j++)
138       MPI_Comm_create(MPI_COMM_WORLD, intracomm_group, &groups[j]);
139 
140     t_comm = MPI_Wtime() - t_start;
141     MPI_Group_free(&intracomm_group);
142 
143     for (j = 0; j < NITER; j++)
144       pgroup_free(&groups[j]);
145 
146     if (me == 0)
147       printf("%6d\t%0.9f\t%0.9f\n", i, t_pg/NITER, t_comm/NITER);
148 
149   }
150 
151   free(glist);
152   MPI_Group_free(&world_group);
153 
154   MPI_Finalize();
155   return 0;
156 }
157