KMR
mrmpi-crmat.c
1 /* ----------------------------------------------------------------------
2  MR-MPI = MapReduce-MPI library
3  http://www.cs.sandia.gov/~sjplimp/mapreduce.html
4  Steve Plimpton, sjplimp@sandia.gov, Sandia National Laboratories
5 
6  Copyright (2009) Sandia Corporation. Under the terms of Contract
7  DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
8  certain rights in this software. This software is distributed under
9  the modified Berkeley Software Distribution (BSD) License.
10 
11  See the README file in the top-level MapReduce directory.
12 ------------------------------------------------------------------------- */
13 
14 /* This is an example in "examples" directory in MR-MPI (11Mar13),
15  converted for KMR (2013-04-24). */
16 
17 /*
18 MapReduce random RMAT matrix generation example in C++
19 Syntax: rmat N Nz a b c d frac seed {outfile}
20  2^N = # of rows in RMAT matrix
21  Nz = non-zeroes per row
22  a,b,c,d = RMAT params (must sum to 1.0)
23  frac = RMAT randomization param (frac < 1, 0 = no randomization)
24  seed = RNG seed (positive int)
25  outfile = output RMAT matrix to this filename (optional)
26 */
27 
28 /* KMR MEMO: Run parameter from the paper (command line arguments):
29  PARAMETER=(20 4 0.57 0.19 0.19 0.05 0.1 1234) */
30 
31 #include <mpi.h>
32 #include <math.h>
33 #include <stdio.h>
34 #include <stdlib.h>
35 #include <string.h>
36 #include <assert.h>
37 #if 0
38 #include "cmapreduce.h"
39 #endif
40 #include "kmr.h"
41 
42 int generate(const struct kmr_kv_box kv,
43  const KMR_KVS *kvi, KMR_KVS *kvo, void *p,
44  const long index);
45 int cull(const struct kmr_kv_box kv[], const long n,
46  const KMR_KVS *kvi, KMR_KVS *kvo, void *p);
47 void output(char *, int, char *, int, int *, void *, void *);
48 void nonzero(char *, int, char *, int, int *, void *, void *);
49 void degree(char *, int, char *, int, int *, void *, void *);
50 void histo(char *, int, char *, int, int *, void *, void *);
51 int ncompare(char *, int, char *, int);
52 void stats(uint64_t, char *, int, char *, int, void *, void *);
53 
54 typedef struct { // RMAT params
55  int nlevels, order;
56  int nnonzero;
57  int ngenerate;
58  double a, b, c, d, fraction;
59  char *outfile;
60  FILE *fp;
61 } RMAT;
62 
63 typedef int VERTEX; // vertex ID
64 
65 typedef struct { // edge = 2 vertices
66  VERTEX vi,vj;
67 } EDGE;
68 
69 /* ---------------------------------------------------------------------- */
70 
71 int main(int argc, char **argv)
72 {
73  int me, nprocs, thlv;
74 
75  MPI_Init_thread(&argc, &argv, MPI_THREAD_SERIALIZED, &thlv);
76  MPI_Comm_rank(MPI_COMM_WORLD, &me);
77  MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
78  kmr_init();
79 
80  // parse command-line argv
81 
82  if (argc != 9 && argc != 10) {
83  if (me == 0) printf("Syntax: rmat N Nz a b c d frac seed {outfile}\n");
84  MPI_Abort(MPI_COMM_WORLD, 1);
85  }
86 
87  RMAT rmat;
88  rmat.nlevels = atoi(argv[1]);
89  rmat.nnonzero = atoi(argv[2]);
90  rmat.a = atof(argv[3]);
91  rmat.b = atof(argv[4]);
92  rmat.c = atof(argv[5]);
93  rmat.d = atof(argv[6]);
94  rmat.fraction = atof(argv[7]);
95  int seed = atoi(argv[8]);
96  if (argc == 10) {
97  int n = (int)(strlen(argv[9]) + 1);
98  rmat.outfile = malloc(sizeof(char) * (size_t)n);
99  strcpy(rmat.outfile, argv[9]);
100  } else {
101  rmat.outfile = NULL;
102  }
103 
104  if (rmat.a + rmat.b + rmat.c + rmat.d != 1.0) {
105  if (me == 0) printf("ERROR: a,b,c,d must sum to 1\n");
106  MPI_Abort(MPI_COMM_WORLD, 1);
107  }
108 
109  if (rmat.fraction >= 1.0) {
110  if (me == 0) printf("ERROR: fraction must be < 1\n");
111  MPI_Abort(MPI_COMM_WORLD, 1);
112  }
113 
114  srand48(seed+me);
115  rmat.order = (1 << rmat.nlevels);
116 
117  KMR *mr = kmr_create_context(MPI_COMM_WORLD, MPI_INFO_NULL, 0);
118  /*void *mr = MR_create(MPI_COMM_WORLD);*/
119  /*MR_set_verbosity(mr, 2);*/
120  /*MR_set_timer(mr, 1);*/
121 
122  // loop until desired number of unique nonzero entries
123 
124  struct kmr_option keepopen = {.keep_open = 1};
125 
126  MPI_Barrier(MPI_COMM_WORLD);
127  double tstart = MPI_Wtime();
128 
129  KMR_KVS *mat0 = kmr_create_kvs(mr, KMR_KV_OPAQUE, KMR_KV_OPAQUE);
130 
131  int niterate = 0;
132  int ntotal = (1 << rmat.nlevels) * rmat.nnonzero;
133  int nremain = ntotal;
134  while (nremain > 0) {
135  niterate++;
136  rmat.ngenerate = nremain/nprocs;
137  if (me < (nremain % nprocs))
138  rmat.ngenerate++;
139  kmr_map_once(mat0, &rmat, kmr_noopt, 0, generate);
140  /*MR_map_add(mr, nprocs, &generate, &rmat, 1);*/
141  /*int nunique = MR_collate(mr, NULL);*/
142  /*if (nunique == ntotal)*/
143  /*break;*/
144  /*MR_reduce(mr, &cull, &rmat);*/
145  KMR_KVS *mat1 = kmr_create_kvs(mr, KMR_KV_OPAQUE, KMR_KV_OPAQUE);
146  kmr_shuffle(mat0, mat1, kmr_noopt);
147  KMR_KVS *mat2 = kmr_create_kvs(mr, KMR_KV_OPAQUE, KMR_KV_OPAQUE);
148  kmr_reduce(mat1, mat2, 0, keepopen, cull);
149  long nunique;
150  mat2->c.stowed = 1;
151  kmr_get_element_count(mat2, &nunique);
152  mat2->c.stowed = 0;
153  nremain = (ntotal - (int)nunique);
154  mat0 = mat2;
155  }
156 
157  MPI_Barrier(MPI_COMM_WORLD);
158  double tstop = MPI_Wtime();
159 
160 #if 0
161  // output matrix if requested
162 
163  if (rmat.outfile) {
164  char fname[128];
165  sprintf(fname, "%s.%d", rmat.outfile, me);
166  rmat.fp = fopen(fname, "w");
167  if (rmat.fp == NULL) {
168  printf("ERROR: Could not open output file\n");
169  MPI_Abort(MPI_COMM_WORLD, 1);
170  }
171  void *mr2 = MR_copy(mr);
172  MR_reduce(mr2, &output, &rmat);
173  fclose(rmat.fp);
174  MR_destroy(mr2);
175  }
176 
177  // stats to screen
178  // include stats on number of nonzeroes per row
179 
180  if (me == 0) {
181  printf("%d rows in matrix\n", rmat.order);
182  printf("%d nonzeroes in matrix\n", ntotal);
183  }
184 
185  MR_reduce(mr, &nonzero, NULL);
186  MR_collate(mr, NULL);
187  MR_reduce(mr, &degree, NULL);
188  MR_collate(mr, NULL);
189  MR_reduce(mr, &histo, NULL);
190  MR_gather(mr, 1);
191  MR_sort_keys(mr, &ncompare);
192  int total = 0;
193  MR_map_mr(mr, mr, &stats, &total);
194  if (me == 0) printf("%d rows with 0 nonzeroes\n", rmat.order-total);
195 #endif
196 
197  if (me == 0) {
198  printf("%g secs to generate matrix on %d procs in %d iterations\n",
199  tstop-tstart, nprocs, niterate);
200  }
201 
202  // clean up
203 
204  kmr_free_kvs(mat0);
205  kmr_free_context(mr);
206  kmr_fin();
207  /*MR_destroy(mr);*/
208  free(rmat.outfile);
209  MPI_Finalize();
210 }
211 
212 /* ----------------------------------------------------------------------
213  generate RMAT matrix entries
214  emit one KV per edge: key = edge, value = NULL
215 ------------------------------------------------------------------------- */
216 
217 /*void generate(int itask, void *kv, void *ptr)*/
218 
219 int generate(const struct kmr_kv_box kv,
220  const KMR_KVS *kvi, KMR_KVS *kvo, void *p,
221  const long index)
222 {
223  RMAT *rmat = (RMAT *)p;
224  int nlevels = rmat->nlevels;
225  int order = rmat->order;
226  int ngenerate = rmat->ngenerate;
227  double a = rmat->a;
228  double b = rmat->b;
229  double c = rmat->c;
230  double d = rmat->d;
231  double fraction = rmat->fraction;
232 
233  int i, j, ilevel, delta, m;
234  double a1, b1, c1, d1, total, rn;
235  EDGE edge;
236 
237  for (m = 0; m < ngenerate; m++) {
238  delta = order >> 1;
239  a1 = a; b1 = b; c1 = c; d1 = d;
240  i = j = 0;
241 
242  for (ilevel = 0; ilevel < nlevels; ilevel++) {
243  rn = drand48();
244  if (rn < a1) {
245  ;
246  } else if (rn < a1+b1) {
247  j += delta;
248  } else if (rn < a1+b1+c1) {
249  i += delta;
250  } else {
251  i += delta;
252  j += delta;
253  }
254 
255  delta /= 2;
256  if (fraction > 0.0) {
257  a1 += a1*fraction * (drand48() - 0.5);
258  b1 += b1*fraction * (drand48() - 0.5);
259  c1 += c1*fraction * (drand48() - 0.5);
260  d1 += d1*fraction * (drand48() - 0.5);
261  total = a1+b1+c1+d1;
262  a1 /= total;
263  b1 /= total;
264  c1 /= total;
265  d1 /= total;
266  }
267  }
268 
269  edge.vi = i;
270  edge.vj = j;
271  struct kmr_kv_box akv = {.klen = sizeof(EDGE),
272  .vlen = 0,
273  .k.p = (char *)&edge,
274  .v.p = 0};
275  kmr_add_kv(kvo, akv);
276  /*MR_kv_add(kv, (char *)&edge, sizeof(EDGE), NULL, 0);*/
277  }
278  return 0;
279 }
280 
281 /* ----------------------------------------------------------------------
282  eliminate duplicate edges
283  input: one KMV per edge, MV has multiple entries if duplicates exist
284  output: one KV per edge: key = edge, value = NULL
285 ------------------------------------------------------------------------- */
286 
287 /*void cull(char *key, int keybytes, char *multivalue,
288  int nvalues, int *valuebytes, void *kv, void *ptr)*/
289 
290 int cull(const struct kmr_kv_box kv[], const long n,
291  const KMR_KVS *kvi, KMR_KVS *kvo, void *p)
292 {
293  struct kmr_kv_box akv = {.klen = sizeof(EDGE),
294  .vlen = 0,
295  .k.p = kv[0].k.p,
296  .v.p = 0};
297  kmr_add_kv(kvo, akv);
298  /*MR_kv_add(kv, key, keybytes, NULL, 0);*/
299  return 0;
300 }
301 
302 /* ----------------------------------------------------------------------
303  write edges to a file unique to this processor
304 ------------------------------------------------------------------------- */
305 
306 void output(char *key, int keybytes, char *multivalue,
307  int nvalues, int *valuebytes, void *kv, void *ptr)
308 {
309  RMAT *rmat = (RMAT *)ptr;
310  EDGE *edge = (EDGE *)key;
311  fprintf(rmat->fp, "%d %d 1\n", edge->vi+1, edge->vj+1);
312 }
313 
314 /* ----------------------------------------------------------------------
315  enumerate nonzeroes in each row
316  input: one KMV per edge
317  output: one KV per edge: key = row I, value = NULL
318 ------------------------------------------------------------------------- */
319 
320 void nonzero(char *key, int keybytes, char *multivalue,
321  int nvalues, int *valuebytes, void *kv, void *ptr)
322 {
323  /*EDGE *edge = (EDGE *)key;*/
324  /*MR_kv_add(kv, (char *)&edge->vi, sizeof(VERTEX), NULL,0);*/
325 }
326 
327 /* ----------------------------------------------------------------------
328  count nonzeroes in each row
329  input: one KMV per row, MV has entry for each nonzero
330  output: one KV: key = # of nonzeroes, value = NULL
331 ------------------------------------------------------------------------- */
332 
333 void degree(char *key, int keybytes, char *multivalue,
334  int nvalues, int *valuebytes, void *kv, void *ptr)
335 {
336  /*MR_kv_add(kv, (char *)&nvalues, sizeof(int), NULL, 0);*/
337 }
338 
339 /* ----------------------------------------------------------------------
340  count rows with same # of nonzeroes
341  input: one KMV per nonzero count, MV has entry for each row
342  output: one KV: key = # of nonzeroes, value = # of rows
343 ------------------------------------------------------------------------- */
344 
345 void histo(char *key, int keybytes, char *multivalue,
346  int nvalues, int *valuebytes, void *kv, void *ptr)
347 {
348  /*MR_kv_add(kv, key, keybytes, (char *)&nvalues, sizeof(int));*/
349 }
350 
351 /* ----------------------------------------------------------------------
352  compare two counts
353  order values by count, largest first
354 ------------------------------------------------------------------------- */
355 
356 int ncompare(char *p1, int len1, char *p2, int len2)
357 {
358  int i1 = *(int *) p1;
359  int i2 = *(int *) p2;
360  if (i1 > i2)
361  return -1;
362  else if (i1 < i2)
363  return 1;
364  else
365  return 0;
366 }
367 
368 /* ----------------------------------------------------------------------
369  print # of rows with a specific # of nonzeroes
370 ------------------------------------------------------------------------- */
371 
372 void stats(uint64_t itask, char *key, int keybytes, char *value,
373  int valuebytes, void *kv, void *ptr)
374 {
375  int *total = (int *) ptr;
376  int nnz = *(int *) key;
377  int ncount = *(int *) value;
378  *total += ncount;
379  printf("%d rows with %d nonzeroes\n",ncount,nnz);
380 }
Key-Value Stream (abstract).
Definition: kmr.h:587
#define kmr_reduce(KVI, KVO, ARG, OPT, R)
Reduces key-value pairs.
Definition: kmr.h:88
Options to Mapping, Shuffling, and Reduction.
Definition: kmr.h:613
int kmr_add_kv(KMR_KVS *kvs, const struct kmr_kv_box kv)
Adds a key-value pair.
Definition: kmrbase.c:751
int kmr_map_once(KMR_KVS *kvo, void *arg, struct kmr_option opt, _Bool rank_zero_only, kmr_mapfn_t m)
Maps once.
Definition: kmrbase.c:1402
#define kmr_create_kvs(MR, KF, VF)
Makes a new key-value stream (of type KMR_KVS) with the specified field datatypes.
Definition: kmr.h:71
int kmr_shuffle(KMR_KVS *kvi, KMR_KVS *kvo, struct kmr_option opt)
Shuffles key-value pairs to the appropriate destination ranks.
Definition: kmrbase.c:2036
KMR Context.
Definition: kmr.h:222
int kmr_free_kvs(KMR_KVS *kvs)
Releases a key-value stream (type KMR_KVS).
Definition: kmrbase.c:621
Handy Copy of a Key-Value Field.
Definition: kmr.h:358
int kmr_fin(void)
Clears the environment.
Definition: kmrbase.c:124
int kmr_get_element_count(KMR_KVS *kvs, long *v)
Gets the total number of key-value pairs.
Definition: kmrmoreops.c:114
#define kmr_init()
Sets up the environment.
Definition: kmr.h:747
int kmr_free_context(KMR *mr)
Releases a context created with kmr_create_context().
Definition: kmrbase.c:326
KMR Interface.
KMR * kmr_create_context(const MPI_Comm comm, const MPI_Info conf, const char *name)
Makes a new KMR context (a context has type KMR).
Definition: kmrbase.c:147