File: mapreduce.html

package info (click to toggle)
vtk 5.8.0-13
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 130,524 kB
  • sloc: cpp: 1,129,256; ansic: 708,203; tcl: 48,526; python: 20,875; xml: 6,779; yacc: 4,208; perl: 3,121; java: 2,788; lex: 931; sh: 660; asm: 471; makefile: 299
file content (2193 lines) | stat: -rw-r--r-- 100,568 bytes parent folder | download | duplicates (8)
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
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193
<HTML>
<HEAD>
<TITLE>MapReduce-MPI Library Users Manual</TITLE>
<META NAME="docnumber" CONTENT="http://www.sandia.gov/~sjplimp/mapreduce.html">
<META NAME="author" CONTENT="Sandia National Laboratories, Copyright (2009) Sandia Corporation">
<META NAME="copyright" CONTENT="This software and manual is distributed under the modified Berkeley Software Distribution (BSD) License.">
</HEAD>

<BODY>

<CENTER><A HREF = "http://www.sandia.gov/~sjplimp/mapreduce.html">MapReduce-MPI WWW Site</A> 
</CENTER>


<HR>

<H1></H1>

<H2><CENTER>MapReduce-MPI Library 
</CENTER></H2>
<P>This document describes the 13 April 2009 version of the open-source
MapReduce-MPI (MR-MPI) library that implements the <A HREF = "http://en.wikipedia.org/wiki/Mapreduce">MapReduce
operation</A> popularized by Google on top of standard MPI message
passing.  The library is designed for parallel execution on
distributed-memory platforms, but will also operate on a single
processor.  The library is written in C++, is callable from hi-level
langauges (C++, C, Fortran, Python, or other scripting languages), and
requires no additional software except linking with an MPI library if
you wish to perform MapReduces in parallel.
</P>
<P>Similar to the original Google design, a user performs a MapReduce by
writing a small program that invokes the library.  The user typically
provides two application-specific functions, a "map" and a "reduce",
that are called by the library when a MapReduce operation is executed.
"Map" and "reduce" are serial functions, meaning they are invoked
independently on individual processors on portions of your data when
performing a MapReduce operation in parallel.
</P>
<P>The goal of this library is to provide a simple and portable interface
for users to create their own MapReduce programs, which can then be
run on any desktop or large parallel machine using MPI.  See the
Background section for features and limitations of this
implementation.
</P>
<P>Source codes for the library is freely available for download from the
<A HREF = "http://www.cs.sandia.gov/~sjplimp/mapreduce.html">MR-MPI web site</A> and is licensed under the modified <A HREF = "http://en.wikipedia.org/wiki/BSD_license">Berkeley
Software Distribution (BSD) License</A>.  This basically means it can
be used by anyone for any purpose.  See the LICENSE file provided with
the distribution for more details.
</P>
<P>The distrubution includes a few examples of simple programs that
illustrate the use of MapReduce.
</P>
<P>The authors of this library are <A HREF = "http://www.cs.sandia.gov/~sjplimp">Steve Plimpton</A> at Sandia
National Laboratories and <A HREF = "http://www.cs.sandia.gov/~kddevin">Karen Devine</A> who can be contacted via
email: sjplimp at sandia.gov, kddevin at sandia.gov.
</P>
<UL><LI><A HREF = "#back">Background</A> 

<LI><A HREF = "#whatis">What is a MapReduce?</A> 

<LI><A HREF = "#start">Getting Started</A> 

<LI><A HREF = "#program">Writing a MapReduce program</A> 

<LI><A HREF = "#lib">C++ Interface to the Library</A> 

<UL><LI>  <A HREF = "#construct">Instantiate a MapReduce object</A> 

<LI>  <A HREF = "#copy">Copy a MapReduce object</A> 

<LI>  <A HREF = "#destruct">Destroy a MapReduce object</A> 

<LI>  <A HREF = "#aggregate">MapReduce::aggregate()</A> 

<LI>  <A HREF = "#clone">MapReduce::clone()</A> 

<LI>  <A HREF = "#collapse">MapReduce::collapse()</A> 

<LI>  <A HREF = "#collate">MapReduce::collate()</A> 

<LI>  <A HREF = "#compress">MapReduce::compress()</A> 

<LI>  <A HREF = "#convert">MapReduce::convert()</A> 

<LI>  <A HREF = "#gather">MapReduce::gather()</A> 

<LI>  <A HREF = "#map">MapReduce::map()</A> 

<LI>  <A HREF = "#reduce">MapReduce::reduce()</A> 

<LI>  <A HREF = "#scrunch">MapReduce::scrunch()</A> 

<LI>  <A HREF = "#sort_keys">MapReduce::sort_keys()</A> 

<LI>  <A HREF = "#sort_values">MapReduce::sort_values()</A> 

<LI>  <A HREF = "#sort_multivalues">MapReduce::sort_multivalues()</A> 

<LI>  <A HREF = "#kvstats">MapReduce::kv_stats()</A> 

<LI>  <A HREF = "#kmvstats">MapReduce::kmv_stats()</A> 

<LI>  <A HREF = "#add">KeyValue::add()</A> 

<LI>  <A HREF = "#settings">Settings and defaults</A> 
</UL>
<LI><A HREF = "#cinterface">C interface to the Library</A> 

<LI><A HREF = "#pyinterface">Python interface to the Library</A> 

<LI><A HREF = "#tech">Technical Details</A> 

<LI><A HREF = "#examples">Examples</A> 

<UL><LI>  <A HREF = "#word">Word frequency</A> 

<LI>  <A HREF = "#rmat">R-MAT matrices</A> 
</UL>

</UL>
<HR>

<A NAME = "back"></A><H3>Background 
</H3>
<P>MapReduce is the programming paradigm popularized by Google
researchers <A HREF = "#Dean">Dean and Ghemawat</A>.  Their motivation was to enable
analysis programs to be rapidly developed and deployed within Google
to operate on the massive data sets residing on their large
distributed clusters.  Their paper introduced a novel way of thinking
about certain kinds of large-scale computations as "map" operations
followed by "reduces".  The power of the paradigm is that when cast in
this way, a traditionally serial algorithm now becomes two highly
parallel application-specific operations (requiring no communication)
sandwiched around an intermediate operation that requires parallel
communication, but which can be encapsulated in a library since the
operation is independent of the application.
</P>
<P>The Google implementation of MapReduce was a C++ library with
communication between networked machines via remote procedure calls.
They allow for fault tolerance when large numbers of machines are
used, and can use disks as out-of-core memory to process huge data
sets.  Thousands of MapReduce programs have since been written by
Google researchers and are part of the daily compute tasks run by the
company.
</P>
<P>While I had heard about MapReduce, I didn't appreciate its power for
scientific computing on a monolithic distributed-memory parallel
machine, until reading a SC08 paper by <A HREF = "#Tu">Tu, et al</A> of the
D.E. Shaw company.  They showed how to think about tasks such as the
post-processing of simulation output as MapReduce operations.  In this
context it can be useful for computations that would normally be
thought of as serial, such as reading in a large data set and scanning
it for events of a desired kind.  As before, the computation can be
formulated as a highly parallel "map" followed by a "reduce".  The
encapsulated parallel operation in the middle requires all-to-all
communication to reorgnanize the data, a familiar MPI operation.
</P>
<P>Tu's implementation of MapReduce was in parallel Python with
communication between processors via MPI, again allowing disks to be
used for out-of-core operations, since their Linux cluster has one
disk per processing node.
</P>
<P>This MapReduce-MPI (MR-MPI) library is a very simple and lightweight
implementation of the basic MapReduce functionality, borrowing ideas
from both the <A HREF = "#Dean">Dean and Sanjay</A> and <A HREF = "#Tu">Tu, et al</A> papers.  It
has the following features:
</P>
<UL><LI>C++ library using MPI for inter-processor communication.  This allows
precise control over the memory allocated during a large-scale
MapReduce. 

<LI>C++ and C and Python interfaces provided.  A C++ interface means that
one or more MapReduce objects can be instantiated and invoked by the
user's program.  A C interface means that the library can also be
called from C or other hi-level languages such as Fortran.  A Python
interface means the library can be called from a Python script,
allowing you to write serial map() and reduce() functions in Python.
If your machine can run Python in parallel, you can also run a
parallel MapReduce in that manner. 

<LI>Small, portable.  The entire library is a few thousand lines of C++
code in a handful of C++ files which can be built on any machine with
a C++ compiler.  For parallel operation, you link with MPI, a standard
message passing library available on all distributed memory machines.
For serial operation, a dummy MPI library can be substituted, which is
provided.  The Python wrapper can be installed on any machine with a
version of Python that includes the ctypes module, typically Python
2.5 or later. 
</UL>
<P>This library also has the following limitations, which may be overcome
in future releases:
</P>
<UL><LI>No fault toleranace.  Current MPI implementations do not enable easy
detection of a dead processor.  So like most MPI programs, a MapReduce
operation will hang or crash if a processor goes away. 

<LI>No out-of-core operations.  Most of the large parallel machines at
Sandia do not have one disk per processor or node.  Rather they have a
few I/O nodes shared by 1000s of processors.  This makes out-of-core
processing via disk access by all processors less effective and less
portable.  While these machines do have huge aggregate memory, it does
mean the library is limited to processing data sets that will fit in
that memory.  This can be a limitation, particularly when the
intermediate data set (between the map and reduce operations) is
large. 
</UL>
<P>Finally, I call attention to <A HREF = "#Gray">recent work</A> by Alexander Gray and
colleagues at Georgia Tech.  They show that various kinds of
scientific computations such as N-body forces via multipole
expansions, k-means clustering, and machine learning algorithms, can
be formulated as MapReduce operations.  Thus there is an expanding set
of data-intense or compute-intense problems that may be amenable to
solution using a MapReduce library such as this.
</P>
<HR>

<A NAME = "whatis"></A><H3>What is a MapReduce? 
</H3>
<P>The canonical example of a MapReduce operation, described in both the
<A HREF = "#Dean">Dean and Sanjay</A> and <A HREF = "#Tu">Tu, et al</A> papers, is counting the
frequency of words in a collection of text files.  Imagine a large
corpus of text comprising Gbytes or Tbytes of data.  To count how
often each word appears, the following algorithm would work, written
in Python:
</P>
<PRE>dict = {}
for file in sys.argv[1:]:
 text = open(file,'r').read()
 words = text.split()
 for word in words:
   if word not in dict: dict[word] = 1
   else: dict[word] += 1
unique = dict.keys()
for word in unique:
 print dict[word],word 
</PRE>
<P>Dict is a "dictionary" or associative array which is a collection of
key/value pairs where the keys are unique.  In this case, the key is a
word and its value is the number of times it appears in any text file.
The program loops over files, and splits the contents into words
(separated by whitespace).  For each word, it either adds it to the
dictionary or increments its associated value.  Finally, the resulting
dictionary of unique words and their counts is printed.
</P>
<P>The drawback of this implementation is that it is inherently serial.
The files are read one by one.  More importantly the dictionary data
structure is updated one word at a time.
</P>
<P>A MapReduce formulation of the same task is as follows:
</P>
<PRE>array = []
for file in sys.argv[1:]:
 array += map(file)
newarray = collate(array)
unique = [] 
for entry in newarray:
 unique += reduce(entry)
for entry in unique:
 print entry[1],entry[0] 
