File: pgroup_bench.c

package info (click to toggle)
armci-mpi 0.0~git20160222-2
  • links: PTS
  • area: main
  • in suites: stretch
  • size: 1,756 kB
  • sloc: ansic: 12,698; sh: 229; makefile: 53; fortran: 44
file content (156 lines) | stat: -rw-r--r-- 3,599 bytes parent folder | download | duplicates (6)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
/*
 * Copyright (C) 2010. See COPYRIGHT in top-level directory.
 */

#include <stdio.h>
#include <stdlib.h>
#include <assert.h>

#include <mpi.h>

#define INTERCOMM_TAG 0
#define NITER 1024


void pgroup_create(int grp_size, int *pid_list, MPI_Comm *group_out);
void pgroup_free(MPI_Comm *group);


/** Free a pgroup
  */
void pgroup_free(MPI_Comm *group) {
  /* Note: It's ok to compare predefined handles */
  if (*group == MPI_COMM_NULL || *group == MPI_COMM_SELF)
    return;

  MPI_Comm_free(group);
}


/* Create a processor group containing the processes in pid_list.
 *
 * NOTE: pid_list list must be identical and sorted on all processes
 */
void pgroup_create(int grp_size, int *pid_list, MPI_Comm *group_out) {
  int       i, grp_me, me, nproc, merge_size;
  MPI_Comm  pgroup, inter_pgroup;

  MPI_Comm_rank(MPI_COMM_WORLD, &me);
  MPI_Comm_size(MPI_COMM_WORLD, &nproc);

  /* CASE: pgroup size 0 */
  if (grp_size == 0) {
    *group_out = MPI_COMM_NULL;
    return;
  }

  /* CASE: pgroup size 1 */
  else if (grp_size == 1 && pid_list[0] == me) {
    *group_out = MPI_COMM_SELF;
    return;
  }

  /* CHECK: If I'm not a member, return COMM_NULL */
  grp_me = -1;
  for (i = 0; i < grp_size; i++) {
    if (pid_list[i] == me) {
      grp_me = i;
      break;
    }
  }

  if (grp_me < 0) {
    *group_out = MPI_COMM_NULL;
    return;
  }

  pgroup = MPI_COMM_SELF;

  for (merge_size = 1; merge_size < grp_size; merge_size *= 2) {
    int      gid        = grp_me / merge_size;
    MPI_Comm pgroup_old = pgroup;

    if (gid % 2 == 0) {
      /* Check if right partner doesn't exist */
      if ((gid+1)*merge_size >= grp_size)
        continue;

      MPI_Intercomm_create(pgroup, 0, MPI_COMM_WORLD, pid_list[(gid+1)*merge_size], INTERCOMM_TAG, &inter_pgroup);
      MPI_Intercomm_merge(inter_pgroup, 0 /* LOW */, &pgroup);
    } else {
      MPI_Intercomm_create(pgroup, 0, MPI_COMM_WORLD, pid_list[(gid-1)*merge_size], INTERCOMM_TAG, &inter_pgroup);
      MPI_Intercomm_merge(inter_pgroup, 1 /* HIGH */, &pgroup);
    }

    MPI_Comm_free(&inter_pgroup);
    if (pgroup_old != MPI_COMM_SELF) MPI_Comm_free(&pgroup_old);
  }

  *group_out = pgroup;
}


int main(int argc, char **argv) {
  int me, nproc, i, j, *glist;
  MPI_Comm groups[NITER];
  MPI_Group world_group;

  MPI_Init(&argc, &argv);

  MPI_Comm_rank(MPI_COMM_WORLD, &me);
  MPI_Comm_size(MPI_COMM_WORLD, &nproc);
  MPI_Comm_group(MPI_COMM_WORLD, &world_group);

  glist = (int*) malloc(nproc*sizeof(int));

  for (i = 0; i < nproc; i++)
    glist[i] = i;

  if (me == 0)
    printf("Gsize\tPGgroup (sec)\tComm (sec)\n");

  MPI_Barrier(MPI_COMM_WORLD);

  for (i = 1; i <= nproc; i*= 2) {
    double t_start, t_pg, t_comm;
    MPI_Group intracomm_group;

    /** Benchmark pgroup creation cost **/

    MPI_Barrier(MPI_COMM_WORLD);
    t_start = MPI_Wtime();

    for (j = 0; j < NITER; j++)
      pgroup_create(i, glist, &groups[j]);

    t_pg = MPI_Wtime() - t_start;

    for (j = 0; j < NITER; j++)
      pgroup_free(&groups[j]);

    /** Benchmark intracommunicator creation cost **/

    MPI_Group_incl(world_group, i, glist, &intracomm_group);
    MPI_Barrier(MPI_COMM_WORLD);
    t_start = MPI_Wtime();

    for (j = 0; j < NITER; j++)
      MPI_Comm_create(MPI_COMM_WORLD, intracomm_group, &groups[j]);

    t_comm = MPI_Wtime() - t_start;
    MPI_Group_free(&intracomm_group);

    for (j = 0; j < NITER; j++)
      pgroup_free(&groups[j]);

    if (me == 0)
      printf("%6d\t%0.9f\t%0.9f\n", i, t_pg/NITER, t_comm/NITER);

  }

  free(glist);
  MPI_Group_free(&world_group);

  MPI_Finalize();
  return 0;
}