</PRE>
<P>Array is now a linear list of key/value pairs where a key may appear
many times (not a dictionary).  The map() function reads a file,
splits it into words, and generates a key/value pair for each word in
the file.  The key is the word itself and the value is the integer 1.
The collate() function reorganizes the (potentially very large) list
of key/value pairs into a new array of key/value pairs where each
unique key appears exactly once and the associated value is a
concatenated list of all the values associated with the same key in
the original array.  Thus, a key/value pair in the new array would be
("dog",[1,1,1,1,1]) if the word "dog" appeared 5 times in the text
corpus.  The reduce() function takes a single key/value entry from the
new array and returns a key/value pair that has the word as its key
and the count as its value, ("dog",5) in this case.  Finally, the
elements of the unique array are printed.
</P>
<P>As written, the MapReduce algorithm could be executed on a single
processor.  However, there is now evident parallelism.  The map()
function calls are independent of each other and can be executed on
different processors simultaneously.  Ditto for the reduce() function
calls.  In this scenario, each processor would accumulate its own
local "array" and "unique" lists of key/value pairs.
</P>
<P>Also note that if the map and reduce functions are viewed as black
boxes that produce a list of key/value pairs (in the case of map) or
convert a single key/value pair into a new key/value pair (in the case
of reduce), then they are the only part of the above algorithm that is
application-specific.  The remaining portions (the collate function,
assignment of map or reduce tasks to processors, combining of the
map/reduce output across processors) can be handled behind the scenes
in an application-independent fashion.  That is the portion of the
code that is handled by the MR-MPI or other MapReduce libraries.  The
user only needs to provide a small driving program to call the library
and serial functions for performing the desired map() and reduce()
operations.
</P>
<HR>

<A NAME = "start"></A><H3>Getting Started 
</H3>
<P>Once you have
<A HREF = "http://www.sandia.gov/~sjplimp/download.html">downloaded</A> the
MapReduce MPI (MR-MPI) library, you should have the tarball
mapreduce.tar.gz on your machine.  Unpack it with the following
commands:
</P>
<PRE>gunzip mapreduce.tar.gz
tar xvf mapreduce.tar 
</PRE>
<P>which should create a mapreduce directory containing the following:
</P>
<UL><LI>README
<LI>LICENSE
<LI>doc
<LI>examples
<LI>mpistubs
<LI>python
<LI>src
<LI>user 
</UL>
<P>The doc directory contains this documentation.  The examples directory
contains a few simple MapReduce programs which call the MR-MPI
library.  These are documented by a README file in that directory and
are discussed below.  The mpistubs directory contains a dummy MPI
library which can be used to build a MapReduce program on a serial
machine.  The python contains the Python wrapper files needed to call
the MR-MPI library from Python.  The src directory contains the files
that comprise the MR-MPI library.  The user directory contains
user-contributed MapReduce programs.  See the README in that directory
for further details.
</P>
<P>To build the library for use by a C++ or C program, go to the src
directory and type
</P>
<PRE>make -f Makefile.foo 
</PRE>
<P>where you should create and use a Makefile.foo appropriate for your
machine, using one of the provided Makefiles as a template.  Note that
Makefile.serial builds the library for a serial machine, using the
dummy MPI library in mpistubs.  Other Makefiles build it for parallel
execution using a MPI library installed on your machine.  You may need
to edit one of the Makefiles to be compatible with your compilers and
MPI installation.  If you use the dummy MPI library, you will need to
build it first, by typing "make" from within mpistubs.  Again, you may
need to edit mpistubs/Makefile for your machine.
</P>
<P>If you successfully build the MR-MPI library, you should produce the
file "libmrmpi.a" which can be linked by other programs.  As discussed
below, both a C++ and C interface are part of the library, so the
library should be usable from any hi-level language.
</P>
<P>To use the MR-MPI library from Python, you don't need to build it from
the src directory.  Instead, you build it as a dynamic library from
the python directory.  Instructions are given below in the <A HREF = "#pyinterface">Python
interface</A> section.
</P>
<P>The MapReduce programs in the examples directory can be built by
typing
</P>
<PRE>make -f Makefile.foo 
</PRE>
<P>from within the examples directory.  Again, one of the provided
Makefiles may need to be modified for your platform.  Some of the
example programs are provided as a C++ program, a C program, and as a
Python script.
</P>
<HR>

<A NAME = "program"></A><H3>Writing a MapReduce program 
</H3>
<P>The usual way to use the MR-MPI library is to write a small main
program that calls the library.  In C++, your program includes two
library header files and uses the MapReduce namespace:
</P>
<PRE>#include "mapreduce.h"
#include "keyvalue.h"
using namespace MAPREDUCE_NS 
</PRE>
<P>Follow these links for info on using the library from a <A HREF = "#cinterface">C
program</A> or from a <A HREF = "#pyinterface">Python program</A>.
</P>
<P>Arguments to the library's "map" and "reduce" methods include function
pointers to serial "mymap" and "myreduce" functions in your code
(named anything you wish), which will be "called back" from the
library as it performs the parallel map and reduce operations.
</P>
<P>A typical MapReduce program involves these steps:
</P>
<PRE>MapReduce *mr = new MapReduce(MPI_COMM_WORLD);   // instantiate an MR object
mr->map(nfiles,&mymap);                          // parallel map
mr->collate()                                    // collate keys
mr->reduce(&myreduce);                           // parallel reduce
delete mr;                                       // delete the MR object 
</PRE>
<P>The main program you write may be no more complicated than this.  The
API for the MR-MPI library is a handful of methods which are
components of a MapReduce operation.  They can be combined in more
complex sequences of calls than listed above.  For example, one map()
may be followed by several reduce() operations to massage your data in
a desired way.  Output of final results is typically performed as part
of a myreduce() function you write which executes on one or more
processors and writes to a file(s) or the screen.
</P>
<P>The MR-MPI library operates on "keys" and "values" which are generated
and manipulated by your mymap() and myreduce() functions.  A key and a
value are simply byte strings of arbitrary length which are logically
associated with each other, and can thus represent anything you wish.
For example, a key can be a text string or a particle or grid cell ID.
A value can be one or more numeric values or a text string or a
composite data structure that you create.
</P>
<HR>

<A NAME = "lib"></A><H3>C++ Interface to the Library 
</H3>
<P>This section discusses how to call the MR-MPI library from a C++
program and gives a description of all its methods and variable
settings.  Use of the library from a <A HREF = "#cinterface">C program</A> (or other
hi-level language) or from <A HREF = "#pyinterface">Python</A> is discussed at the
end of the section.
</P>
<P>All the library methods operate on two basic data structures stored
within the MapReduce object, a KeyValue object (KV) and a
KeyMultiValue object (KMV).  When running in parallel, these objects
are stored in a distributed fashion across multiple processors.
</P>
<P>A KV is a collection of key/value pairs.  The same key may appear many
times in the collection, associated with values which may or may not
be the same.
</P>
<P>A KMV is also a collection of key/value pairs.  But each key in the
KMV is unique, meaning it appears exactly once (see the clone() method
for a possible exception).  The value associated with a KMV key is a
concatenated list (a multi-value) of all the values associated with
the same key in the KV.
</P>
<P>More details about how KV and KMV objects are stored are given in the
<A HREF = "#tech">Technical Details</A> section.
</P>
<P>Here is an overview of how the various library methods operate on KV
and KMV objects.  This is useful to understand, since this determines
how the various operations can be chained together in your program:
</P>
<DIV ALIGN=center><TABLE  BORDER=1 >
<TR><TD >aggregate()</TD><TD > KV -> KV</TD><TD > keys are aggregated onto procs</TD><TD > parallel</TD></TR>
<TR><TD >clone()</TD><TD > KV -> KMV</TD><TD > each KV key becomes a KMV key</TD><TD > serial</TD></TR>
<TR><TD >collapse()</TD><TD > KV -> KMV</TD><TD > all KV keys become one KMV key</TD><TD > serial</TD></TR>
<TR><TD >collate()</TD><TD > KV -> KMV</TD><TD > aggregate + convert</TD><TD > parallel</TD></TR>
<TR><TD >compress()</TD><TD > KV -> KV</TD><TD > calls back to user program to compress duplicate keys</TD><TD > serial</TD></TR>
<TR><TD >convert()</TD><TD > KV -> KMV</TD><TD > duplicate KV keys become one KMV key</TD><TD > serial</TD></TR>
<TR><TD >gather()</TD><TD > KV -> KV</TD><TD > collect keys on many procs to few procs</TD><TD > parallel</TD></TR>
<TR><TD >map()</TD><TD > create or add to a KV</TD><TD > calls back to user program to generate keys</TD><TD > serial</TD></TR>
<TR><TD >reduce()</TD><TD > KMV -> KV</TD><TD > calls back to user program to process multi-values</TD><TD > serial</TD></TR>
<TR><TD >scrunch()</TD><TD > KV -> KMV</TD><TD > gather + collapse</TD><TD > parallel</TD></TR>
<TR><TD >sort_keys()</TD><TD > KV -> KV</TD><TD > calls back to user program to sort keys</TD><TD > serial</TD></TR>
<TR><TD >sort_values()</TD><TD > KV -> KV</TD><TD > calls back to user program to sort values</TD><TD > serial</TD></TR>
<TR><TD >sort_multivalues()</TD><TD > KMV -> KMV</TD><TD > calls back to user program to sort multi-values</TD><TD > serial 
</TD></TR></TABLE></DIV>

<P>If a method creates a new KV or KMV, the old one is deleted, if it
existed.  The methods flagged as "serial" perform their operation on
the portion of a KV or KMV owned by an individual processor.  They
involve only local computation (performed simultaneously on all
processors) and no parallel comuunication.  The methods flagged as
"parallel" involve communication between processors.
</P>
<HR>

<A NAME = "construct"></A><H4>Instantiate a MapReduce object 
</H4>
<PRE>MapReduce::MapReduce(MPI_Comm comm)
MapReduce::MapReduce()
MapReduce::MapReduce(double dummy) 
</PRE>
<P>You can create a MapReduce object in any of the three ways shown, as
well as by copying from an existing MapReduce object (see
<A HREF = "#copy">below</A>).  The three creation methods differ slightly in how MPI
is initialized and finalized.
</P>
<P>In the first case, you pass an MPI communicator to the constructor.
This means your program should initialize (and finalize) MPI, which
creates the MPI_COMM_WORLD communicator (all the processors you are
running on).  Normally this is what you pass to the MapReduce
constructor, but you can pass a communicator for a subset of your
processors if desired.  You can also instantiate multiple MapReduce
objects, giving them each a communicator for all the processors or
communicators for a subset of processors.
</P>
<P>The second case can be used if your program does not use MPI at all.
The library will initialize MPI if it has not already been
initialized.  It will not finalize MPI, but this should be fine.
Worst case, your program may complain when it exits if MPI has not
been finalized.
</P>
<P>The third case is the same as the second except that the library will
finalize MPI when the last instance of a MapReduce object is
destructed.  Note that this means your program cannot delete all its
MapReduce objects in a early phase of the program and then instantiate
more MapReduce objects later.  This limitation is why the second case
is provided.  The third case is invoked by passing a double to the
constructor.  If this is done for any instantiated MapReduce object,
then the library will finalize MPI.  The value of the double doesn't
matter as it isn't used.  The use of a double is simply to make it
different than the first case, since MPI_Comm is often implemented by
MPI libraries as a type cast to an integer.
</P>
<P>As examples, any of these lines of code will create a MapReduce object:
</P>
<PRE>MapReduce *mr = new MapReduce(MPI_COMM_WORLD);
MapReduce *mr = new MapReduce();
MapReduce *mr = new MapReduce(0.0);
MapReduce mr(MPI_COMM_WORLD);
MapReduce mr();
MapReduce mr;
MapReduce mr(0.0); 
</PRE>
<HR>

<A NAME = "copy"></A><H4>Copy a MapReduce object 
</H4>
<PRE>MapReduce::MapReduce(MapReduce &) 
</PRE>
<P>Any of these segments of code will create a new MapReduce object mr2
with contents that are a deep copy from the existing MapReduce object
mr1:
</P>
<PRE>MapReduce *mr = new MapReduce(MPI_COMM_WORLD);
MapReduce *mr2 = new MapReduce(*mr); 
</PRE>
<PRE>MapReduce mr(MPI_COMM_WORLD);
MapReduce mr2 = mr; 
</PRE>
<PRE>MapReduce mr(MPI_COMM_WORLD);
MapReduce mr2(mr); 
</PRE>
<P>A deep copy means that all the key/value and/or key/multivalue pairs
contained in the first MapReduce object are copied into the new
MapReduce object.  Thus the first MapReduce object could be deleted
without affecting the new MapReduce object.
</P>
<P>This is useful if you wish to retain a copy of a set of key/value
pairs before processing it further.  See the <A HREF = "#add">KeyValue::add(KeyValue
*kv)</A> method for how to merge the key/value pairs from two
MapReduce objects into one.
</P>
<HR>

<A NAME = "destruct"></A><H4>Destroy a MapReduce object 
</H4>
<PRE>MapReduce::~MapReduce() 
</PRE>
<P>This destroys a previously created MapReduce object, freeing all the
memory it allocated internally to store keys and values.
</P>
<P>If you created the MapReduce object in this manner
</P>
<PRE>MapReduce *mr = new MapReduce(MPI_COMM_WORLD); 
</PRE>
<P>then you should destroy it with
</P>
<PRE>delete mr 
</PRE>
<HR>

<A NAME = "aggregate"></A><H4>MapReduce aggregate() method 
</H4>
<PRE>int MapReduce::aggregate(int (*myhash)(char *, int)) 
</PRE>
<P>This calls the aggregate() method of a MapReduce object, which
reorganizes a KeyValue object across processors into a new KeyValue
object.  In the original object, duplicates of the same key may be
stored on many processors.  In the new object, all duplicates of a key
are stored by the same processor.  The method returns the total number
of key/value pairs in the new KeyValue object, which will be the same
as the number in the original object.
</P>
<P>A hashing function is used to assign keys to processors.  Typically
you will not care how this is done, in which case you can specify a
NULL, i.e. mr->aggregate(NULL), and the MR-MPI library will use its
own internal hash function, which will distribute them randomly and
hopefully evenly across processors.
</P>
<P>On the other had, if you know the best way to do this for your data,
then you should provide the hashing function.  For example, if your
keys are integer IDs for particles or grid cells, you might want to
use the ID (modulo the processor count) to choose the processor it is
assigned to.  Ideally, you want a hash function that will distribute
keys to processors in a load-balanced fashion.
</P>
<P>In this example the user function is called myhash() and it must have
the following interface:
</P>
<PRE>int iproc = myhash(char *key, int keybytes) 
</PRE>
<P>Your function will be passed a key (byte string) and its length in
bytes.  Typically you want to return an integer such that 0 <= iproc <
P, where P is the number of processors.  But you can return any
integer, since the MR-MPI library uses the result in this manner to
assign the key to a processor:
</P>
<PRE>int iproc = myhash(key,keybytes) % P; 
</PRE>
<P>Because the aggregate() method will, in general, reassign all
key/value pairs to new processors, it incurs a large volume of
all-to-all communication.  However, this is performed concurrently,
taking advantage of the large bisection bandwidth most large parallel
machines provide.
</P>
<P>The aggregate() method should load-balance key/value pairs across
processors if they are initially imbalanced.
</P>
<HR>

<A NAME = "clone"></A><H4>MapReduce clone() method 
</H4>
<PRE>int MapReduce::clone() 
</PRE>
<P>This calls the clone() method of a MapReduce object, which converts a
KeyValue object directly into a KeyMultiValue object.  It simply turns
each key in KeyValue object into a key in the new KeyMultiValue
object, with the same value.  The method returns the total number of
key/value pairs in the KeyMultiValue object, which will be the same as
the number in the KeyValue object.
</P>
<P>This method essentially enables a KeyValue object to be passed
directly to a reduce operation, which requires a KeyMultiValue object
as input.  Typically you would only do this if the keys in the
KeyValue object are already unique, to avoid the extra overhead of an
aggregate() or convert() or collate(), but this is not required.  If
they are not, then there will also be duplicate keys in the
KeyMultiValue object.
</P>
<P>Note that one of the map() methods allows an existing KeyValue object
to be passed as input to a user map() function, generating a new
Keyvalue object in the process.  Thus there is typically no
need to invoke clone() followed by reduce().
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, the key/value pairs of the new KeyMultiValue
object are stored on the same processor which owns the corresponding
KeyValue pairs.
</P>
<HR>

<A NAME = "collapse"></A><H4>MapReduce collapse() method 
</H4>
<PRE>int MapReduce::collapse(char *key, int keybytes) 
</PRE>
<P>This calls the collapse() method of a MapReduce object, which
collapses a KeyValue object into a KeyMultiValue object with a single
new key, given as an argument with its length in bytes.  The single
new value in the KeyMultiValue object is a concatentated list of all
the keys and values in the KeyValue object.  The method returns the
total number of key/value pairs in the KeyMultiValue object, which
will be 1 for each processor owning pairs.
</P>
<P>For example, if the KeyValue object contains these key/value pairs:
</P>
<PRE>("dog",3), ("me",45), ("parallel",1) 
</PRE>
<P>then the new KeyMultiValue object will contain a single key/value pair:
</P>
<PRE>(key,["dog",3,"me",45,"parallel",1]) 
</PRE>
<P>This method can be used to collect a set of key/value pairs to use in
a reduce() method so that it can all be passed to a single invocation
of your myreduce() function for output.  See the <A HREF = "#tech">Technical
Details</A> section for details on how the clone() method affects
the alignment of keys and values that may eventually be passed
to your myreduce() function via the reduce() method.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor collapses the key/value pairs it
owns into a single key/value pair.  Thus each processor will assign
the same key to its new pair.  See the gather() and scrunch() methods
for ways to collect all key/value pairs on to one or a few processors.
</P>
<HR>

<A NAME = "collate"></A><H4>MapReduce collate() method 
</H4>
<PRE>int MapReduce::collate(int (*myhash)(char *, int)) 
</PRE>
<P>This calls the collate() method of a MapReduce object, which
aggregates a KeyValue object across processors and converts it into a
KeyMultiValue object.  This method is exactly the same as performing
an aggregate() followed by a convert().  The method returns the total
number of unique key/value pairs in the KeyMultiValue object.
</P>
<P>The hash argument is used by the aggregate() portion of the operation.
See the <A HREF = "#aggregate">description of that method</A> for details.
</P>
<P>Note that if your map operation does not produce duplicate keys, you
do not typically need to perform a collate().  Instead you can convert
a KeyValue object into a KeyMultiValue object directly via the clone()
method, which requires no communication.  One exception would be if
your map operation produces a KeyValue object which is highly
imbalanced across processors.  The aggregate() method should
redistribute the key/value pairs more evenly.
</P>
<P>This method is a parallel operation (aggregate), followed by an
on-processor operation (convert).
</P>
<HR>

<A NAME = "compress"></A><H4>MapReduce compress() method 
</H4>
<PRE>int MapReduce::compress(void (*mycompress)(char *, int, char *, int, int *, KeyValue *, void *), void *ptr) 
</PRE>
<P>This calls the compress() method of a MapReduce object, which
compresses a KeyValue object with duplicate keys into a new KeyValue
object, where each key appears once (on that processor) and has a
single new value.  The new value is a combination of the values
associated with that key in the original KeyValue object.  The
mycompress() function you provide generates the new value.  The method
returns the total number of key/value pairs in the new KeyValue
object.
</P>
<P>This method is used to compress a large set of key/value pairs
produced by the map() method into a smaller set before proceeding with
the rest of a MapReduce operation, e.g. with a collate() and reduce().
</P>
<P>You can give this method a pointer (void *ptr) which will be returned
to your mycompress() function.  See the <A HREF = "#tech">Technical Details</A>
section for why this can be useful.  Just specify a NULL if you don't
need this.
</P>
<P>In this example the user function is called mycompress() and it must
have the following interface, which is the same as that used
by the reduce() method:
</P>
<PRE>void mycompress(char *key, int keybytes, char *multivalue, int nvalues, int *valuebytes, KeyValue *kv, void *ptr) 
</PRE>
<P>A single key/multi-value pair is passed to your function from a
temporary KeyMultiValue object created by the library.  That object
creates a multi-value for each unique key in the KeyValue object which
contains a list of the nvalues associated with that key.  Note that
this is only the values on this processor, not across all processors.
The char *multivalue argument is a pointer to the beginning of the
multi-value which contains all nvalues, packed one after the other.
The int *valuebytes argument is an array which stores the length of
each value in bytes.  If needed, it can be used by your function to
compute an offset into char *values for where each individual value
begins.  Your function is also passed a kv pointer to a new KeyValue
object created and stored internally by the MapReduce object.
</P>
<P>Your mycompress() function should typicaly produce a single key/value
pair which it registers with the MapReduce object by calling the add()
method of the KeyValue object.  The syntax for registration is
described below with the KeyValue::add() method.  For example, if the
set of nvalues were integers, the compressed value might be the sum of
those integers.
</P>
<P>See the <A HREF = "#tech">Technical Details</A> section for details ont the
byte-alignment of keys and values that are passed to your mycompress()
function and on how to byte-align keys and values you register with
the KeyValue::add() method.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor operates only on the key/value
pairs it stores.  Thus you are NOT compressing all values associated
with a particular key across all processors, but only those currently
owned by one processor.
</P>
<HR>

<A NAME = "convert"></A><H4>MapReduce convert() method 
</H4>
<PRE>int MapReduce::convert() 
</PRE>
<P>This calls the convert() method of a MapReduce object, which converts
a KeyValue object into a KeyMultiValue object.  It does this by
finding duplicate keys (stored only by this processor) and
concatenating their values into a list of values which it associates
with the key in the KeyMultiValue object.  The method returns the
total number of key/value pairs in the KeyMultiValue object, which
will be the number of unique keys in the KeyValue object.
</P>
<P>This operation creates a hash table to find duplicate keys
efficiently.  More details are given in the <A HREF = "#tech">Technical Details</A>
section.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor converts only the key/value pairs
it owns into key/multi-value pairs.  Thus, this operation is typically
performed only after the aggregate() method has collected all
duplicate keys to the same processor.  The collate() method performs
an aggregate() followed by a convert().
</P>
<HR>

<A NAME = "gather"></A><H4>MapReduce gather() method 
</H4>
<PRE>int MapReduce::gather(int nprocs) 
</PRE>
<P>This calls the gather() method of a MapReduce object, which collects
the key/value pairs of a KeyValue object spread across all processors
to form a new KeyValue object on a subset (nprocs) of processors.
Nprocs can be 1 or any number smaller than P, the total number of
processors.  The gathering is done to the lowest ID processors, from 0
to nprocs-1.  Processors with ID >= nprocs end up with an empty
KeyValue object containing no key/value pairs.  The method returns the
total number of key/value pairs in the new KeyValue object, which will
be the same as in the original KeyValue object.
</P>
<P>This method can be used to collect the results of a reduce() to a
single processor for output.  See the collapse() and scrunch() methods
for related ways to collect key/value pairs for output.  A gather()
may also be useful before a reduce() if the number of unique key/value
pairs is small enough that you wish to perform the reduce tasks on
fewer processors.
</P>
<P>This method requires parallel point-to-point communication as
processors send their key/value pairs to other processors.
</P>
<HR>

<A NAME = "map"></A><H4>MapReduce map() method 
</H4>
<PRE>int MapReduce::map(int nmap, void (*mymap)(int, KeyValue *, void *), void *ptr)
int MapReduce::map(int nmap, void (*mymap)(int, KeyValue *, void *), void *ptr, int addflag) 
</PRE>
<PRE>int MapReduce::map(char *file, void (*mymap)(int, char *, KeyValue *, void *), void *ptr)
int MapReduce::map(char *file, void (*mymap)(int, char *, KeyValue *, void *), void *ptr, int addflag) 
</PRE>
<PRE>int MapReduce::map(int nmap, int nfiles, char **files, char sepchar, int delta, void (*mymap)(int, char *, int, KeyValue *, void *), void *ptr)
int MapReduce::map(int nmap, int nfiles, char **files, char sepchar, int delta, void (*mymap)(int, char *, int, KeyValue *, void *), void *ptr, int addflag) 
</PRE>
<PRE>int MapReduce::map(int nmap, int nfiles, char **files, char *sepstr, int delta, void (*mymap)(int, char *, int, KeyValue *, void *), void *ptr)
int MapReduce::map(int nmap, int nfiles, char **files, char *sepstr, int delta, void (*mymap)(int, char *, int, KeyValue *, void *), void *ptr, int addflag) 
</PRE>
<PRE>int MapReduce::map(KeyValue *, void (*mymap)(int, char *, int, char *, int, KeyValue *, void *), void *ptr)
int MapReduce::map(KeyValue *, void (*mymap)(int, char *, int, char *, int, KeyValue *, void *), void *ptr, int addflag) 
</PRE>
<P>This calls the map() method of a MapReduce object.  A function pointer
to a mapping function you write is specified as an argument.  This
method either creates a new KeyValue object to store all the key/value
pairs generated by your mymap function, or adds them to an existing
KeyValue object.  The method returns the total number of key/value
pairs in the KeyValue object.
</P>
<P>For the first set of variants (with and without addflag) you specify a
total number of map tasks (nmap) to perform across all processors.
The index of a map task is passed back to your mymap() function.
</P>
<P>For the second set of variants you specify a master file that contains
a list of filenames.  A filename is passed back to your mymap()
function.  The master file should list one file per line.  Blank lines
are not allowed.  Leading and trailing whitespace around the filename
is OK.
</P>
<P>For the third set of variants you specify an array of one or more file
names and a separation character (sepchar).  For the fourth set of
variants, you specify an array of one or more files names and a
separation string (sepstr).  The file(s) are split into nmap chunks
with roughly equal numbers of bytes in each chunk.  One chunk from one
file is read and passed back to your mymap() function, so your code
does not read the file.  See details below about the splitting
methodology and the delta input parameter.
</P>
<P>For the fifth set of variants, you specify an existing set of
key/value pairs, either from this MapReduce object or another one.
One key and its value are passed back to your mymap() function,
allowing you to generate new key/value pairs from an existing set.
</P>
<P>You can give any of the map() methods a pointer (void *ptr) which will
be returned to your mymap() function.  See the <A HREF = "#tech">Technical
Details</A> section for why this can be useful.  Just specify a
NULL if you don't need this.
</P>
<P>If the last argument <I>addflag</I> is omitted or is specified as 0, then
map() will create a new KeyValue object, deleting any existing
KeyValue object.  If addflag is non-zero, then key/value pairs
generated by your mymap() function are added to an existing KeyValue
object, which is created if needed.
</P>
<P>If the fifth map() variant is called using the KeyValue object stored
by the MapReduce object itself as an argument, and if addflag is 0,
then the existing KeyValue object is effectively replaced by the newly
generated key/value pairs.  If addflag is non-zero, then the newly
generated key/value pairs are added to the existing KeyValue object.
</P>
<P>In this example the user function is called mymap() and it has one of
four interfaces depending on which variant of the map() method is
invoked:
</P>
<PRE>void mymap(int itask, KeyValue *kv, void *ptr)
void mymap(int itask, char *file, KeyValue *kv, void *ptr)
void mymap(int itask, char *str, int size, KeyValue *kv, void *ptr)
void mymap(int itask, char *key, int keybytes, char *value, int valuebytes, KeyValue *kv, void *ptr) 
</PRE>
<P>In all cases, the final 2 arguments passed to your function are a
pointer to a KeyValue object (kv) stored internally by the MapReduce
object, and the original pointer you specified as an argument to the
map() method, as void *ptr.
</P>
<P>In the first case, itask is passed to your function with a value 0 <=
itask < nmap, where nmap was specified in the map() call.  For
example, you could use itask to select a file from a list stored by
your application.  Your mymap() function could open and read the file
or perform some other operation.
</P>
<P>In the second case, itask will have a value 0 <= itask < nfiles, where
nfiles is the number of filenames in the master file you specified.
Your function is also passed a single filename, which it will
presumably open and read.
</P>
<P>In the third case, itask will have a value from 0 <= itask < nmap,
where nmap was specified in the map() call and is the number of file
segments generated.  It is also passed a string of bytes (str) of
length size from one of the files.  Size includes a trailing '\0' that
is appended to the string.
</P>
<P>For map() methods that take files and a separation criterion as
arguments, you must specify nmap >= nfiles, so that there is one or
more map tasks per file.  For files that are split into multiple
chunks, the split is done at occurrences of the separation character
or string.  You specify a delta of how many extra bytes to read with
each chunk that will guarantee the splitting character or string is
found within that many bytes.  For example if the files are lines of
text, you could choose a newline character '\n' as the sepchar, and a
delta of 80 (if the longest line in your files is 80 characters).  If
the files are snapshots of simulation data where each snapshot is 1000
lines (no more than 80 characters per line), you could choose the
first line of each snapshot (e.g. "Snapshot") as the sepstr, and a
delta of 80000.  Note that if the separation character or string is
not found within delta bytes, an error will be generated.  Also note
that there is no harm in choosing a large delta so long as it is not
larger than the chunk size for a particular file.
</P>
<P>If the separation criterion is a character (sepchar), the chunk of
bytes passed to your mymap() function will start with the character
after a sepchar, and will end with a sepchar (followed by a '\0').  If
the separation criterion is a string (sepstr), the chunk of bytes
passed to your mymap() function will start with sepstr, and will end
with the character immediately preceeding a sepstr (followed by a
'\0').  Note that this means your mymap() function will be passed
different byte strings if you specify sepchar = 'A' vs sepstr = "A".
</P>
<P>In the fourth case, itask will have a value from 0 <= itask < nkey,
where nkey is the number of key/value pairs in the specified KeyValue
object.  Key and value are the byte strings for a single key/value
pair and are of length keybytes and valuebytes respectively.
</P>
<P>The MapReduce map() method assigns the map tasks to processors.
Options for how it does this can be controlled by MapReduce settings,
described below.  Basically, nmap/P tasks are assigned to each
processor, where P is the number of processors in the MPI communicator
you instantiated the MapReduce object with.
</P>
<P>Typically, your mymap() function will produce key/value pairs which it
registers with the MapReduce object by calling the add() method of the
KeyValue object.  The syntax for registration is described below with
the KeyValue::add() method.
</P>
<P>See the <A HREF = "#tech">Technical Details</A> section for details on the
length and byte-alignment of keys and values you register with the
KeyValue::add() method.
</P>
<P>Aside from the assignment of tasks to processors, this method is
really an on-processor operation, requiring no communication.  When
run in parallel, each processor generates key/value pairs and stores
them, independently of other processors.
</P>
<HR>

<A NAME = "reduce"></A><H4>MapReduce reduce() method 
</H4>
<PRE>int MapReduce::reduce(void (*myreduce)(char *, int, char *, int, int *, KeyValue *, void *), void *ptr) 
</PRE>
<P>This calls the reduce() method of a MapReduce object, passing it a
function pointer to a reduce function you write.  It operates on a
KeyMultiValue object, calling your myreduce function once for each
unique key/multi-value pair owned by that processor.  A new KeyValue
object is created which stores all the key/value pairs generated by
your myreduce() function.  The method returns the total number of new
key/value pairs stored by all processors.
</P>
<P>You can give this method a pointer (void *ptr) which will be returned
to your myreduce() function.  See the <A HREF = "#tech">Technical Details</A>
section for why this can be useful.  Just specify a NULL if you don't
need this.
</P>
<P>In this example the user function is called myreduce() and it must
have the following interface:
</P>
<PRE>void myreduce(char *key, int keybytes, char *multivalue, int nvalues, int *valuebytes, KeyValue *kv, void *ptr) 
</PRE>
<P>A single key/multi-value pair is passed to your function from the
KeyMultiValue object stored by the MapReduce object.  The key is
typically unique to this reduce task and the multi-value is a list of
the nvalues associated with that key in the KeyMultiValue object.  The
char *multivalue argument is a pointer to the beginning of the
multi-value which contains all nvalues, packed one after the other.
The int *valuebytes argument is an array which stores the length of
each value in bytes.  If needed, it can be used by your function to
compute an offset into char *values for where each individual value
begins.  Your function is also passed a kv pointer to a new KeyValue
object created and stored internally by the MapReduce object.
</P>
<P>Your myreduce() function can produce key/value pairs (though this is
not required) which it registers with the MapReduce object by calling
the add() method of the KeyValue object.  The syntax for registration
is described below with the KeyValue::add() method.  Alternatively,
your myreduce() function can write information to an output file.
</P>
<P>See the <A HREF = "#tech">Technical Details</A> section for details on the
byte-alignment of keys and values that are passed to your myreduce()
function and on how to byte-align keys and values you register with
the KeyValue::add() method.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor performs a myreduce() on each of
the key/value pairs it owns and stores any new key/value pairs it
generates.
</P>
<HR>

<A NAME = "scrunch"></A><H4>MapReduce scrunch() method 
</H4>
<PRE>int MapReduce::scrunch(int nprocs, char *key, int keybytes) 
</PRE>
<P>This calls the scrunch() method of a MapReduce object, which gathers a
KeyValue object onto nprocs and collapses it into a KeyMultiValue
object.  This method is exactly the same as performing a gather()
followed by a collapse().  The method returns the total number of
key/value pairs in the KeyMultiValue object which should be one for
each of the nprocs.
</P>
<P>The nprocs argument is used by the gather() portion of the operation.
See the <A HREF = "#gather">description of that method</A> for details.  The key and
keybytes arguments are used by the callapse() portion of the operation.
See the <A HREF = "#collapse">description of that method</A> for details.
</P>
<P>Note that if nprocs > 1, then the same key will be assigned to the
collapsed key/multi-value pairs on each processor.
</P>
<P>This method can be used to collect a set of key/value pairs to use in
a reduce() method so that it can all be passed to a single invocation
of your myreduce() function for output.
</P>
<P>This method is a parallel operation (gather), followed by an
on-processor operation (collapse).
</P>
<HR>

<A NAME = "sort_keys"></A><H4>MapReduce sort_keys() method 
</H4>
<PRE>int MapReduce::sort_keys(int (*mycompare)(char *, int, char *, int)) 
</PRE>
<P>This calls the sort_keys() method of a MapReduce object, which sorts a
KeyValue object by its keys to produce a new KeyValue object.  The
mycompare() function you provide compares pairs of keys for the sort,
since the MapReduce object does not know how to interpret the content
of your keys.  The method returns the total number of key/value pairs
in the new KeyValue object which will be the same as in the original.
</P>
<P>This method is used to sort key/value pairs by key before a KeyValue
object is transformed into a KeyMultiValue object, e.g. via the
clone(), collapse(), or convert() methods.  Note that these operations
preserve the order of paires in the KeyValue object when creating a
KeyMultiValue object, which can then be passed to your application for
output, e.g. via the reduce() method.  Note however, that sort_keys()
does NOT sort keys across all processors but only sorts the keys on
each processor within the KeyValue object.  Thus if you gather() or
aggregate() after performing a sort_keys(), the sorted order will be
lost, since those methods move key/value pairs to new processors.
</P>
<P>In this example the user function is called mycompare() and it must
have the following interface
</P>
<PRE>int mycompare(char *key1, int len1, char *key2, int len2) 
</PRE>
<P>Key1 and key2 are pointers to the byte strings for 2 keys, each of
length len1 and len2.  Your function should compare them and return a
-1, 0, or 1 if key1 is less than, equal to, or greater than key2,
respectively.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor operates only on the key/value
pairs it stores.
</P>
<HR>

<A NAME = "sort_values"></A><H4>MapReduce sort_values() method 
</H4>
<PRE>int MapReduce::sort_values(int (*mycompare)(char *, int, char *, int)) 
</PRE>
<P>This calls the sort_values() method of a MapReduce object, which sorts
a KeyValue object by its values to produce a new KeyValue object.  The
mycompare() function you provide compares pairs of values for the
sort, since the MapReduce object does not know how to interpret the
content of your values.  The method returns the total number of
key/value pairs in the new KeyValue object which will be the same as
in the original.
</P>
<P>This method is used to sort key/value pairs by value before a KeyValue
object is transformed into a KeyMultiValue object, e.g. via the
clone(), collapse(), or convert() methods.  Note that these operations
preserve the order of pairs in the KeyValue object when creating a
KeyMultiValue object, which can then be passed to your application for
output, e.g. via the reduce() method.  Note however, that
sort_values() does NOT sort values across all processors but only
sorts the values on each processor within the KeyValue object.  Thus
if you gather() or aggregate() after performing a sort_values(), the
sorted order will be lost, since those methods move key/value pairs to
new processors.
</P>
<P>In this example the user function is called mycompare() and it must
have the following interface
</P>
<PRE>int mycompare(char *value1, int len1, char *value2, int len2) 
</PRE>
<P>Value1 and value2 are pointers to the byte strings for 2 values, each
of length len1 and len2.  Your function should compare them and return
a -1, 0, or 1 if value1 is less than, equal to, or greater than
value2, respectively.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor operates only on the key/value
pairs it stores.
</P>
<HR>

<A NAME = "sort_multivalues"></A><H4>MapReduce sort_multivalues() method 
</H4>
<PRE>int MapReduce::sort_multivalues(int (*mycompare)(char *, int, char *, int)) 
</PRE>
<P>This calls the sort_multivalues() method of a MapReduce object, which
sorts the values for each key within a KeyMultiValue object to produce
a new KeyMultiValue object.  The mycompare() function you provide
compares pairs of values for the sort, since the MapReduce object does
not know how to interpret the content of your values.  The method
returns the total number of key/value pairs in the new KeyMultiValue
object which will be the same as in the original.
</P>
<P>This method can be used to sort a set of multi-values within a key
before they are passed to your application, e.g. via the reduce()
method.  Note that it typically only makes sense to use
sort_multivalues() for a KeyMultiValue object created by the convert()
or collate() methods, not KeyMultiValue objects created by the clone()
or collapse() or scrunch() methods.
</P>
<P>In this example the user function is called mycompare() and it must
have the following interface
</P>
<PRE>int mycompare(char *value1, int len1, char *value2, int len2) 
</PRE>
<P>Value1 and value2 are pointers to the byte strings for 2 values, each
of length len1 and len2.  Your function should compare them and return
a -1, 0, or 1 if value1 is less than, equal to, or greater than
value2, respectively.
</P>
<P>This method is an on-processor operation, requiring no communication.
When run in parallel, each processor operates only on the
key/multi-value pairs it stores.
</P>
<HR>

<A NAME = "kvstats"></A><H4>MapReduce kv_stats() method 
</H4>
<PRE>void MapReduce::kv_stats(int level) 
</PRE>
<P>Calling this method prints statistics about the KeyValue object stored
within the MapReduce object.  If level = 1 is specified, a one-line
summary is printed for all the key/value pairs across all processors.
If a level = 2 is specified, per-processor information is also printed
in a one-line histogram format.
</P>
<HR>

<A NAME = "kmvstats"></A><H4>MapReduce kmv_stats() method 
</H4>
<PRE>void MapReduce::kmv_stats(int level) 
</PRE>
<P>Calling this method prints statistics about the KeyMultiValue object
stored within the MapReduce object.  If level = 1 is specified, a
one-line summary is printed for all the key/multi-value pairs across
all processors.  If a level = 2 is specified, per-processor
information is also printed in a one-line histogram format.
</P>
<HR>

<A NAME = "add"></A><H4>KeyValue add() method 
</H4>
<PRE>void KeyValue::add(char *key, int keybytes, char *value, int valuebytes)
void KeyValue::add(int n, char *keys, int keybytes, char *values, int valuebytes) 
</PRE>
<PRE>void KeyValue::add(int n, char *keys, int *keybytes, char *values, int *valuebytes)
void KeyValue::add(KeyValue *kv) 
</PRE>
<P>The first three of these add methods are called by the mymap(),
mycompress(), and myreduce() functions in your program to register
key/value pairs with the KeyValue object stored by the MapReduce
object whose map(), compress(), or reduce() method was invoked.  The
first version registers a single key/value pair.  The second version
registers N key/value pairs, where the keys are all the same length
and the values are all the same length.  The third version registers a
set of N key/value pairs where the length of each key and of each
value is specified.
</P>
<P>As explained above, from the perspective of the MR-MPI library, keys
and values are variable-length byte strings.  To register such
strings, you must specify their length in bytes.  This is done via the
keybytes and valuebytes arguments, either as a single length or as a
vectors of lengths.  Note that if your key or value is a text string,
it should typically include a trailing "0" to terminate the string.
See the <A HREF = "#tech">Technical Details</A> section for details on the length
and byte-alignment of keys and values you register with the add()
method.
</P>
<P>The fourth of these add methods can be called directly to add the
key/value pairs from one MapReduce object to those of another.  For
example, these lines
</P>
<PRE>MapReduce *mr1 = new MapReduce(MPI_COMM_WORLD);
mr1->map(ntasks,&mymap,NULL);
MapReduce *mr2 = new MapReduce(*mr1);
mr2->collate();
mr2->reduce(&myreduce2,NULL);
mr1->kv->add(mr2->kv);
delete mr2;
mr1->reduce(&myreduce1,NULL); 
</PRE>
<P>would generate one set of key/value pairs from the initial map()
operation, then make a copy of them, which are then collated and
reduced to a new set of key/value pairs.  The new set of key/value
pairs are added to the original set produced by the map() operation to
form an augmented set of key/value pairs, which could be further
processed.
</P>
<HR>

<A NAME = "settings"></A><H4>Settings and defaults 
</H4>
<P>These are the library variables that can be set by your program:
</P>
<PRE>mapstyle = 0 (chunk) or 1 (stride) or 2 (master/slave)
verbosity = 0 (none) or 1 (summary) or 2 (histogrammed)
timer = 0 (none) or 1 (summary) or 2 (histogrammed) 
</PRE>
<P>Your program can set or change these values at any time, e.g.
</P>
<PRE>MapReduce *mr = new MapReduce(MPI_COMM_WORLD);
mr->verbosity = 1; 
</PRE>
<P>The <I>mapstyle</I> setting determines how the N map tasks are assigned to
the P processors.
</P>
<P>A value of 0 means split the tasks into "chunks" so that processor 0
is given tasks from 0 to N/P, proc 1 is given tasks from N/P to 2N/P,
etc.  Proc P-1 is given tasks from N - N/P to N.
</P>
<P>A value of 1 means "strided" assignment, so proc 0 is given tasks
0,P,2P,etc and proc 1 is given tasks 1,P+1,2P+1,etc and so forth.
</P>
<P>A value of 2 uses a "master/slave" paradigm for assigning tasks.  Proc
0 becomes the "master"; the remaining processors are "slaves".  Each
is given an initial task by the master and reports back when it is
finished.  It is then assigned the next available task which continues
until all tasks are completed.  This is a good choice if the CPU time
required by various mapping tasks varies greatly, since it will tend
to load-balance the work across processors.  Note however that proc 0
performs no mapping tasks.
</P>
<P>The default value for <I>mapstyle</I> is 0.
</P>
<P>The <I>verbosity</I> setting determines how much diagnostic output each
library call prints to the screen.  A value of 0 means "none".  A
value of 1 means a "summary" of the results across all processors is
printed, typically a count of total key/value pairs and the memory
required to store them.  A value of 2 prints the summary results and
also a "histogram" of these quantities by processor, so that you can
detect memory usage imbalance.
</P>
<P>The default value for <I>verbosity</I> is 0.
</P>
<P>The <I>timer</I> setting prints out timing information for each call to the
library.  A value of 0 means "none".  A value of 1 invokes an
MPI_Barrier() at the beginning and end of the operation and prints the
elapsed time, which will be the same on all processors.  A value of 2
invokes no MPI_Barrier() calls and prints a one-line summary of timing
results across all processors and also a "histogram" of the time on
each processor, so that you can detect computational imbalance.
</P>
<P>The default value for <I>timer</I> is 0.
</P>
<HR>

<A NAME = "cinterface"></A><H4>C interface to the Library 
</H4>
<P>The MR-MPI library can be called from a C program, using the interface
defined in src/cmapreduce.h.  This is a C file which should be
included in your C program to define the API to the library:
</P>
<PRE>#include "cmapreduce.h" 
</PRE>
<P>Note that the C interface should also be usable to call the MapReduce
MPI library from Fortran or other hi-level languages, including
scripting languages.  See information below on how to do this from
<A HREF = "#pyinterface">Python</A>.
</P>
<P>The C interface consists of the following functions.  Their
functionality and arguments are described in the <A HREF = "#lib">C++ interface
section</A>.
</P>
<PRE>void *MR_create(MPI_Comm comm);
void *MR_create_mpi();
void *MR_create_mpi_finalize();
void *MR_copy(void *MRptr);
void MR_destroy(void *MRptr); 
</PRE>
<PRE>int MR_aggregate(void *MRptr, int (*myhash)(char *, int));
int MR_clone(void *MRptr);
int MR_collapse(void *MRptr, char *key, int keybytes);
int MR_collate(void *MRptr, int (*myhash)(char *, int));
int MR_compress(void *MRptr, 
    void (*mycompress)(char *, int, char *,
           int, int *, void *KVptr, void *APPptr),
    void *APPptr);
int MR_convert(void *MRptr);
int MR_gather(void *MRptr, int numprocs); 
</PRE>
<PRE>int MR_map(void *MRptr, int nmap,
     void (*mymap)(int, void *KVptr, void *APPptr),
     void *APPptr);
int MR_map_add(void *MRptr, int nmap,
         void (*mymap)(int, void *KVptr, void *APPptr),
         void *APPptr, int addflag);
int MR_map_file_list(void *MRptr, char *file,
               void (*mymap)(int, char *, void *KVptr, void *APPptr),
         void *APPptr);
int MR_map_file_list_add(void *MRptr, char *file,
       void (*mymap)(int, char *, void *KVptr, void *APPptr),
       void *APPptr, int addflag);
int MR_map_file_char(void *MRptr, int nmap, int files, char **files,
         char sepchar, int delta,
         void (*mymap)(int, char *, int, void *KVptr, void *APPptr),
         void *APPptr);
int MR_map_file_char_add(void *MRptr, int nmap, int files, char **files,
       char sepchar, int delta,
       void (*mymap)(int, char *, int, void *KVptr, void *APPptr),
       void *APPptr, int addflag);
int MR_map_file_str(void *MRptr, int nmap, int files, char **files,
        char *sepstr, int delta,
        void (*mymap)(int, char *, int, void *KVptr, void *APPptr),
        void *APPptr);
int MR_map_file_str_add(void *MRptr, int nmap, int files, char **files,
      char *sepstr, int delta,
      void (*mymap)(int, char *, int, void *KVptr, void *APPptr),
      void *APPptr, int addflag);
int MR_map_kv(void *MRptr, void *MRptr2,
        void (*mymap)(int, char *, int, char *, int *, void *KVptr, void *APPptr),
        void *APPptr);
int MR_map_kv_add(void *MRptr, void *MRptr2,
      void (*mymap)(int, char *, int, char *, int *, void *KVptr, void *APPptr),
      void *APPptr, int addflag); 
</PRE>
<PRE>int MR_reduce(void *MRptr,
        void (*myreduce)(char *, int, char *,
             int, int *, void *KVptr, void *APPptr),
        void *APPptr);
int MR_scrunch(void *MRptr, int numprocs, char *key, int keybytes); 
</PRE>
<PRE>int MR_sort_keys(void *MRptr, 
     int (*mycompare)(char *, int, char *, int));
int MR_sort_values(void *MRptr,
       int (*mycompare)(char *, int, char *, int));
int MR_sort_multivalues(void *MRptr,
      int (*mycompare)(char *, int, char *, int)); 
</PRE>
<PRE>void MR_kv_stats(void *MRptr, int level);
void MR_kmv_stats(void *MRptr, int level); 
</PRE>
<PRE>void MR_set_mapstyle(void *MRptr, int value);
void MR_set_verbosity(void *MRptr, int value);
void MR_set_timer(void *MRptr, int value); 
</PRE>
<PRE>void MR_kv_add(void *KVptr, char *key, int keybytes, 
         char *value, int valuebytes);
void MR_kv_add_multi_static(void *KVptr, int n,
          char *key, int keybytes,
          char *value, int valuebytes);
void MR_kv_add_multi_dynamic(void *KVptr, int n,
          char *key, int *keybytes,
          char *value, int *valuebytes);
void MR_kv_add_kv(void *MRptr, void *MRptr2); 
</PRE>
<P>These functions correspond one-to-one with the C++ methods described
above, except that for C++ methods with multiple interfaces
(e.g. create() or map()), there are multiple C functions, with
slightly different names.  The MR_set() functions are added to the C
interface to enable the corresponding library variables to be set.
</P>
<P>Note that when you call MR_create() or MR_copy(), they return a "void
*MRptr" which is a pointer to the MapReduce object created by the
library.  This pointer is used as the first argument of all the other
MR calls.  This means a C program can effectively instantiate multiple
MapReduce objects by simply keeping track of the pointers returned to
it.
</P>
<P>The remaining arguments of each function call are the same as those
used with the C++ methods.  The only exceptions are several of the
MR_kv_add() functions which take a KVptr as their first argument.
This is a pointer to a KeyValue object.  These calls are made from
your program's mymap(), myreduce(), and mycompress() functions to
register key/value pairs with the MR-MPI library.  The KVptr is passed
as an argument to your functions when they are called back from the
MR-MPI library.
</P>
<P>See the C programs in the examples directory for examples of how these
calls are made from a C program.  They are conceptually identical to
the C++ programs in the same directory.
</P>
<HR>

<A NAME = "pyinterface"></A><H4>Python interface to the Library 
</H4>
<P>A Python wrapper for the MR-MPI library is included in the
distribution.  The advantage of using Python is how concise the
language is, enabling rapid development and debugging of MapReduce
programs.  The disadvantage is speed, since Python is slower than a
compiled language.  Using the MR-MPI library from Python incurs two
additional overheads, discussed in the <A HREF = "#tech">Technical Details</A>
section.
</P>
<P>Before using the MR-MPI library in a Python script, the Python on your
machine must be "extended" to include an interface to the MR-MPI
library.  If your Python script will invoke MPI operations, you will
also need to extend your Python with an interface to MPI itself.
</P>
<P>Thus you should first decide how you intend to use the MR-MPI library
from Python.  There are 3 options:
</P>
<P>(1) Use the library on a single processor running Python.
</P>
<P>(2) Use the library in parallel, where each processor runs Python, but
your application script does not use MPI.
</P>
<P>(3) Use the library in parallel, where each processor runs Python, and
your application also makes MPI calls through a Python/MPI interface.
</P>
<P>Note that for (2) and (3) you will not be able to interact with Python
interactively by typing commands and getting a response.  This is
because when you have multiple instances of Python running (e.g. on a
parallel machine) they cannot all read what you type.
</P>
<P>Working in mode (1) does not require your machine to have MPI
installed.  You should extend your Python with a serial version of the
MR-MPI library and its dummy MPI library.  See instructions below on
how to do this.
</P>
<P>Working in mode (2) requires your machine to have an MPI library
installed, but your Python does not need to be extended with MPI
itself.  The MPI library must be a shared library (e.g. a *.so file on
Linux) which is not typically created when MPI is built/installed.
See instruction below on how to do this.  You should extend your
Python with the parallel MR-MPI library which will use the shared MPI
system library.  See instructions below on how to do this.
</P>
<P>Working in mode (3) requires your machine to have MPI installed (as a
shared library as in (2)).  You must also extend your Python with the
parallel MR-MPI library (same as in (2)) and with MPI itself, via one
of several available Python/MPI packages.  See instructions below on
how to do the latter task.
</P>
<P>The following sub-sections cover the rest of the Python discussion:
</P>
<UL><LI>Extending Python with a serial version of the MR-MPI library
<LI>Creating a shared MPI library
<LI>Extending Python with a parallel version of the MR-MPI library
<LI>Extending Python with MPI itself
<LI>Testing the MR-MPI library from Python
<LI>Using the MR-MPI library from Python 
</UL>
<HR>

<P><B>Extending Python with a serial version of the MR-MPI library</B>
</P>
<P>From the python directory, type
</P>
<PRE>python setup_serial.py build 
</PRE>
<P>and then one of these commands:
</P>
<PRE>sudo python setup_serial.py install
python setup_serial.py install --home=~/foo 
</PRE>
<P>The "build" command should compile all the needed MR-MPI C++ files,
including the dummy MPI library.  The first "install" command will put
the needed files in your Python's site-packages sub-directory, so that
Python can load them.  For example, if you installed Python yourself
on a Linux machine, it would typically be somewhere like
/usr/local/lib/python2.5/site-packages.  Installing Python packages
this way often requires you to be able to write to the Python
directories, which may require root priveleges, hence the "sudo"
prefix.  If this is not the case, you can drop the "sudo".
</P>
<P>Alternatively, you can install the MR-MPI files (or any other Python
packages) in your own user space.  The second "install" command does
this, where you should replace "foo" with your directory of choice.
</P>
<P>If these commands are successful, an <I>mrmpi.py</I> and <I>_mrmpi_serial.so</I>
file will be put in the appropriate directory.
</P>
<HR>

<P><B>Creating a shared MPI library</B>
</P>
<P>A shared library is one that is dynamically loadable, which is what
Python requires.  On Linux this is a library file that ends in ".so",
not ".a".  Such a shared library is normally not built if you
installed MPI yourself, but it is easy to do.  Here is how to do it
for <A HREF = "http://www-unix.mcs.anl.gov/mpi">MPICH</A>, a popular open-source version of MPI, distributed
by Argonne National Labs.  From within the mpich directory, type
</P>
<PRE>./configure --enable-sharedlib=gcc
make
make install 
</PRE>
<P>You may need to use "sudo make install" in place of the last line.
The end result should be the file libmpich.so in /usr/local/lib.  Note
that if the file libmpich.a already existed in /usr/local/lib, you
will now have both a static and shared MPICH library.  This will be
fine for Python MR-MPI since it only uses the shared library.  But if
you build other codes with libmpich.a, then those builds may fail if
the linker uses libmpich.so instead, unless other dynamic libraries
are also linked to.
</P>
<HR>

<P><B>Extending Python with a parallel version of the MR-MPI library</B>
</P>
<P>From the python directory, type
</P>
<PRE>python setup.py build 
</PRE>
<P>and then one of these commands:
</P>
<PRE>sudo python setup.py install
python setup.py install --home=~/foo 
</PRE>
<P>The "build" command should compile all the needed MR-MPI C++ files,
which will require MPI to be installed on your system.  This means it
must find both the header file mpi.h and a shared library file,
e.g. libmpich.so if the MPICH version of MPI is installed.  See the
preceding section for how to create a build MPI as a shared library if
it does not exist.
</P>
<P>The first "install" command will put the needed files in your Python's
site-packages sub-directory, so that Python can load them.  For
example, if you installed Python yourself on a Linux machine, it would
typically be somewhere like /usr/local/lib/python2.5/site-packages.
Installing Python packages this way often requires you to be able to
write to the Python directories, which may require root priveleges,
hence the "sudo" prefix.  If this is not the case, you can drop the
"sudo".
</P>
<P>Alternatively, you can install the MR-MPI files (or any other Python
packages) in your own user space.  The second "install" command does
this, where you should replace "foo" with your directory of choice.
</P>
<P>If these commands are successful, an <I>mrmpi.py</I> and <I>_mrmpi.so</I> file
will be put in the appropriate directory.
</P>
<HR>

<P><B>Extending Python with MPI itself</B>
</P>
<P>There are several Python packages available that purport to wrap MPI
and allow its functions to be called from Python.
</P>
<P>These include
</P>
<UL><LI><A HREF = "http://pympi.sourceforge.net/">pyMPI</A>
<LI><A HREF = "http://code.google.com/p/maroonmpi/">maroonmpi</A>
<LI><A HREF = "http://code.google.com/p/mpi4py/">mpi4py</A>
<LI><A HREF = "http://nbcr.sdsc.edu/forum/viewtopic.php?t=89&sid=c997fefc3933bd66204875b436940f16">myMPI</A>
<LI><A HREF = "http://datamining.anu.edu.au/~ole/pypar/">Pypar</A> 
</UL>
<P>All of these except pyMPI work by wrapping the MPI library (which must
be available on your system as a shared library, as discussed above),
and exposing (some portion of) its interface to your Python script.
This means they cannot be used interactively in parallel, since they
do not address the issue of interactive input to multiple instances of
Python running on different processors.  The one exception is pyMPI,
which alters the Python interpreter to address this issue, and (I
believe) creates a new alternate executable (in place of python
itself) as a result.
</P>
<P>In principle any of these Python/MPI packages should work with the
MR-MPI library.  However, when I downloaded and looked at a few of
them, their docuemtation was incomplete and I had trouble with their
installation.  It's not clear if some of the packages are still being
actively developed and supported.
</P>
<P>The one I recommend, since I have successfully used it with the MR-MPI
library, is Pypar.  Pypar requires the ubiqtuitous <A HREF = "http://numpy.scipy.org">Numpy
package</A> be installed in your Python.  After
launching python, type
</P>
<PRE>>>> import numpy 
</PRE>
<P>to see if it is installed.  If not, here is how to install it (version
1.3.0b1 as of April 2009).  Unpack the numpy tarball and from its
top-level directory, type
</P>
<PRE>python setup.py build
sudo python setup.py install 
</PRE>
<P>The "sudo" is only needed if required to copy Numpy files into your
Python distribution's site-packages directory.
</P>
<P>To install PyPar (version pypar-2.1.0_66 as of April 2009), unpack it
and from its "source" directory, type
</P>
<PRE>python setup.py build
sudo python setup.py install 
</PRE>
<P>Again, the "sudo" is only needed if required to copy PyPar files into
your Python distribution's site-packages directory.
</P>
<P>If you have successully installed Pypar, you should be able to run
python serially and type
</P>
<PRE>>>> import pypar 
</PRE>
<P>without error.  You should also be able to run python in parallel
on a simple test script
</P>
<PRE>% mpirun -np 4 python test.script 
</PRE>
<P>where test.script contains the lines
</P>
<PRE>import pypar
print "Proc %d out of %d procs" % (pypar.rank(),pypar.size()) 
</PRE>
<P>and see one line of output for each processor you ran on.
</P>
<HR>

<P><B>Testing the MR-MPI library from Python</B>
</P>
<P>Before importing the MR-MPI library in a Python program, one more step
is needed.  The interface to the library is via Python ctypes, which
loads the shared MR-MPI library via a CDLL() call, which in turn is a
wrapper on the C-library dlopen().  This command is different than a
normal Python "import" and needs to be able to find the MR-MPI shared
library, which is either in the Python site-packages directory or in a
local directory you specified in the "python setup.py install"
command, as described above.
</P>
<P>The simplest way to do this is add a line like this to your
.cshrc or other shell start-up file.
</P>
<PRE>setenv LD_LIBRARY_PATH $<I>LD_LIBRARY_PATH</I>:/usr/local/lib/python2.5/site-packages 
</PRE>
<P>and then execute the file to insure the path has been updated.  This
will extend the path that dlopen() uses to look for shared libraries.
</P>
<P>To test if the MR-MPI library has been successfully installed,
launch python in serial and type
</P>
<PRE>>>> from mrmpi import mrmpi
>>> mr = mrmpi() 
</PRE>
<P>If you get no errors, you're ready to use the library, as described
below.
</P>
<P>If you built the MR-MPI library for parallel use, launch python 
in parallel
</P>
<PRE>% mpirun -np 4 python test.script 
</PRE>
<P>where test.script contains the lines
</P>
<PRE>import pypar
from mrmpi import mrmpi
mr = mrmpi()
print "Proc %d out of %d procs has" % (pypar.rank(),pypar.size()), mr 
</PRE>
<P>Again, if you get no errors, you're good to go.
</P>
<HR>

<P><B>Using the MR-MPI library from Python</B>
</P>
<P>The Python interface to the MR-MPI library consists of an "mrmpi"
class which creates a "mrmpi" object, with a set of methods that can
be invoked on that object.  The sample code lines below assume you
have first imported the "mrmpi" module as follows:
</P>
<PRE>from mrmpi import mrmpi 
</PRE>
<P>Note that when your script imports the Pypar package (same with some
other Python/MPI packages), it initializes MPI for you.  Pypar does
not, however, make the global MPI communicator (MPI_COMM_WORLD)
visible to your program, so you can't pass it to the MR-MPI library.
When using Pypar, the last line of your input script should thus be
pypar.finalize(), to insure MPI is shut down correctly.
</P>
<P>Some of the methods defined by the mrmpi class take callback functions
as arguments, e.g. map() and reduce().  These are Python functions you
define elsewhere in your script.  When you register "keys" and
"values" with the library, they can be simple quantities like strings
or ints or floats.  Or they can be Python data structures like lists
or tuples.
</P>
<P>These are the class methods defined by the mrmpi module.  Their
functionality and arguments are described in the <A HREF = "#lib">C++ interface
section</A>.
</P>
<PRE>mr = mrmpi()                # create an mrmpi object
mr = mrmpi(mpi_comm)        # ditto, but with a specified MPI communicator
mr = mrmpi(0.0)             # ditto, and the library will finalize MPI 
</PRE>
<PRE>mr2 = mr.copy()             # copy mr to create mr2 
</PRE>
<PRE>mr.destroy()                # destroy an mrmpi object, freeing its memory
                            # this will also occur if Python garbage collects 
</PRE>
<PRE>mr.aggregate()
mr.aggregate(myhash)        # if specified, myhash is a hash function
          #   called back from the library as myhash(key)
          # myhash() should return an integer (a proc ID)
mr.clone()
mr.collapse(key)
mr.collate()
mr.collate(myhash)          # if specified, myhash is the same function
          #   as for aggregate() 
</PRE>
<PRE>mr.compress(mycompress)     # mycompress is a function called back from the
          #   library as mycompress(key,mvalue,mr,ptr)
          #   where mvalue is a list of values associated
          #   with the key, mr is the MapReduce object,
          #   and you (optionally) provide ptr (see below)
          # your mycompress function should typically
          #   make calls like mr->add(key,value)
mr.compress(mycompress,ptr) # if specified, ptr is any Python datum
                #    and is passed back to your mycompress()
          # if not specified, ptr = None 
</PRE>
<PRE>mr.convert()
mr.gather(nprocs) 
</PRE>
<PRE>mr.map(nmap,mymap)          # mymap is a function called back from the
          #   library as mymap(itask,mr,ptr)
          #   where mr is the MapReduce object,
          #   and you (optionally) provide ptr (see below)
          # your mymap function should typically
          #   make calls like mr->add(key,value)
mr.map(nmap,mymap,ptr)      # if specified, ptr is any Python datum
          #    and is passed back to your mymap()
          # if not specified, ptr = None
mr.map(nmap,mymap,ptr,addflag) # if addflag is specfied as a non-zero int,
             #   new key/value pairs will be added to the
             #   existing key/value pairs 
</PRE>
<PRE>mr.map_file_list(file,mymap) # file is a file containing a list of filenames
           # mymap is a function called back from the
           #   library as mymap(itask,filename,mr,ptr)
           # as above, ptr and addflag are optional args
mr.map_file_char(nmap,files,sepchar,delta,mymap)
                             # files is a list of filenames
           # mymap is a function called back from the
           #   library as mymap(itask,str,mr,ptr)
           # as above, ptr and addflag are optional args
mr.map_file_str(nmap,files,sepstr,delta,mymap)
                             # files is a list of filenames
           # mymap is a function called back from the
           #   library as mymap(itask,str,mr,ptr)
           # as above, ptr and addflag are optional args
mr.map_kv(mr2,mymap)         # pass key/values in mr2 to mymap
                             # mymap is a function called back from the
           #   library as mymap(itask,key,value,mr,ptr)
           # as above, ptr and addflag are optional args 
</PRE>
<PRE>mr.reduce(myreduce)         # myreduce is a function called back from the
          #   library as myreduce(key,mvalue,mr,ptr)
          #   where mvalue is a list of values associated
          #   with the key, mr is the MapReduce object,
          #   and you (optionally) provide ptr (see below)
          # your myreduce function should typically
          #   make calls like mr->add(key,value)
mr.reduce(myreduce,ptr)     # if specified, ptr is any Python datum
          #    and is passed back to your myreduce()
          # if not specified, ptr = None
mr.scrunch(nprocs,key) 
</PRE>
<PRE>mr.sort_keys(mycompare)
mr.sort_values(mycompare)
mr.sort_multivalues(mycompare) # compare is a function called back from the
             #   library as mycompare(a,b) where
             #   a and b are two keys or two values
             # your mycompare() should compare them
             #   and return a -1, 0, or 1 
             #   if a < b, or a == b, or a > b 
</PRE>
<PRE>mr.kv_stats(level)
mr.kmv_stats(level) 
</PRE>
<PRE>mr.mapstyle(value)             # set the mapstyle to value
mr.verbosity(value)            # set the verbosity to value
mr.timer(value)                # set the timer to value 
</PRE>
<PRE>mr.add(key,value)                 # add single key and value
mr.add_multi_static(keys,values)  # add list of keys and values
          # all keys are assumed to be same length
          # all values are assumed to be same length
mr.add_multi_dynamic(keys,values) # add list of keys and values
          # each key may be different length
          # each value may be different length
mr.add_kv(mr2)                    # add key/values in mr2 to those in mr 
</PRE>
<P>These class methods correspond one-to-one with the C++ methods
described above, except that for C++ methods with multiple interfaces
(e.g. create() or map()), there are multiple Python methods with
slightly different names, similar to the <A HREF = "#cinterface">C interface</A>.
</P>
<P>See the Python scripts in the examples directory for examples of how
these calls are made from a Python program.  They are conceptually
identical to the C++ and C programs in the same directory.
</P>
<HR>

<A NAME = "tech"></A><H3>Technical Details 
</H3>
<P>This section provides additional details about using the MapReduce
library and how it is implemented.  These topics are covered:
</P>
<UL><LI>Length and byte-alignment of keys and values
<LI>Memory requirements for KeyValue and KeyMultiValue objects
<LI>Hash functions
<LI>Callback functions
<LI>Python overhead
<LI>Error messages 
</UL>
<H4>Length and byte-alignment of keys and values 
</H4>
<P>As explained in <A HREF = "#lib">this section</A>, keys and values are
variable-length strings of bytes.  The MR-MPI library knows nothing of
their contents and simply treats them as contiguous chunks of bytes.
</P>
<P>When you register a key and value in your <A HREF = "#map">mymap()</A> or
<A HREF = "#compress">mycompress()</A> or <A HREF = "#reduce">myreduce()</A> function via the
<A HREF = "#add">KeyValue::add() method</A>, you specify their lengths in bytes.
Keys and values are typically returned to your program for further
processing or output, e.g. as arguments passed to your myreduce()
function by the <A HREF = "#reduce">MapReduce reduce() operation</A>, and their
lengths are also returned to your function.
</P>
<P>Keys and values are passed as character pointers to your functions
which may need to convert the pointer to an appropriate data type and
then correctly interpret the byte string.  If the key or value is a
variable-length text string, it should typically be terminated by a
"0", so that C-library style string functions can be invoked on it.
If a key or value is a complex data structure, your function must be
able to decode it.
</P>
<P>A related issue with keys and values is the byte-alignment of integer
or floating point values they include.  For example, it is usually a
bad idea to store an 8-byte double in memory such that it is
mis-aligned with respect to 8-byte boundaries.  The reason is that
accessing a mis-aligned double for computaion may be slower.
</P>
<P>As an example, say your "value" is a 4-byte integer followed by an
8-byte double.  You might think it can be stored and registered as 12
contiguous bytes.  However, this would likely mean the double is
mis-aligned.  One solution is to convert the int to a double before
storing both quantities in a 16-byte value string.  Another solution
is to create a struct to store the int and double and use the sizeof()
function to determine the length of the struct and use that as the
length of your "value".  The compiler should then guarantee proper
alignment of each structure member.
</P>
<P>Special care should also be taken if your keys and values are
heterogeneous.  This is because the MR-MPI library packs keys one
after the other into one long byte string, and similarly for values;
see below for a discussion of memory usage.  When values from a
KeyValue object are orgnized into a multi-value in a KeyMultiValue
object, the values are likewise packed one after the other.  Thus even
if the components of your keys and values were aligned when you
registered them with the <A HREF = "#add">KeyValue::add() method</A>, they may not
be if they later get mixed in with other keys and values of different
lengths.  For example, the <A HREF = "#collapse">collapse</A> operation creates a
multi-value that is sequence of key,value,key,value,etc from a KV.  If
the keys are variable-length text strings and the values are ints,
then the values will not be aligned on 4-byte boundaries.
</P>
<P>The solution in this case is for your callback function to copy the
bytes of a key or value into a local data structure (e.g. using the C
memcpy() function) that has the proper alignment.  E.g. in the collapse
example above, these lines of code:
</P>
<PRE>int myvalue;
memcpy(&myvalue,&values<B>offsets<B>i</B></B>,sizeof(int)); 
</PRE>
<P>would load the 4 bytes of the Ith value in the multi-value into the
local integer "myvalue", where it can safely be used for computation.
</P>
<H4>Memory requirements for KeyValue and KeyMultiValue objects 
</H4>
<P>KeyValue and KeyMultiValue objects were described in <A HREF = "#lib">this
section</A>.  An instance of a MapReduce object contains a single
KeyValue object (KV) and a single KeyMultiValue object (KMV),
depending on which methods you have invoked.
</P>
<P>The memory cost for storing a KV is as follows.  The key and value
strings are packed into contiguous arrays of bytes.  For each key and
for each value an integer is also stored, which is the offset into the
byte arrays of where that key or value starts.  Thus the total memory
of a KV is the memory for the key/value data itself plus 2 ints per
pair.
</P>
<P>Recall that a KMV contains key/multi-value pairs where the number of
pairs is typically the number of unique keys in the original KV.  The
memory cost for storing a KMV is as follows.  The key and multi-value
strings are packed into contiguous arrays of bytes.  For each key and
for each multi-value an integer is also stored, which is the offset
into the byte arrays of where that key or multi-value starts.
For each key an additional integer stores the number of values in the
associated multi-value.  An array of offsets for individual values in 
each multi-value is also stored.
</P>
<P>Thus the total memory of a KMV is the memory for the key/multi-value
data itself plus 3 ints per pair plus 1 int per value in the original
KV.  Note that memory for key data is typically less than for the
original KV, since the KMV only stores unique keys.  The memory for
multi-value data is exactly the same as the value data in the original
KV, since all the KV values are in the multi-values.
</P>
<P>When the KMV is formed from a KV via a <A HREF = "#convert">convert</A> operation,
additional memory is temporarily required for a hash table used to
find duplicate keys.  This requires approximately 5 ints per unique
key in the KMV.
</P>
<P>Note that in parallel, for a KV or KMV, each processor stores the
above data for only the fraction of key/value pairs it generated
during a map operation or acquired during other operations.  If this
is imbalanced, one processor may run out of memory before others do.
</P>
<P>If you run out of memory when performing a MapReduce operation, the
library should print an error message and die.  One solution is to run
on more processors.  Another way to save memory, if possible, is to
write your mymap() function so that it combines multiple key/value
pairs into fewer or smaller key/value pairs before registering them
with the MapReduce object.  Calling the compress() method before
performing an aggregate() or collate() may also save memory.
</P>
<H4>Hash functions 
</H4>
<P>The convert() and collate() methods use a hash function to organize
keys and find duplicates.  The MR-MPI library uses the hashlittle()
function from lookup3.c, written by Bob Jenkins and available freely
on the WWW.  It operates on arbitrary-length byte strings (a key) and
produces a 32-bit integer hash value, a portion of which is used as a
bucket index into a hash table.
</P>
<H4>Callback functions 
</H4>
<P>Several of the library methods take a callback function as an
argument, meaning that function is called back to from the library
when the method is invoked.  These functions are part of your
MapReduce program and can perform any operation you wish on your data
(or on no data), so long as they produce the appropriate information.
E.g. they generate key/value pairs in the case of map() or compress()
or reduce(), or they hash a key to a processor in the case of
aggregate() or collate(), or they compare two keys or values in the
case of sort_keys() or sort_values().
</P>
<P>The mymap() and myreduce() functions can perform simple operations or
very complex, compute-intensive operations.  For example, if your
parallel machine supports it, they could invoke another program or
script to read/parse an input file or calculate some result.
</P>
<P>Note that in your program, a callback function CANNOT be a class
method unless it is declared to be "static".  It can also be a
non-class method, i.e. just a stand-alone function.  In either case,
such a function cannot access class data.
</P>
<P>One way to get around this restriction is to define global variables
that allow your function to access information it needs.
</P>
<P>Another way around this restriction is to use the feature provided by
several of the library methods with callback function arguments which
allow you to pass in a pointer to whatever data you wish.  This
pointer is returned as an argument when the callback is made.  This
pointer should be cast to (void *) when passed in, and your callback
function can later cast it back to the appropriate data type.  For
example, a class could set the pointer to an array or an internal data
structure or the class itself as "(void *) this".  Specify a NULL if
your function doesn't need the pointer.
</P>
<H4>Python overhead 
</H4>
<P>Using the MR-MPI library from Python incurs two not-so-obvious
overheads beyond the usual slowdown due to using an interpreted
language.  First, Python objects used as keys and values are "pickled"
and "unpickled" using the cPickle Python library when passed into and
out of the C++ library.  This is because the library stores them as
byte strings.  The pickling process serializes a Python object
(e.g. an integer, a string, a tuple, or a list) into a byte stream in
a way that it can be unpickled into the same Python object.
</P>
<P>The second overhead is due to the complexity of making a double
callbacks between the library and your Python script.  I.e. the
library calls back once to the user program which then calls back into
the library.  Consider what happens during a map() operation when the
library is called from a C++ program.
</P>
<UL><LI>the program calls the library map() method
<LI>the library map() calls back to the user map() callback function
<LI>the user map() calls the library add() method to register a key/value pair 
</UL>
<P>When doing this from Python there are 3 additional layers between the
Python program and the library, the Python mrmpi class, an invisible C
layer (created by ctypes), and the C interface on the C++ library
itself.  Thus the callback operation proceeds as follows:
</P>
<UL><LI>the program calls the mrmpi class map() method
<LI>the mrmpi class map() calls the invisible C map() function
<LI>the invisible map() calls the C interface map() function
<LI>the C interface map() calls the library map() method
<LI>the library map() calls back to the invisible C callback function
<LI>the invisible callback calls the mrmpi class callback method
<LI>the mrmpi callback calls the user map() callback function
<LI>the user map() calls the mrmpi class add() method to register a key/value pair
<LI>the mrmpi class add() calls the invisible C add() function
<LI>the invisible add() calls the C interface add() function
<LI>the C interface add() calls the library add() method 
</UL>
<P>Thus 3 calls have become 11 due to the 3 additional layers data must
pass through.  Some of these pass throughs are very simple, but others
require massaging and copying of data, like the pickling/unpickling
described above, which occurs in the mrmpi class methods.  I was
somewhat surprised this double-callback sequence works as well and as
transparently as it does - Python ctypes is amazing!
</P>
<H4>Error messages 
</H4>
<P>The error messages printed out by the MR-MPI library are hopefully
self-explanatory.  At some point they will be listed in these doc
pages.
</P>
<HR>

<A NAME = "examples"></A><H3>Examples 
</H3>
<P>This section describes the MapReduce programs provided in the examples
directory of the distribution:
</P>
<UL><LI>wordfreq
<LI>rmat 
</UL>
<P>Each are provided in 3 formats: as a C++ program, C program, and
Python script.  Note that the Python scripts use the PyPar package
which provides a Python/MPI interface, as discussed above in the
<A HREF = "#pyinterface">Python Interface</A> section, so you must have
<A HREF = "http://datamining.anu.edu.au/~ole/pypar">PyPar</A> installed in your Python to run them.
</P>
<P>The C++ and C programs can be built (assuming you have already built
the MR-MPI library) by typing
</P>
<PRE>make -f Makefile.foo 
</PRE>
<P>from within the examples directory, using one of the provided
Makefiles.  As with the library itself, you may need to edit one of
the Makefiles to create a new version appropriate to your machine.
</P>
<A NAME = "word"></A><H4>Word frequency example 
</H4>
<P>The wordfreq programs implement the word frequency counting algorithm
described above in <A HREF = "#whatis">this section</A>.  The wordfreq programs are
run by giving a list of text files as arguments, e.g.
</P>
<PRE>wordfreq ~/mydir/*.cpp
mpirun -np 8 cwordfreq ~/mydir/*.cpp
python wordfreq.py ~/mydir/*.cpp
mpirun -np 8 python wordfreq.py ~/mydir/*.cpp 
</PRE>
<P>Total word counts and a list of the top 10 words should be printed to
the screen, along with the time to perform the operation.
</P>
<P>The 3 different versions of the wordfreq program should give the same
answers, although if non-text files are used, the parsing of the
contents into words can be done differently by the C library strtok()
function and the Python string "split" method.
</P>
<A NAME = "rmat"></A><H4>R-MAT matrices example 
</H4>
<P>The rmat programs generate a particular form of randomized sparse
matrix known as an <A HREF = "#RMAT">R-MAT matrix</A>.  Depending on the parameters
chosen, the sparsity pattern in the resulting matrix can be highly
non-uniform, and a good model for irregular graphs, such as ones
representing a network of computers or WWW page links.
</P>
<P>The rmat programs are run by specifying a few parameters, e.g.
</P>
<PRE>rmat N Nz a b c d frac outfile
mpirun -np 8 crmat N Nz a b c d frac outfile
python rmat.py N Nz a b c d frac outfile
mpirun -np 8 python rmat.py N Nz a b c d frac outfile 
</PRE>
<P>The meaning of the parameters is as follows.  Note that only matrices
with a power-of-2 number of rows can be generated, so specifying N=20
creates a matrix with over a million rows.
</P>
<UL><LI>2^N = # of rows in matrix
<LI>Nz = average # of non-zeroes per row
<LI>a,b,c,d = generation params for matrix entries, must sum to 1
<LI>frac = randomization parameter between 0 and 1
<LI>seed = random # seed, positive integer
<LI>outfile = optional output file 
</UL>
<P>A full description of the R-MAT generation algorithm is beyond the
scope of this doc page, but here's the brief version.  The <I>a,b,c,d</I>
parameters are effectively weights on the 4 quadrants of the matrix.
To generate a single new matrix element, one quadrant is chosen, with
a probability proportional to its weight.  This operation is repeated
recursively within the chosen quadrant, applying the <I>frac</I> parameter
to randomize the weights a bit.  After N iterations, a single I,J
matrix location has been identified and its value is set (to 1 in this
case).
</P>
<P>The total number of matrix entries generated is Nx * 2^N.  This
procedure can generate duplicates, so those are removed, and new
elements generated until the desired number is reached.
</P>
<P>When completed, the matrix statistics are printed to the screen, along
with the time to generate the matrix.  If the optional <I>outfile</I>
parameter is specified, then the matrix entries are written to files
(one per processor).  Each line of any file has the form
</P>
<PRE>I J value 
</PRE>
<P>where I,J are the matrix row,column and value is the matrix entry (all
are 1 in this case).  If the files are concatenated together,
the full set of matrix entries should result.
</P>
<P>The 3 different versions of the rmat programs should give the same
answers in a statistical sense.  The answers will not be identical
because the same random number generation scheme is not used in all 3
programs.
</P>
<HR>

<H4>Citations 
</H4>
<A NAME = "Dean"></A>

<P><B>(Dean)</B> J. Dean and S. Ghemawat, "MapReduce: Simplified Data
Processing on Large Clusters", OSDI'04 conference (2004); J. Dean and
S. Ghemawat, "MapReduce: Simplified Data Processing on Large
Clusters", Communications of the ACM, 51, p 107-113 (2008).
</P>
<A NAME = "Tu"></A>

<P><B>(Tu)</B> T. Tu, C. A. Rendleman, D. W. Borhani, R. O. Dror,
J. Gullingsrud, M. O. Jensen, J. L. Kelpeis, P. Maragakis, P. Miller,
K. A. Stafford, D. E. Shaw, "A Scalable Parallel Framework for
Analyzing Terascale Molecular Dynamics Trajectories", SC08 proceedings
(2008).
</P>
<A NAME = "Gray"></A>

<P><B>(Gray)</B> A. Gray, Georgia Tech, http://www.cc.gatech.edu/~agray
</P>
<A NAME = "RMAT"></A>

<P><B>(RMAT)</B> D. Chakrabarti, Y. Zhan, C. Faloutsos, R-MAT: A Recursive
Model for Graph Mining", if Proceedings of the SIAM Conference on Data
Mining (2004), available at
http://www.cs.cmu.edu/~deepay/mywww/papers/siam04.pdf.
</P>
<P>Alexandaer Gray, Georgia Tech, http://www.cc.gatech.edu/~agray
</P>














</BODY>

</HTML>