This file is indexed.

/usr/include/deal.II/base/data_out_base.h is in libdeal.ii-dev 6.3.1-1.1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

   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
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
2529
2530
2531
2532
2533
2534
2535
2536
2537
2538
2539
2540
2541
2542
2543
2544
2545
2546
2547
2548
2549
2550
2551
2552
2553
2554
2555
2556
2557
2558
2559
2560
2561
2562
2563
2564
2565
2566
2567
2568
2569
2570
2571
2572
2573
2574
2575
2576
2577
2578
2579
2580
2581
2582
2583
2584
2585
2586
2587
2588
2589
2590
2591
2592
2593
2594
2595
2596
2597
2598
2599
2600
2601
2602
2603
2604
2605
2606
2607
2608
2609
2610
2611
2612
2613
2614
2615
2616
2617
2618
2619
2620
2621
2622
2623
2624
2625
2626
2627
2628
2629
2630
2631
2632
2633
2634
2635
2636
2637
2638
2639
2640
2641
2642
2643
2644
2645
2646
2647
2648
2649
2650
2651
2652
2653
2654
2655
2656
2657
2658
2659
2660
2661
2662
2663
2664
2665
2666
2667
2668
2669
2670
2671
2672
2673
2674
2675
2676
2677
2678
2679
2680
2681
2682
2683
2684
2685
2686
2687
2688
2689
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701
2702
2703
2704
2705
2706
2707
2708
2709
2710
2711
2712
2713
2714
2715
2716
2717
2718
2719
2720
2721
2722
2723
2724
2725
2726
2727
2728
2729
2730
2731
2732
2733
2734
2735
2736
2737
2738
2739
2740
2741
2742
2743
2744
2745
2746
2747
2748
2749
2750
2751
2752
2753
2754
2755
2756
2757
2758
2759
2760
2761
2762
2763
2764
2765
2766
2767
2768
2769
2770
2771
2772
2773
2774
2775
2776
2777
2778
2779
2780
2781
2782
2783
2784
2785
2786
2787
2788
2789
2790
2791
2792
2793
2794
2795
2796
2797
2798
2799
2800
2801
2802
2803
2804
2805
2806
2807
2808
2809
2810
2811
2812
2813
2814
2815
2816
2817
2818
2819
2820
2821
2822
2823
2824
2825
2826
2827
2828
2829
2830
2831
2832
2833
2834
2835
2836
2837
2838
2839
2840
2841
2842
2843
2844
2845
2846
2847
2848
2849
2850
2851
2852
2853
2854
2855
2856
2857
2858
2859
2860
2861
2862
2863
2864
2865
2866
2867
2868
2869
2870
2871
2872
2873
2874
2875
2876
2877
2878
2879
2880
2881
2882
2883
2884
2885
2886
2887
2888
2889
2890
2891
2892
2893
2894
2895
//---------------------------------------------------------------------------
//    $Id: data_out_base.h 21203 2010-06-11 19:48:01Z bangerth $
//    Version: $Name$
//
//    Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors
//
//    This file is subject to QPL and may not be  distributed
//    without copyright and license information. Please refer
//    to the file deal.II/doc/license.html for the  text  and
//    further information on this license.
//
//---------------------------------------------------------------------------
#ifndef __deal2__data_out_base_h
#define __deal2__data_out_base_h


#include <base/config.h>
#include <base/point.h>
#include <base/table.h>
#include <base/geometry_info.h>
#include <base/std_cxx1x/tuple.h>

#include <vector>
#include <string>

// Only include the Tecplot API header if the appropriate files
// were detected by configure
#ifdef DEAL_II_HAVE_TECPLOT
#  include "TECIO.h"
#  include <string.h>
#endif

// we only need output streams, but older compilers did not provide
// them in a separate include file
#ifdef HAVE_STD_OSTREAM_HEADER
#  include <ostream>
#else
#  include <iostream>
#endif

DEAL_II_NAMESPACE_OPEN


class ParameterHandler;


/**
 * This is a base class for output of data on meshes of very general
 * form. Output data is expected as a set of <tt>patches</tt> and
 * written to the output stream in the format expected by the
 * visualization tool. For a list of output formats, check the
 * enumeration #OutputFormat. For each format listed there, this class
 * contains a function <tt>write_format</tt>, writing the
 * output. Refer to the documentation of those functions for details
 * on a certain format.
 *
 * <h3>Structure of the output data</h3>
 *
 * Data is not written with the deal.II mesh structure. Instead, it
 * relies on a set of <tt>patches</tt> created by a derived class (for
 * example the DataOut, DataOutStack, DataOutFaces, DataOutRotation,
 * or MatrixOut classes).  Each Patch describes a single logical cell
 * of a mesh, possibly subdivided a number of times to represent
 * higher order polynomials defined on this cell. To this end, a patch
 * consists of a <tt>dim</tt>-dimensional regular grid with the same
 * number of grid points in each direction. In the simplest case it
 * may consist of the corner points of a single mesh cell.  For each
 * point of this local grid, the Patch contains an arbitrary number of
 * data values, though the number of data sets must be the same for
 * each point on each patch.
 *
 * By offering this interface to the different output formats, it is simple
 * to extend this class to new formats without depending on such things
 * as actual triangulations and handling of data vectors. These things shall
 * be provided by derived class which have a user callable interface then.
 *
 * Inside each patch, the data is organized in the usual
 * lexicographical order, <i>x</i> running fastest, then <i>y</i> and
 * <i>z</i>. Nodes are stored in this order and cells as well. Each
 * cell in 3D is stored such that the front face is in the
 * <i>xz</i>-plane. In order to enhance intellegibility of this
 * concept, the following two sections are kept from a previous
 * version of this documentation.
 *
 *
 * <h4>Patches</h4>
 *
 * Grids can be thought of as a collection of cells; if you want to write out
 * data on such a grid, you can do so by writing them one cell at a time.
 * The functions in this class therefore take a list of objects describing the
 * data on one cell each. This data for each cell usually consists of a list
 * of vertices for this cell, and a list of data values (for example solution
 * data, error information, etc) at each of these vertices.
 *
 * In some cases, this interface to a cell is too restricted, however. For
 * example, you may have higher order elements and printing the values at
 * the vertices only is not enough. For this reason, we not only provide
 * writing the data on the vertices only, but the data is organizes as a
 * tensor product grid on each cell. The parameter <tt>n_subdivision</tt>, which is
 * given for each patch separately, denotes how often the cell is to be
 * divided for output; for example, <tt>n_subdivisions==1</tt> yields no subdivision
 * of the cell, <tt>n_subdivisions==2</tt> will produce a grid of 3 times 3 points
 * in two spatial dimensions and 3 times 3 times 3 points in three dimensions,
 * <tt>n_subdivisions==3</tt> will yield 4 times 4 (times 4) points, etc. The actual
 * location of these points on the patch will be computed by a multilinear
 * transformation from the vertices given for this patch.

 * For cells at the boundary, a mapping might be used to calculate the position
 * of the inner points. In that case the coordinates are stored inside the
 * Patch, as they cannot be easily recovered otherwise.
 *
 * Given these comments, the actual data to be printed on this patch of
 * points consists of several data sets each of which has a value at each
 * of the patch points. For example with <tt>n_subdivisions==2</tt> in two space
 * dimensions, each data set has to provide nine values, and since the
 * patch is to be printed as a tensor product (or its transformation to the
 * real space cell), its values are to be ordered like
 * <i>(x0,y0) (x0,y1) (x0,y2) (x1,y0) (x1,y1) (x1,y2) (x2,y0) (x2,y1) (x2,y2)</i>,
 * i.e. the z-coordinate runs fastest, then the y-coordinate, then x (if there
 * are that many space directions).
 *
 *
 * <h4>Generalized patches</h4>
 *
 * In general, the patches as explained above might be too
 * restricted. For example, one might want to draw only the outer
 * faces of a domain in a three-dimensional computation, if one is not
 * interested in what happens inside. Then, the objects that should be
 * drawn are two-dimensional in a three-dimensional world. The
 * Patch class and associated output functions handle these
 * cases. The Patch class therefore takes two template parameters,
 * the first, named <tt>dim</tt> denoting the dimension of the object (in
 * the above example, this would be two), while the second, named
 * <tt>spacedim</tt>, denotes the dimension of the embedding space (this
 * would be three). The corner points of a patch have the dimension of
 * the space, while their number is determined by the dimension of the
 * patch. By default, the second template parameter has the same value
 * as the first, which would correspond to outputting a cell, rather
 * than a face or something else.
 *
 * <h3>DataOutBaseInterface</h3>
 *
 * This class has an interface that is not usually called by a user
 * directly; also, it consists of <tt>static</tt> functions
 * only. Usually, derived classes will inherit this class
 * <tt>protected</tt> to hide this interface to the users of thes
 * classes.
 *
 * The interface of this class basically consists of the declaration of a data
 * type describing a patch and a bunch of functions taking a list of patches
 * and writing them in one format or other to the stream. It is in the
 * responsibility of the derived classes to provide this list of patches.
 * In addition to the list of patches, a name for each data set may be given.
 *
 *
 * <h3>Querying interface</h3>
 *
 * This class also provides a few functions (parse_output_format(),
 * get_output_format_names(), default_suffix()) that can be used to query
 * which output formats this class supports. The provide a list of names for
 * all the formats we can output, parse a string and return an enum indicating
 * each format, and provide a way to convert a value of this enum into the
 * usual suffix used for files of that name. Using these functions, one can
 * entirely free applications from knowledge which formats the library
 * presently allows to output; several of the example programs show how to do
 * this.
 *
 * <h3>Output parameters</h3>
 *
 * All functions take a parameter which is a structure of type
 * <tt>XFlags</tt>, where <tt>X</tt> is the name of the output
 * format. To find out what flags are presently supported, read the
 * documentation of the different structures.
 *
 * Note that usually the output formats used for scientific visualization
 * programs have no or very few parameters (apart from some compatibility flags)
 * because there the actual appearance of output is determined using the
 * visualization program and the files produced by this class store more or less
 * only raw data.
 *
 * The direct output formats, like Postscript or Povray need to be given a lot
 * more parameters, though, since there the output file has to contain all
 * details of the viewpoint, light source, etc.
 *
 * <h3>Writing backends</h3>
 *
 * An abstraction layer has been introduced to facilitate coding
 * backends for additional visualization tools. It is applicable for
 * data formats separating the information into a field of vertices, a
 * field of connection information for the grid cells and data fields.
 *
 * For each of these fields, output functions are implemented, namely
 * write_nodes(), write_cells() and write_data(). In order to use
 * these functions, a format specific output stream must be written,
 * following the examples of DXStream, GmvStream, VtkStream and so on,
 * implemented in the .cc file.
 *
 * In this framework, the implementation of a new output format is
 * reduced to writing the section headers and the new output stream
 * class for writing a single mesh object.
 *
 * <h3>Credits</h3>
 * <ul>
 *
 * <li>EPS output based on an earlier implementation by Stefan Nauber
 * for the old DataOut class
 *
 * <li>@ref SoftwarePovray output by Thomas Richter
 *
 * <li>@ref SoftwareTecplot output by Benjamin Shelton Kirk
 *
 * </ul>
 *
 * @ingroup output
 * @author Wolfgang Bangerth, Guido Kanschat 1999, 2000, 2001, 2002, 2005, 2006.
 */
class DataOutBase
{
  public:
				     /**
				      * Data structure describing a patch of
				      * data in <tt>dim</tt> space
				      * dimensions.
				      *
				      * A patch consists of the following data:
				      * <ul>
				      * <li>the corner #vertices,
				      * <li> the number
				      * #n_subdivisions of the number
				      * of cells the Patch has in each
				      * space direction,
				      * <li> the #data attached to
				      * each vertex, in the usual
				      * lexicographic ordering,
				      * <li> Information on #neighbors.
				      * </ul>
				      *
				      * See the general
				      * documentation of the
				      * <tt>DataOutBase</tt> class for more
				      * information on its contents and
				      * purposes.  In the case of two
				      * dimensions, the next picture ist an
				      * example of <tt>n_subdivision</tt> = 4
				      * because the number of (sub)cells
				      * within each patch is equal to
				      * <tt>2^dim</tt>.
				      *
				      * @ingroup output
				      *
				      * @author Wolfgang Bangerth, Guido Kanschat
				      */
    template <int dim, int spacedim=dim>
    struct Patch
    {
					 /**
					  * Make the <tt>spacedim</tt> template
					  * parameter available.
					  */
	static const unsigned int space_dim=spacedim;

					 /**
					  * Corner points of a patch.
					  * Inner points are computed by
					  * a multilinear transform of
					  * the unit cell to the cell
					  * specified by these corner
					  * points. The order of points
					  * is the same as for cells
					  * in the triangulation.
					  */
	Point<spacedim> vertices[GeometryInfo<dim>::vertices_per_cell];

					 /**
					  * Numbers of neighbors of a patch.
					  * OpenDX format requires
					  * neighbor information for
					  * advanced output. Here the
					  * neighborship relationship
					  * of patches is
					  * stored. During output,
					  * this must be transformed
					  * into neighborship of
					  * sub-grid cells.
					  */
	unsigned int neighbors[GeometryInfo<dim>::faces_per_cell];

					 /**
					  * Number of this
					  * patch. Since we are not
					  * sure patches are handled
					  * in the same order, always,
					  * we better store this.
					  */
	unsigned int patch_index;

					 /**
					  * Number of subdivisions with
					  * which this patch is to be
					  * written. <tt>1</tt> means no
					  * subdivision, <tt>2</tt> means
					  * bisection, <tt>3</tt> trisection,
					  * etc.
					  */
	unsigned int n_subdivisions;

					 /**
					  * Data vectors. The format is
					  * as follows:
					  * <tt>data(i,.)</tt> denotes the data
					  * belonging to the <tt>i</tt>th data
					  * vector. <tt>data.n()</tt>
					  * therefore equals the number
					  * of output points; this
					  * number is <tt>(subdivisions+1)^{dim</tt>}.
					  * <tt>data.m()</tt> equals the number of
					  * data vectors.
					  *
					  * Within each column,
					  * <tt>data(.,j)</tt> are the
					  * data values at the output
					  * point <tt>j</tt>, where
					  * <tt>j</tt> denotes the
					  * usual lexicographic
					  * ordering in deal.II. This
					  * is also the order of
					  * points as provided by the
					  * <tt>QIterated</tt> class
					  * when used with the
					  * <tt>QTrapez</tt> class as
					  * subquadrature.
					  *
					  * Since the number of data vectors
					  * is usually the same for all
					  * patches to be printed,
					  * <tt>data.size()</tt> should yield
					  * the same value for all patches
					  * provided. The exception are
					  * patches for which
					  * points_are_available are set,
					  * where the actual coordinates of
					  * the point are appended to the
					  * 'data' field, see the
					  * documentation of the
					  * points_are_available flag.
					  */
	Table<2,float> data;

					 /**
					  * Bool flag indicating, whether the
					  * coordinates of the inner patch
					  * points are appended to the @p data
					  * table (@ true) or not (@ false),
					  * where the second is the standard and
					  * can be found for all cells in the
					  * interior of a domain.
					  *
					  * On the boundary of a domain, patch
					  * points are evaluated using a
					  * Mapping and therefore have to be
					  * stored inside the patch, as the
					  * Mapping and the corresponding
					  * boundary information are no longer
					  * available later on when we
					  * actually write the patch out to an
					  * output stream.
					  */
	bool points_are_available;

					 /**
					  * Default constructor. Sets
					  * <tt>n_subdivisions</tt> to one.
					  */
	Patch ();

					 /**
					  * Compare the present patch
					  * for equality with another
					  * one. This is used in a few
					  * of the automated tests in
					  * our testsuite.
					  */
	bool operator == (const Patch &patch) const;

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;

					 /**
					  * Value to be used if this
					  * patch has no neighbor on
					  * one side.
					  */
	static const unsigned int no_neighbor = numbers::invalid_unsigned_int;
					 /** @addtogroup Exceptions
					  * @{ */

					 /**
					  * Exception
					  */
	DeclException2 (ExcInvalidCombinationOfDimensions,
			int, int,
			<< "It is not possible to have a structural dimension of " << arg1
			<< " to be larger than the space dimension of the surrounding"
			<< " space " << arg2);
					 //@}
    };

    				     /**
				      * Flags controlling the details of
				      * output in @ref SoftwareOpenDX format.
				      *
				      * @ingroup output
				      */
    struct DXFlags
    {
					 /**
					  * Write neighbor
					  * information. This
					  * information is necessary
					  * for instance, if OpenDX is
					  * supposed to compute
					  * integral curves
					  * (streamlines). If it is
					  * not present, streamlines
					  * end at cell boundaries.
					  */
	bool write_neighbors;
					 /**
					  * Write integer values of
					  * the Triangulation in
					  * binary format.
					  */
	bool int_binary;
					 /**
					  * Write coordinate vectors in
					  * binary format.
					  */
	bool coordinates_binary;

					 /**
					  * Write data vectors in
					  * binary format.
					  */
	bool data_binary;

					 /**
					  * Write binary coordinate
					  * vectors as double (64 bit)
					  * numbers instead of float (32 bit).
					  */
	bool data_double;

					 /**
					  * Constructor.
					  */
	DXFlags (const bool write_neighbors = false,
		 const bool int_binary = false,
		 const bool coordinates_binary = false,
		 const bool data_binary = false);

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm);

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object.
					  */
	unsigned int memory_consumption () const;
    };

				     /**
				      * Flags controlling the details
				      * of output in UCD format for
				      * @ref SoftwareAVS.
				      *
				      * @ingroup output
				      */
    struct UcdFlags
    {
					 /**
					  * Write a comment at the
					  * beginning of the file
					  * stating the date of
					  * creation and some other
					  * data.  While this is
					  * supported by the UCD
					  * format and @ref
					  * SoftwareAVS, some other
					  * programs get confused by
					  * this, so the default is to
					  * not write a
					  * preamble. However, a
					  * preamble can be written
					  * using this flag.
					  *
					  * Default: <code>false</code>.
					  */
	bool write_preamble;

					 /**
					  * Constructor.
					  */
	UcdFlags (const bool write_preamble = false);

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm);

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };

				     /**
				      * Flags controlling the details of
				      * output in Gnuplot format. At
				      * present no flags are implemented.
				      *
				      * @ingroup output
				      */
    struct GnuplotFlags
    {
      private:
					 /**
					  * Dummy entry to suppress compiler
					  * warnings when copying an empty
					  * structure. Remove this member
					  * when adding the first flag to
					  * this structure (and remove the
					  * <tt>private</tt> as well).
					  */
	int dummy;

      public:
					 /**
					  * Constructor.
					  */
	GnuplotFlags ();

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm) const;

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };

    				     /**
				      * Flags controlling the details
				      * of output in @ref SoftwarePovray
				      * format. Several flags are
				      * implemented, see their
				      * respective documentation.
				      *
				      * @ingroup output
				      */
    struct PovrayFlags
    {
					 /**
					  * Normal vector interpolation,
					  * if set to true
					  *
					  * default = false
					  */
	bool smooth;

					 /**
					  * Use bicubic patches (b-splines)
					  * instead of triangles.
					  *
					  * default = false
					  */
	bool bicubic_patch;

					 /**
					  * include external "data.inc"
					  * with camera, light and
					  * texture definition for the
					  * scene.
					  *
					  * default = false
					  */
	bool external_data;

					 /**
					  * Constructor.
					  */
	PovrayFlags (const bool smooth = false,
		     const bool bicubic_patch = false,
		     const bool external_data = false);

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm);

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };


				     /**
				      * Flags controlling the details of
				      * output in encapsulated postscript
				      * format.
				      *
				      * @ingroup output
				      */
    struct EpsFlags
    {
					 /**
					  * This denotes the number of the
					  * data vector which shall be used
					  * for generating the height
					  * information. By default, the
					  * first data vector is taken,
					  * i.e. <tt>height_vector==0</tt>, if
					  * there is any data vector. If there
					  * is no data vector, no height
					  * information is generated.
					  */
	unsigned int height_vector;

					 /**
					  * Number of the vector which is
					  * to be taken to colorize cells.
					  * The same applies as for
					  * <tt>height_vector</tt>.
					  */
	unsigned int color_vector;

					 /**
					  * Enum denoting the possibilities
					  * whether the scaling should be done
					  * such that the given <tt>size</tt> equals
					  * the width or the height of
					  * the resulting picture.
					  */
	enum SizeType {
					       /// Scale to given width
	      width,
					       /// Scale to given height
	      height
	};

					 /**
					  * See above. Default is <tt>width</tt>.
					  */
	SizeType size_type;

					 /**
					  * Width or height of the output
					  * as given in postscript units
					  * This usually is given by the
					  * strange unit 1/72 inch. Whether
					  * this is height or width is
					  * specified by the flag
					  * <tt>size_type</tt>.
					  *
					  * Default is 300, which represents
					  * a size of roughly 10 cm.
					  */
	unsigned int size;

					 /**
					  * Width of a line in postscript
					  * units. Default is 0.5.
					  */
	double line_width;

					 /**
					  * Angle of the line origin-viewer
					  * against the z-axis in degrees.
					  *
					  * Default is the Gnuplot-default
					  * of 60.
					  */
	double azimut_angle;

					 /**
					  * Angle by which the viewers
					  * position projected onto the
					  * x-y-plane is rotated around
					  * the z-axis, in positive sense
					  * when viewed from above. The
					  * unit are degrees, and zero
					  * equals a position above or below
					  * the negative y-axis.
					  *
					  * Default is the
					  * Gnuplot-default of 30.
					  * An exemple of a
					  * Gnuplot-default of 0 is
					  * the following:
					  *
					  * @verbatim
					  *
					  *          3________7
					  *          /       /|
					  *         /       / |
					  *       2/______6/  |
					  *       |   |   |   |
					  * O-->  |   0___|___4
					  *       |  /    |  /
					  *       | /     | /
					  *      1|/______5/
					  *
					  * @endverbatim
					  */
	double turn_angle;

					 /**
					  * Factor by which the z-axis is to
					  * be stretched as compared to the
					  * x- and y-axes. This is to compensate
					  * for the different sizes that
					  * coordinate and solution values may
					  * have and to prevent that the plot
					  * looks to much out-of-place (no
					  * elevation at all if solution values
					  * are much smaller than coordinate
					  * values, or the common "extremely
					  * mountainous area" in the opposite
					  * case.
					  *
					  * Default is <tt>1.0</tt>.
					  */
	double z_scaling;

					 /**
					  * Flag the determines whether the
					  * lines bounding the cells (or the
					  * parts of each patch) are to be
					  * plotted.
					  *
					  * Default: <tt>true</tt>.
					  */
	bool   draw_mesh;

					 /**
					  * Flag whether to fill the regions
					  * between the lines bounding the cells
					  * or not. If not, no hidden line removal
					  * is performed, which in this crude
					  * implementation is done through
					  * writing the cells in a back-to-front
					  * order, thereby hiding the cells in
					  * the background by cells in the
					  * foreground.
					  *
					  * If this flag is <tt>false</tt> and <tt>draw_mesh</tt>
					  * is <tt>false</tt> as well, nothing will be
					  * printed.
					  *
					  * If this flag is <tt>true</tt>, then the cells
					  * will be drawn either colored by one
					  * of the data sets (if <tt>shade_cells</tt> is
					  * <tt>true</tt>), or pure white (if
					  * <tt>shade_cells</tt> is false or if there are
					  * no data sets).
					  *
					  * Default is <tt>true</tt>.
					  */
	bool   draw_cells;

					 /**
					  * Flag to determine whether the cells
					  * shall be colorized by the data
					  * set denoted by <tt>color_vector</tt>, or
					  * simply be painted in white. This
					  * flag only makes sense if
					  * <tt>draw_cells==true</tt>. Colorization is
					  * done through the <tt>color_function</tt>.
					  *
					  * Default is <tt>true</tt>.
					  */
	bool   shade_cells;

					 /**
					  * Structure keeping the three color
					  * values in the RGB system.
					  */
	struct RgbValues
	{
	    float red;
	    float green;
	    float blue;

					     /**
					      * Return <tt>true</tt> if the
					      * color represented by
					      * the three color values
					      * is a grey scale,
					      * i.e. all components
					      * are equal.
					      */
	    bool is_grey () const;
	};

					 /**
					  * Definition of a function pointer
					  * type taking a value and returning
					  * a triple of color values in RGB
					  * values.
					  *
					  * Besides the actual value by which
					  * the color is to be computed, min
					  * and max values of the data to
					  * be colorized are given as well.
					  */
	typedef RgbValues (*ColorFunction) (const double value,
					    const double min_value,
					    const double max_value);

					 /**
					  * This is a pointer to the function
					  * which is used to colorize the cells.
					  * By default, it points to the
					  * static function <tt>default_color_function</tt>
					  * which is a member of this class.
					  */
	ColorFunction color_function;


					 /**
					  * Default colorization function. This
					  * one does what one usually wants:
					  * It shifts colors from black (lowest
					  * value) through blue, green and red
					  * to white (highest value). For the
					  * exact defition of the color scale
					  * refer to the implementation.
					  *
					  * This function was originally written
					  * by Stefan Nauber.
					  */
	static RgbValues
        default_color_function (const double value,
                                const double min_value,
                                const double max_value);

					 /**
					  * This is an alternative color
					  * function producing a grey scale
					  * between black (lowest values)
					  * and white (highest values). You
					  * may use it by setting the
					  * <tt>color_function</tt> variable to the
					  * address of this function.
					  */
	static RgbValues
        grey_scale_color_function (const double value,
                                   const double min_value,
                                   const double max_value);

					 /**
					  * This is one more
					  * alternative color function
					  * producing a grey scale
					  * between white (lowest
					  * values) and black (highest
					  * values), i.e. the scale is
					  * reversed to the previous
					  * one. You may use it by
					  * setting the
					  * <tt>color_function</tt>
					  * variable to the address of
					  * this function.
					  */
	static RgbValues
        reverse_grey_scale_color_function (const double value,
                                           const double min_value,
                                           const double max_value);

					 /**
					  * Constructor.
					  */
	EpsFlags (const unsigned int  height_vector = 0,
		  const unsigned int  color_vector  = 0,
		  const SizeType      size_type     = width,
		  const unsigned int  size          = 300,
		  const double        line_width    = 0.5,
		  const double        azimut_angle  = 60,
		  const double        turn_angle    = 30,
		  const double        z_scaling     = 1.0,
		  const bool          draw_mesh     = true,
		  const bool          draw_cells    = true,
		  const bool          shade_cells   = true,
		  const ColorFunction color_function= &default_color_function);

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  *
					  * For coloring, only the color
					  * functions declared in this
					  * class are offered.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm);

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };

    				     /**
				      * Flags controlling the details
				      * of output in @ref SoftwareGMV
				      * format. At present no flags
				      * are implemented.
				      *
				      * @ingroup output
				      */
    struct GmvFlags
    {
      private:
					 /**
					  * Dummy entry to suppress compiler
					  * warnings when copying an empty
					  * structure. Remove this member
					  * when adding the first flag to
					  * this structure (and remove the
					  * <tt>private</tt> as well).
					  */
	int dummy;

      public:
					 /**
					  * Default constructor.
					  */
	GmvFlags ();

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm) const;

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };

    				     /**
				      * Flags controlling the details
				      * of output in @ref
				      * SoftwareTecplot format.
				      *
				      * @ingroup output
				      */
    struct TecplotFlags
    {

      public:

					 /**
					  * This variable is needed to hold the
					  * output file name when using the
					  * Tecplot API to write binary files.
					  * If the user doesn't set the file
					  * name with this variable only
					  * ASCII Tecplot output will be
					  * produced.
					  */
        const char* tecplot_binary_file_name;

					 /**
					  * Tecplot allows to assign
					  * names to zones. This
					  * variable stores this name.
					  */
	const char* zone_name;

	                                 /**
	                                  * Constructor
	                                  **/
	TecplotFlags (const char* tecplot_binary_file_name = NULL,
		      const char* zone_name = NULL);

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm) const;

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };

    				     /**
				      * Flags controlling the details
				      * of output in @ref SoftwareVTK
				      * format. At present no flags
				      * are implemented.
				      *
				      * @ingroup output
				      */
    struct VtkFlags
    {
      private:
					 /**
					  * Dummy entry to suppress compiler
					  * warnings when copying an empty
					  * structure. Remove this member
					  * when adding the first flag to
					  * this structure (and remove the
					  * <tt>private</tt> as well).
					  */
	int dummy;

      public:
					 /**
					  * Default constructor.
					  */
	VtkFlags ();

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm) const;

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };
    				     /**
				      * Flags controlling the details
				      * of output in deal.II
				      * intermediate format. At
				      * present no flags are
				      * implemented.
				      *
				      * @ingroup output
				      */
    struct Deal_II_IntermediateFlags
    {
					 /**
					  * An indicator of the
					  * currect file format
					  * version used to write
					  * intermediate format. We do
					  * not attempt to be backward
					  * compatible, so this number
					  * is used only to verify
					  * that the format we are
					  * writing is what the
					  * current readers and
					  * writers understand.
					  */
        static const unsigned int format_version = 3;
      private:
					 /**
					  * Dummy entry to suppress compiler
					  * warnings when copying an empty
					  * structure. Remove this member
					  * when adding the first flag to
					  * this structure (and remove the
					  * <tt>private</tt> as well).
					  */
	int dummy;

      public:
					 /**
					  * Constructor.
					  */
	Deal_II_IntermediateFlags ();

					 /**
					  * Declare all flags with name
					  * and type as offered by this
					  * class, for use in input files.
					  */
	static void declare_parameters (ParameterHandler &prm);

					 /**
					  * Read the parameters declared in
					  * <tt>declare_parameters</tt> and set the
					  * flags for this output format
					  * accordingly.
					  *
					  * The flags thus obtained overwrite
					  * all previous contents of this object.
					  */
	void parse_parameters (const ParameterHandler &prm) const;

					 /**
					  * Determine an estimate for
					  * the memory consumption (in
					  * bytes) of this
					  * object. Since sometimes
					  * the size of objects can
					  * not be determined exactly
					  * (for example: what is the
					  * memory consumption of an
					  * STL <tt>std::map</tt> type with a
					  * certain number of
					  * elements?), this is only
					  * an estimate. however often
					  * quite close to the true
					  * value.
					  */
	unsigned int memory_consumption () const;
    };

				     /**
				      * Provide a data type specifying
				      * the presently supported output
				      * formats.
				      */
    enum OutputFormat
    {
					   /**
					    * Use the format already
					    * stored in the object.
					    */
	  default_format,
					   /**
					    * Do not write any output.
					    */
	  none,
					   /**
					    * Output for @ref SoftwareOpenDX.
					    */
	  dx,
					   /**
					    * Output in the UCD format
					    * for @ref SoftwareAVS.
					    */
	  ucd,
					   /**
					    * Output for the @ref
					    * SoftwareGnuplot tool.
					    */
	  gnuplot,
					   /**
					    * Output for the @ref
					    * SoftwarePovray
					    * raytracer.
					    */
	  povray,
					   /**
					    * Output in encapsulated
					    * PostScript.
					    */
	  eps,
					   /**
					    * Output for @ref SoftwareGMV.
					    */
	  gmv,
					   /**
					    * Output for @ref
					    * SoftwareTecplot in text
					    * format.
					    */

	  tecplot,
					   /**
					    * Output for @ref
					    * SoftwareTecplot in
					    * binary format. Faster
					    * and smaller than text
					    * format.
					    */
	  tecplot_binary,

					   /**
					    * Output in @ref
					    * SoftwareVTK format.
					    */
	  vtk,

					   /**
					    * Output in @ref
					    * SoftwareVTK format.
					    */
	  vtu,

					   /**
					    * Output in deal.II
					    * intermediate format.
					    */
	  deal_II_intermediate
    };


/**
 * Write the given list of patches to the output stream in @ref
 * SoftwareOpenDX format.
 *
 * Since @ref SoftwareOpenDX uses some kind of visual data flow oriented
 * programming language, some of these programs are provided in
 * <tt>contrib/dx</tt>.
 */
    template <int dim, int spacedim>
    static void write_dx (const std::vector<Patch<dim,spacedim> > &patches,
			  const std::vector<std::string>          &data_names,
			  const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			  const DXFlags                           &flags,
			  std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in eps format.
 *
 * Output in this format circumvents the use of auxiliary graphic
 * programs converting some output format into a graphics format. This
 * has the advantage that output is easy and fast, and the
 * disadvantage that you have to give a whole bunch of parameters
 * which determine the direction of sight, the mode of colorization,
 * the scaling of the height axis, etc. (Of course, all these
 * parameters have reasonable default values, which you may want to
 * change from time to time.) At present, this format only supports
 * output for two-dimensional data, with values in the third direction
 * taken from a data vector.
 *
 * Basically, output consists of the mesh and the cells in between
 * them. You can draw either of these, or both, or none if you are
 * really interested in an empty picture. If written, the mesh uses
 * black lines. The cells in between the mesh are either not printed
 * (this will result in a loss of hidden line removal, i.e.  you can
 * "see through" the cells to lines behind), printed in white (which
 * does nothing apart from the hidden line removal), or colorized
 * using one of the data vectors (which need not be the same as the
 * one used for computing the height information) and a customizable
 * color function. The default color functions chooses the color
 * between black, blue, green, red and white, with growing values of
 * the data field chosen for colorization. At present, cells are
 * displayed with one color per cell only, which is taken from the
 * value of the data field at the center of the cell; bilinear
 * interpolation of the color on a cell is not used.
 *
 * By default, the viewpoint is chosen like the default viewpoint in
 * GNUPLOT, i.e.  with an angle of 60 degrees with respect to the
 * positive z-axis and rotated 30 degrees in positive sense (as seen
 * from above) away from the negative y-axis.  Of course you can
 * change these settings.
 *
 * EPS output is written without a border around the picture, i.e. the
 * bounding box is close to the output on all four sides. Coordinates
 * are written using at most five digits, to keep picture size at a
 * reasonable size.
 *
 * All parameters along with their default values are listed in the
 * documentation of the <tt>EpsFlags</tt> member class of this
 * class. See there for more and detailed information.
 */
    template <int dim, int spacedim>
    static void write_eps (const std::vector<Patch<dim,spacedim> > &patches,
			   const std::vector<std::string>          &data_names,
			   const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			   const EpsFlags                          &flags,
			   std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in @ref
 * SoftwareGMV format.
 *
 * Data is written in the following format: nodes are considered the
 * points of the patches. In spatial dimensions less than three,
 * zeroes are inserted for the missing coordinates. The data vectors
 * are written as node or cell data, where for the first the data
 * space is interpolated to (bi-,tri-)linear elements.
 */
    template <int dim, int spacedim>
    static void write_gmv (const std::vector<Patch<dim,spacedim> > &patches,
			   const std::vector<std::string>          &data_names,
			   const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			   const GmvFlags                          &flags,
			   std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in gnuplot
 * format. Visualization of two-dimensional data can then be achieved by
 * starting <tt>gnuplot</tt> and endtering the commands
 *
 * @verbatim
 * set data style lines
 * splot "filename" using 1:2:n
 * @endverbatim
 * This example assumes that the number of the data vector displayed
 * is <b>n-2</b>.
 *
 * The GNUPLOT format is not able to handle data on unstructured grids
 * directly. Directly would mean that you only give the vertices and
 * the solution values thereon and the program constructs its own grid
 * to represent the data. This is only possible for a structured
 * tensor product grid in two dimensions. However, it is possible to
 * give several such patches within one file, which is exactly what
 * the respective function of this class does: writing each cell's
 * data as a patch of data, at least if the patches as passed from
 * derived classes represent cells. Note that the functions on patches
 * need not be continuous at interfaces between patches, so this
 * method also works for discontinuous elements. Note also, that
 * GNUPLOT can do hidden line removal for patched data.
 *
 * While this discussion applies to two spatial dimensions, it is more
 * complicated in 3d. The reason is that we could still use patches,
 * but it is difficult when trying to visualize them, since if we use
 * a cut through the data (by, for example, using x- and
 * z-coordinates, a fixed y-value and plot function values in
 * z-direction, then the patched data is not a patch in the sense
 * GNUPLOT wants it any more. Therefore, we use another approach,
 * namely writing the data on the 3d grid as a sequence of lines,
 * i.e. two points each associated with one or more data sets.  There
 * are therefore 12 lines for each subcells of a patch.
 *
 * Given the lines as described above, a cut through this data in
 * Gnuplot can then be achieved like this (& stands for the dollar
 * sign in the following):
 * @verbatim
 *   set data style lines
 *   splot [:][:][0:] "T" using 1:2:(&3==.5 ? &4 : -1)
 * @endverbatim
 *
 * This command plots data in x- and y-direction unbounded, but in
 * z-direction only those data points which are above the x-y-plane
 * (we assume here a positive solution, if it has negative values, you
 * might want to decrease the lower bound). Furthermore, it only takes
 * the data points with z-values (<tt>&3</tt>) equal to 0.5, i.e. a
 * cut through the domain at <tt>z=0.5</tt>. For the data points on
 * this plane, the data values of the first data set (<tt>&4</tt>) are
 * raised in z-direction above the x-y-plane; all other points are
 * denoted the value <tt>-1</tt> instead of the value of the data
 * vector and are not plotted due to the lower bound in z plotting
 * direction, given in the third pair of brackets.
 *
 * More complex cuts are possible, including nonlinear ones. Note
 * however, that only those points which are actually on the
 * cut-surface are plotted.
 */
    template <int dim, int spacedim>
    static void write_gnuplot (const std::vector<Patch<dim,spacedim> > &patches,
			       const std::vector<std::string>          &data_names,
			       const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			       const GnuplotFlags                      &flags,
			       std::ostream                            &out);

/**
 * Write the given list of patches to the output stream for the @ref
 * SoftwarePovray raytracer.
 *
 * Output in this format creates a povray source file, include
 * standard camera and light source definition for rendering with
 * povray 3.1 At present, this format only supports output for
 * two-dimensional data, with values in the third direction taken from
 * a data vector.
 *
 * The output uses two different povray-objects:
 *
 * <ul>
 * <li> <tt>BICUBIC_PATCH</tt>
 * A <tt>bicubic_patch</tt> is a 3-dimensional Bezier patch. It consists of 16 Points
 * describing the surface. The 4 corner points are touched by the object,
 * while the other 12 points pull and stretch the patch into shape.
 * One <tt>bicubic_patch</tt> is generated on each patch. Therefor the number of
 * subdivisions has to be 3 to provide the patch with 16 points.
 * A bicubic patch is not exact but generates very smooth images.
 *
 * <li> <tt>MESH</tt>
 * The mesh object is used to store large number of triangles.
 * Every square of the patch data is split into one upper-left and one
 * lower-right triangle. If the number of subdivisions is three, 32 triangle
 * are generated for every patch.
 *
 * Using the smooth flag povray interpolates the normals on the triangles,
 * imitating a curved surface
 * </ul>
 *
 * All objects get one texture definition called Tex. This texture has to be
 * declared somewhere before the object data. This may be in an external
 * data file or at the beginning of the output file.
 * Setting the <tt>external_data</tt> flag to false, an standard camera, light and
 * texture (scaled to fit the scene) is added to the outputfile. Set to true
 * an include file "data.inc" is included. This file is not generated by deal
 * and has to include camera, light and the texture definition Tex.
 *
 * You need povray (>=3.0) to render the scene. The minimum options for povray
 * are:
 * @verbatim
 *   povray +I<inputfile> +W<horiz. size> +H<ver. size> +L<include path>
 * @endverbatim
 * If the external file "data.inc" is used, the path to this file has to be
 * included in the povray options.
 */
    template <int dim, int spacedim>
    static void write_povray (const std::vector<Patch<dim,spacedim> > &patches,
			      const std::vector<std::string>          &data_names,
			      const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			      const PovrayFlags                       &flags,
			      std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in @ref
 * SoftwareTecplot ASCII format (FEBLOCK).
 *
 * For more information consult the Tecplot Users and Reference
 * manuals.
 */
    template <int dim, int spacedim>
    static void write_tecplot (const std::vector<Patch<dim,spacedim> > &patches,
			       const std::vector<std::string>          &data_names,
			       const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			       const TecplotFlags                      &flags,
			       std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in @ref
 * SoftwareTecplot binary format.
 *
 * For this to work properly <tt>./configure</tt> checks for the
 * Tecplot API at build time. To write Tecplot binary files directly
 * make sure that the TECHOME environment variable points to the
 * Tecplot installation directory, and that the files
 * @$TECHOME/include/TECIO.h and @$TECHOME/lib/tecio.a are readable.
 * If these files are not availabe (or in the case of 1D) this
 * function will simply call write_tecplot() and thus larger ASCII
 * data files will be produced rather than more efficient Tecplot
 * binary files.
 *
 * @warning TecplotFlags::tecplot_binary_file_name indicates the name
 * of the file to be written.  If the file name is not set ASCII
 * output is produced.
 *
 * For more information consult the Tecplot Users and Reference
 * manuals.
 */
    template <int dim, int spacedim>
    static void write_tecplot_binary (
      const std::vector<Patch<dim,spacedim> > &patches,
      const std::vector<std::string>          &data_names,
      const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
      const TecplotFlags                      &flags,
      std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in UCD format
 * described in the AVS developer's guide (now @ref SoftwareAVS). Due
 * to limitations in the present format, only node based data can be
 * output, which in one reason why we invented the patch concept. In
 * order to write higher order elements, you may split them up into
 * several subdivisions of each cell. These subcells will then,
 * however, also appear as different cells by programs which
 * understand the UCD format.
 *
 * No use is made of the possibility to give model data since these
 * are not supported by all UCD aware programs. You may give cell data
 * in derived classes by setting all values of a given data set on a
 * patch to the same value.
 */
    template <int dim, int spacedim>
    static void write_ucd (const std::vector<Patch<dim,spacedim> > &patches,
			   const std::vector<std::string>          &data_names,
			   const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			   const UcdFlags                          &flags,
			   std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in @ref SoftwareVTK
 * format. The data is written in the traditional VTK format as opposed to the
 * XML-based format that write_vtu() produces.
 *
 * The vector_data_ranges argument denotes ranges of components in the
 * output that are considered a vector, rather than simply a
 * collection of scalar fields. The VTK output format has special
 * provisions that allow these components to be output by a single
 * name rather than having to group several scalar fields into a
 * vector later on in the visualization program.
 */
    template <int dim, int spacedim>
    static void write_vtk (const std::vector<Patch<dim,spacedim> > &patches,
			   const std::vector<std::string>          &data_names,
			   const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			   const VtkFlags                          &flags,
			   std::ostream                            &out);


/**
 * Write the given list of patches to the output stream in @ref SoftwareVTK
 * format. The data is written in the XML-based VTK format as opposed to the
 * traditional format that write_vtk() produces.
 *
 * The vector_data_ranges argument denotes ranges of components in the
 * output that are considered a vector, rather than simply a
 * collection of scalar fields. The VTK output format has special
 * provisions that allow these components to be output by a single
 * name rather than having to group several scalar fields into a
 * vector later on in the visualization program.
 *
 * Some visualization programs, such as ParaView, can read several separate
 * VTU files to parallelize visualization. In that case, you need a
 * <code>.pvtu</code> file that describes which VTU files form a group. The
 * DataOutInterface::write_pvtu_record() function can generate such a master record.
 */
    template <int dim, int spacedim>
    static void write_vtu (const std::vector<Patch<dim,spacedim> > &patches,
			   const std::vector<std::string>          &data_names,
			   const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
			   const VtkFlags                          &flags,
			   std::ostream                            &out);

/**
 * Write the given list of patches to the output stream in deal.II
 * intermediate format. This is not a format understood by any other
 * graphics program, but is rather a direct dump of the intermediate
 * internal format used by deal.II. This internal format is generated
 * by the various classes that can generate output using the
 * DataOutBase class, for example from a finite element solution, and
 * is then converted in the present class to the final graphics
 * format.
 *
 * Note that the intermediate format is what its name suggests: a
 * direct representation of internal data. It isn't standardized and
 * will change whenever we change our internal representation. You can
 * only expect to process files written in this format using the same
 * version of deal.II that was used for writing.
 *
 * The reason why we offer to write out this intermediate format is
 * that it can be read back into a deal.II program using the
 * DataOutReader class, which is helpful in at least two contexts:
 * First, this can be used to later generate graphical output in any
 * other graphics format presently understood; this way, it is not
 * necessary to know at run-time which output format is requested, or
 * if multiple output files in different formats are needed. Secondly,
 * in contrast to almost all other graphics formats, it is possible to
 * merge several files that contain intermediate format data, and
 * generate a single output file from it, which may be again in
 * intermediate format or any of the final formats. This latter option
 * is most helpful for parallel programs: as demonstrated in the
 * step-17 example program, it is possible to let only one processor
 * generate the graphical output for the entire parallel program, but
 * this can become vastly inefficient if many processors are involved,
 * because the load is no longer balanced. The way out is to let each
 * processor generate intermediate graphical output for its chunk of
 * the domain, and the later merge the different files into one, which
 * is an operation that is much cheaper than the generation of the
 * intermediate data.
 *
 * Intermediate format deal.II data is usually stored in files with
 * the ending <tt>.d2</tt>.
 */
    template <int dim, int spacedim>
    static void write_deal_II_intermediate (
      const std::vector<Patch<dim,spacedim> > &patches,
      const std::vector<std::string>          &data_names,
      const std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> > &vector_data_ranges,
      const Deal_II_IntermediateFlags         &flags,
      std::ostream                            &out);


                                     /**
                                      * Given an input stream that contains
                                      * data written by
                                      * write_deal_II_intermediate, determine
                                      * the <tt>dim</tt> and <tt>spacedim</tt>
                                      * template parameters with which that
                                      * function was called, and return them
                                      * as a pair of values.
                                      *
                                      * Note that this function eats a number
                                      * of elements at the present position of
                                      * the stream, and therefore alters
                                      * it. In order to read from it using,
                                      * for example, the DataOutReader class,
                                      * you may wish to either reset the
                                      * stream to its previous position, or
                                      * close and reopen it.
                                      */
    static
    std::pair<unsigned int, unsigned int>
    determine_intermediate_format_dimensions (std::istream &input);

				     /**
				      * Return the <tt>OutputFormat</tt>
				      * value corresponding to the
				      * given string. If the string
				      * does not match any known
				      * format, an exception is
				      * thrown.
				      *
				      * Since this function does not
				      * need data from this object, it
				      * is static and can thus be
				      * called without creating an
				      * object of this class. Its main
				      * purpose is to allow a program
				      * to use any implemented output
				      * format without the need to
				      * extend the program's parser
				      * each time a new format is
				      * implemented.
				      *
				      * To get a list of presently
				      * available format names,
				      * e.g. to give it to the
				      * ParameterHandler class,
				      * use the function
				      * get_output_format_names().
				      */
    static OutputFormat parse_output_format (const std::string &format_name);

				     /**
				      * Return a list of implemented
				      * output formats. The different
				      * names are separated by
				      * vertical bar signs (<tt>`|'</tt>)
				      * as used by the
				      * ParameterHandler classes.
				      */
    static std::string get_output_format_names ();

				     /**
				      * Provide a function which tells us which
				      * suffix a file with a given output format
				      * usually has. At present the following
				      * formats are defined:
				      * <ul>
				      * <li> <tt>dx</tt>: <tt>.dx</tt>
				      * <li> <tt>ucd</tt>: <tt>.inp</tt>
				      * <li> <tt>gnuplot</tt>: <tt>.gnuplot</tt>
				      * <li> <tt>povray</tt>: <tt>.pov</tt>
				      * <li> <tt>eps</tt>: <tt>.eps</tt>
				      * <li> <tt>gmv</tt>: <tt>.gmv</tt>
				      * <li> <tt>vtk</tt>: <tt>.vtk</tt>
				      * <li> <tt>deal_II_intermediate</tt>: <tt>.d2</tt>.
				      * </ul>
				      */
    static std::string default_suffix (const OutputFormat output_format);

				     /**
				      * Determine an estimate for
				      * the memory consumption (in
				      * bytes) of this
				      * object. Since sometimes
				      * the size of objects can
				      * not be determined exactly
				      * (for example: what is the
				      * memory consumption of an
				      * STL <tt>std::map</tt> type with a
				      * certain number of
				      * elements?), this is only
				      * an estimate. however often
				      * quite close to the true
				      * value.
				      */
    static unsigned int memory_consumption ();

				     /** @addtogroup Exceptions
				      * @{ */

				     /**
				      * Exception
				      */
    DeclException2 (ExcInvalidDatasetSize,
		    int, int,
		    << "The number of points in this data set is " << arg1
		    << ", but we expected " << arg2 << " in each space direction.");
				     /**
				      * An output function did not
				      * receive any patches for
				      * writing.
				      */
    DeclException0 (ExcNoPatches);
				     /**
				       * Exception
				       */
    DeclException0 (ExcTecplotAPIError);
 				      /**
				       * Exception
				       */
    DeclException1 (ExcErrorOpeningTecplotFile,
		    char*,
		    << "There was an error opening Tecplot file " << arg1
		    << " for output");

				     //@}
  private:
				     /**
				      * Write the coordinates of nodes
				      * in the desired format.
				      */
    template <int dim, int spacedim, typename STREAM>
    static void write_nodes (const std::vector<Patch<dim,spacedim> >& patches,
			     STREAM& out);

				     /**
				      * Write the node numbers of a
				      * cell in the desired format.
				      */
    template <int dim, int spacedim, typename STREAM>
    static void write_cells (const std::vector<Patch<dim,spacedim> >& patches,
			     STREAM& out);

				     /**
				      * Write data in the desired
				      * format.
				      */
    template <int dim, int spacedim, class STREAM>
    static void write_data (const std::vector<Patch<dim,spacedim> >& patches,
			    const unsigned int n_data_sets,
			    const bool double_precision,
			    STREAM& out);

				     /**
				      * Class holding the data of one
				      * cell of a patch in two space
				      * dimensions for output. It is
				      * the projection of a cell in
				      * three-dimensional space (two
				      * coordinates, one height value)
				      * to the direction of sight.
				      */
    class EpsCell2d
    {
      public:

					 /**
					  * Vector of vertices of this cell.
					  */
	Point<2> vertices[4];

					 /**
					  * Data value from which the actual
					  * colors will be computed by
					  * the colorization function stated
					  * in the <tt>EpsFlags</tt> class.
					  */
	float color_value;

					 /**
					  * Depth into the picture, which
					  * is defined as the distance from
					  * an observer at an the origin in
					  * direction of the line of sight.
					  */
	float depth;

					 /**
					  * Comparison operator for
					  * sorting.
					  */
	bool operator < (const EpsCell2d &) const;
    };


				     /**
				      * This is a helper function for
				      * the <tt>write_gmv</tt>
				      * function. There, the data in
				      * the patches needs to be copied
				      * around as output is one
				      * variable globally at a time,
				      * rather than all data on each
				      * vertex at a time. This copying
				      * around can be done detached
				      * from the main thread, and is
				      * thus moved into this separate
				      * function.
				      *
				      * Note that because of the
				      * similarity of the formats,
				      * this function is also used by
				      * the Vtk and Tecplot output
				      * functions.
				      */
    template <int dim, int spacedim>
    static void
    write_gmv_reorder_data_vectors (const std::vector<Patch<dim,spacedim> > &patches,
				    Table<2,double>                         &data_vectors);

};




/**
 * This class is the interface to the <tt>DataOutBase</tt> class, as already its name
 * might suggest. It does not offer much functionality apart from a way to
 * access the implemented formats and a way to dynamically dispatch what output
 * format to chose.
 *
 * This class is thought as a base class to classes actually
 * generating data for output. It has two abstract virtual functions,
 * <tt>get_patches</tt> and <tt>get_dataset_names</tt> which are to produce the data
 * which is actually needed. These are the only functions that need to
 * be overloaded by a derived class.  In additional to that, it has a
 * function for each output format supported by the underlying base
 * class which gets the output data using these two virtual functions
 * and passes them to the raw output functions.
 *
 * The purpose of this class is mainly two-fold: to support storing flags
 * by which the output in the different output formats are controlled,
 * and means to work with output in a way where output format, flags and
 * other things are determined at run time. In addition to that it offers
 * the abstract interface to derived classes briefly discussed above.
 *
 *
 * <h3>Output flags</h3>
 *
 * The way we treat flags in this class is very similar to that used in
 * the <tt>GridOut</tt> class. For detailed information on the why's and how's,
 * as well as an example of programming, we refer to the documentation
 * of that class.
 *
 * Basically, this class stores a set of flags for each output format
 * supported by the underlying <tt>DataOutBase</tt> class. These are used
 * whenever one of the <tt>write_*</tt> functions is used. By default, the
 * values of these flags are set to reasonable start-ups, but in case
 * you want to change them, you can create a structure holding the flags
 * for one of the output formats and set it using the <tt>set_flags</tt> functions
 * of this class to determine all future output the object might produce
 * by that output format.
 *
 * For information on what parameters are supported by different output
 * functions, please see the documentation of the <tt>DataOutBase</tt> class and
 * its member classes.
 *
 *
 * <h3>Run time selection of output parameters</h3>
 *
 * In the output flags classes, described above, many flags are
 * defined for output in the different formats. In order to make them
 * available to the input file handler class <tt>ParameterHandler</tt>, each
 * of these has a function declaring these flags to the parameter
 * handler and to read them back from an actual input file. In order
 * to avoid that in user programs these functions have to be called
 * for each available output format and the respective flag class, the
 * present <tt>DataOutInterface</tt> class offers a function
 * <tt>declare_parameters</tt> which calls the respective function of all
 * known output format flags classes. The flags of each such format
 * are packed together in a subsection in the input file.
 * Likewise, there is a function <tt>parse_parameters</tt> which reads
 * these parameters and stores them in the flags associated with
 * this object (see above).
 *
 * Using these functions, you do not have to track which formats are
 * presently implemented.
 *
 * Usage is as follows:
 * @code
 *                               // within function declaring parameters:
 *   ...
 *   prm.enter_subsection ("Output format options");
 *     DataOutInterface<dim>::declare_parameters (prm);
 *   prm.leave_subsection ();
 *   ...
 *
 *
 *                               // within function doing the output:
 *   ...
 *   DataOut<dim> out;
 *   prm.enter_subsection ("Output format options");
 *   out.parse_parameters (prm);
 *   prm.leave_subsection ();
 *   ...
 * @endcode
 * Note that in the present example, the class <tt>DataOut</tt> was used. However, any
 * other class derived from <tt>DataOutInterface</tt> would work alike.
 *
 *
 * <h3>Run time selection of formats</h3>
 *
 * This class, much like the <tt>GridOut</tt> class, has a set of functions
 * providing a list of supported output formats, an <tt>enum</tt> denoting all
 * these and a function to parse a string and return the respective
 * <tt>enum</tt> value if it is a valid output format's name (actually, these
 * functions are inherited from the base class). Finally, there is a function
 * <tt>write</tt>, which takes a value of this <tt>enum</tt> and dispatches to
 * one of the actual <tt>write_*</tt> functions depending on the output format
 * selected by this value.
 *
 * The functions offering the different output format names are,
 * respectively, <tt>default_suffix</tt>, <tt>parse_output_format</tt>, and
 * <tt>get_output_format_names</tt>. They make the selection of ouput formats
 * in parameter files much easier, and especially independent of
 * the formats presently implemented. User programs need therefore not
 * be changed whenever a new format is implemented.
 *
 * Additionally, objects of this class have a default format, which
 * can be set by the parameter "Output format" of the parameter
 * file. Within a program, this can be changed by the member function
 * <tt>set_default_format</tt>. Using this default format, it is possible to leave
 * the format selection completely to the parameter file. A suitable
 * suffix for the output file name can be obtained by <tt>default_suffix</tt>
 * without arguments.
 *
 * @ingroup output
 * @author Wolfgang Bangerth, 1999
 */
template <int dim, int spacedim=dim>
class DataOutInterface : private DataOutBase
{
  public:
                                     /*
                                      * Import a few names that were
                                      * previously in this class and have then
                                      * moved to the base class. Since the
                                      * base class is inherited from
                                      * privately, we need to re-import these
                                      * symbols to make sure that references
                                      * to DataOutInterface<dim,spacedim>::XXX
                                      * remain valid.
                                      */
    using DataOutBase::OutputFormat;
    using DataOutBase::default_format;
    using DataOutBase::dx;
    using DataOutBase::gnuplot;
    using DataOutBase::povray;
    using DataOutBase::eps;
    using DataOutBase::tecplot;
    using DataOutBase::tecplot_binary;
    using DataOutBase::vtk;
    using DataOutBase::vtu;
    using DataOutBase::deal_II_intermediate;
    using DataOutBase::parse_output_format;
    using DataOutBase::get_output_format_names;
    using DataOutBase::determine_intermediate_format_dimensions;

				     /**
				      * Constructor.
				      */
    DataOutInterface ();

                                     /**
                                      * Destructor. Does nothing, but is
                                      * declared virtual since this class has
                                      * virtual functions.
                                      */
    virtual ~DataOutInterface ();

				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in OpenDX format. See
				      * DataOutBase::write_dx.
				      */
    void write_dx (std::ostream &out) const;

				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in EPS format. See
				      * DataOutBase::write_eps.
				      */
    void write_eps (std::ostream &out) const;

    				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in GMV format. See
				      * DataOutBase::write_gmv.
				      */
    void write_gmv (std::ostream &out) const;

				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in GNUPLOT format. See
				      * DataOutBase::write_gnuplot.
				      */
    void write_gnuplot (std::ostream &out) const;

    				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in POVRAY format. See
				      * DataOutBase::write_povray.
				      */
    void write_povray (std::ostream &out) const;

    				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in Tecplot format. See
				      * DataOutBase::write_tecplot.
				      */
    void write_tecplot (std::ostream &out) const;

    				     /**
				      * Obtain data through
				      * get_patches() and write it in
				      * the Tecplot binary output
				      * format. Note that the name of
				      * the output file must be
				      * specified through the
				      * TecplotFlags interface.
				      */
    void write_tecplot_binary (std::ostream &out) const;

				     /**
				      * Obtain data through
				      * get_patches() and write it to
				      * <tt>out</tt> in UCD format for
				      * @ref SoftwareAVS. See
				      * DataOutBase::write_ucd.
				      */
    void write_ucd (std::ostream &out) const;

    				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in Vtk format. See
				      * DataOutBase::write_vtk.
				      */
    void write_vtk (std::ostream &out) const;

				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in Vtu (VTK's XML) format. See
				      * DataOutBase::write_vtu.
				      *
				      * Some visualization programs, such as
				      * ParaView, can read several separate
				      * VTU files to parallelize
				      * visualization. In that case, you need
				      * a <code>.pvtu</code> file that
				      * describes which VTU files form a
				      * group. The
				      * DataOutInterface::write_pvtu_record()
				      * function can generate such a master
				      * record.
				      */
    void write_vtu (std::ostream &out) const;

				     /**
				      * Some visualization programs, such as
				      * ParaView, can read several separate
				      * VTU files to parallelize
				      * visualization. In that case, you need
				      * a <code>.pvtu</code> file that
				      * describes which VTU files (written,
				      * for example, through the write_vtu()
				      * function) form a group. The current
				      * function can generate such a master
				      * record.
				      *
				      * The file so written contains a list of
				      * (scalar or vector) fields whose values
				      * are described by the individual files
				      * that comprise the set of parallel VTU
				      * files along with the names of these
				      * files. This function gets the names
				      * and types of fields through the
				      * get_patches() function of this class
				      * like all the other write_xxx()
				      * functions. The second argument to this
				      * function specifies the names of the
				      * files that form the parallel set.
				      *
				      * @note See DataOutBase::write_vtu for
				      * writing each piece. Also note that
				      * only one parallel process needs to
				      * call the current function, listing the
				      * names of the files written by all
				      * parallel processes.
				      */
    void write_pvtu_record (std::ostream &out,
			    const std::vector<std::string> &piece_names) const;


    				     /**
				      * Obtain data through get_patches()
				      * and write it to <tt>out</tt>
				      * in deal.II intermediate
				      * format. See
				      * DataOutBase::write_deal_II_intermediate.
				      *
				      * Note that the intermediate
				      * format is what its name
				      * suggests: a direct
				      * representation of internal
				      * data. It isn't standardized
				      * and will change whenever we
				      * change our internal
				      * representation. You can only
				      * expect to process files
				      * written in this format using
				      * the same version of deal.II
				      * that was used for writing.
				      */
    void write_deal_II_intermediate (std::ostream &out) const;

				     /**
				      * Write data and grid to <tt>out</tt>
				      * according to the given data
				      * format. This function simply
				      * calls the appropriate
				      * <tt>write_*</tt> function. If no
				      * output format is requested,
				      * the <tt>default_format</tt> is
				      * written.
				      *
				      * An error occurs if no format
				      * is provided and the default
				      * format is <tt>default_format</tt>.
				      */
    void write (std::ostream       &out,
		const OutputFormat  output_format = default_format) const;

				     /**
				      * Set the default format. The
				      * value set here is used
				      * anytime, output for format
				      * <tt>default_format</tt> is
				      * requested.
				      */
    void set_default_format (const OutputFormat default_format);

				     /**
				      * Set the flags to be used for
				      * output in OpenDX format.
				      */
    void set_flags (const DXFlags &dx_flags);

				     /**
				      * Set the flags to be used for
				      * output in UCD format.
				      */
    void set_flags (const UcdFlags &ucd_flags);

    				     /**
				      * Set the flags to be used for
				      * output in GNUPLOT format.
				      */
    void set_flags (const GnuplotFlags &gnuplot_flags);

    				     /**
				      * Set the flags to be used for
				      * output in POVRAY format.
				      */
    void set_flags (const PovrayFlags &povray_flags);

    				     /**
				      * Set the flags to be used for
				      * output in EPS output.
				      */
    void set_flags (const EpsFlags &eps_flags);

    				     /**
				      * Set the flags to be used for
				      * output in GMV format.
				      */
    void set_flags (const GmvFlags &gmv_flags);

    				     /**
				      * Set the flags to be used for
				      * output in Tecplot format.
				      */
    void set_flags (const TecplotFlags &tecplot_flags);

    				     /**
				      * Set the flags to be used for
				      * output in VTK format.
				      */
    void set_flags (const VtkFlags &vtk_flags);

    				     /**
				      * Set the flags to be used for output in
				      * deal.II intermediate format.
				      */
    void set_flags (const Deal_II_IntermediateFlags &deal_II_intermediate_flags);

				     /**
				      * A function that returns the same
				      * string as the respective function in
				      * the base class does; the only
				      * exception being that if the parameter
				      * is omitted, then the value for the
				      * present default format is returned.
				      */
    std::string
    default_suffix (const OutputFormat output_format = default_format) const;

				     /**
				      * Declare parameters for all
				      * output formats by declaring
				      * subsections within the
				      * parameter file for each output
				      * format and call the respective
				      * <tt>declare_parameters</tt>
				      * functions of the flag classes
				      * for each output format.
				      *
				      * Some of the declared
				      * subsections may not contain
				      * entries, if the respective
				      * format does not export any
				      * flags.
				      *
				      * Note that the top-level
				      * parameters denoting the number
				      * of subdivisions per patch and
				      * the output format are not
				      * declared, since they are only
				      * passed to virtual functions
				      * and are not stored inside
				      * objects of this type. You have
				      * to declare them yourself.
				      */
    static void declare_parameters (ParameterHandler &prm);

				     /**
				      * Read the parameters declared
				      * in <tt>declare_parameters</tt> and
				      * set the flags for the output
				      * formats accordingly.
				      *
				      * The flags thus obtained
				      * overwrite all previous
				      * contents of the flag objects
				      * as default-constructed or set
				      * by the set_flags() function.
				      */
    void parse_parameters (ParameterHandler &prm);

				     /**
				      * Determine an estimate for
				      * the memory consumption (in
				      * bytes) of this
				      * object. Since sometimes
				      * the size of objects can
				      * not be determined exactly
				      * (for example: what is the
				      * memory consumption of an
				      * STL <tt>std::map</tt> type with a
				      * certain number of
				      * elements?), this is only
				      * an estimate. however often
				      * quite close to the true
				      * value.
				      */
    unsigned int memory_consumption () const;

  protected:
				     /**
				      * This is the abstract function
				      * through which derived classes
				      * propagate preprocessed data in
				      * the form of Patch
				      * structures (declared in the
				      * base class DataOutBase) to
				      * the actual output
				      * function. You need to overload
				      * this function to allow the
				      * output functions to know what
				      * they shall print.
				      */
    virtual
    const std::vector<typename DataOutBase::Patch<dim,spacedim> > &
    get_patches () const = 0;

				     /**
				      * Abstract virtual function
				      * through which the names of
				      * data sets are obtained by the
				      * output functions of the base
				      * class.
				      */
    virtual
    std::vector<std::string>
    get_dataset_names () const = 0;

				     /**
				      * This functions returns
				      * information about how the
				      * individual components of
				      * output files that consist of
				      * more than one data set are to
				      * be interpreted.
				      *
				      * It returns a list of index
				      * pairs and corresponding name
				      * indicating which components of
				      * the output are to be
				      * considered vector-valued
				      * rather than just a collection
				      * of scalar data. The index
				      * pairs are inclusive; for
				      * example, if we have a Stokes
				      * problem in 2d with components
				      * (u,v,p), then the
				      * corresponding vector data
				      * range should be (0,1), and the
				      * returned list would consist of
				      * only a single element with a
				      * tuple such as (0,1,"velocity").
				      *
				      * Since some of the derived
				      * classes do not know about
				      * vector data, this function has
				      * a default implementation that
				      * simply returns an empty
				      * string, meaning that all data
				      * is to be considered a
				      * collection of scalar fields.
				      */
    virtual
    std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> >
    get_vector_data_ranges () const;

				     /**
				      * The default number of
				      * subdivisions for patches. This
				      * is filled by parse_parameters()
				      * and should be obeyed by
				      * build_patches() in derived
				      * classes.
				      */
    unsigned int default_subdivisions;

  private:
				     /**
				      * Standard output format.  Use
				      * this format, if output format
				      * default_format is
				      * requested. It can be changed
				      * by the <tt>set_format</tt> function
				      * or in a parameter file.
				      */
    OutputFormat default_fmt;

				     /**
				      * Flags to be used upon output
				      * of OpenDX data. Can be changed by
				      * using the <tt>set_flags</tt>
				      * function.
				      */
    DXFlags     dx_flags;

				     /**
				      * Flags to be used upon output
				      * of UCD data. Can be changed by
				      * using the <tt>set_flags</tt>
				      * function.
				      */
    UcdFlags     ucd_flags;

				     /**
				      * Flags to be used upon output
				      * of GNUPLOT data. Can be
				      * changed by using the
				      * <tt>set_flags</tt> function.
				      */
    GnuplotFlags gnuplot_flags;

    				     /**
				      * Flags to be used upon output
				      * of POVRAY data. Can be changed
				      * by using the <tt>set_flags</tt>
				      * function.
				      */
    PovrayFlags povray_flags;

				     /**
				      * Flags to be used upon output
				      * of EPS data in one space
				      * dimension. Can be changed by
				      * using the <tt>set_flags</tt>
				      * function.
				      */
    EpsFlags     eps_flags;

				     /**
				      * Flags to be used upon output
				      * of gmv data in one space
				      * dimension. Can be changed by
				      * using the <tt>set_flags</tt>
				      * function.
				      */
    GmvFlags     gmv_flags;

				     /**
				      * Flags to be used upon output
				      * of Tecplot data in one space
				      * dimension. Can be changed by
				      * using the <tt>set_flags</tt>
				      * function.
				      */
    TecplotFlags tecplot_flags;

				     /**
				      * Flags to be used upon output
				      * of vtk data in one space
				      * dimension. Can be changed by
				      * using the <tt>set_flags</tt>
				      * function.
				      */
    VtkFlags     vtk_flags;

				     /**
				      * Flags to be used upon output of
				      * deal.II intermediate data in one space
				      * dimension. Can be changed by using the
				      * <tt>set_flags</tt> function.
				      */
    Deal_II_IntermediateFlags     deal_II_intermediate_flags;
};



/**
 * A class that is used to read data written in deal.II intermediate
 * format back in, so that it can be written out in any of the other
 * supported graphics formats. This class has two main purposes:
 *
 * The first use of this class is so that application programs can
 * defer the decision of which graphics format to use until after the
 * program has been run. The data is written in intermediate format
 * into a file, and later on it can then be converted into any
 * graphics format you wish. This may be useful, for example, if you
 * want to convert it to gnuplot format to get a quick glimpse and
 * later on want to convert it to OpenDX format as well to get a high
 * quality version of the data. The present class allows to read this
 * intermediate format back into the program, and allows it to be
 * written in any other supported format using the relevant functions
 * of the base class.
 *
 * The second use is mostly useful in parallel programs: rather than
 * having one central process generate the graphical output for the
 * entire program, one can let each process generate the graphical
 * data for the cells it owns, and write it into a separate file in
 * intermediate format. Later on, all these intermediate files can
 * then be read back in and merged together, a process that is fast
 * compared to generating the data in the first place. The use of the
 * intermediate format is mostly because it allows separate files to
 * be merged, while this is almost impossible once the data has been
 * written out in any of the supported established graphics formats.
 *
 * This second use scenario is explained in some detail in the step-18
 * example program.
 *
 * Both these applications are implemented in the step-19 example program.
 * There, a slight complication is also explained: in order to read data back
 * into this object, you have to know the template parameters for the space
 * dimension which were used when writing the data. If this knowledge is
 * available at compile time, then this is no problem. However, if it is not
 * (such as in a simple format converter), then it needs to be figured out at
 * run time, even though the compiler already needs it at compile time. A way
 * around using the DataOutBase::determine_intermediate_format_dimensions()
 * function is explained in step-19.
 *
 * Note that the intermediate format is what its name suggests: a
 * direct representation of internal data. It isn't standardized and
 * will change whenever we change our internal representation. You can
 * only expect to process files written in this format using the same
 * version of deal.II that was used for writing.
 *
 * @ingroup input output
 * @author Wolfgang Bangerth, 2005
 */
template <int dim, int spacedim=dim>
class DataOutReader : public DataOutInterface<dim,spacedim>
{
  public:
				     /**
				      * Read a sequence of patches as
				      * written previously by
				      * <tt>DataOutBase::write_deal_II_intermediate</tt>
				      * and store them in the present
				      * object. This overwrites any
				      * previous content.
				      */
    void read (std::istream &in);

                                     /**
                                      * This function can be used to
                                      * merge the patches read by the
                                      * other object into the patches
                                      * that this present object
                                      * stores. This is sometimes
                                      * handy if one has, for example,
                                      * a domain decomposition
                                      * algorithm where each block is
                                      * represented by a DoFHandler of
                                      * its own, but one wants to
                                      * output the solution on all the
                                      * blocks at the same
                                      * time. Alternatively, it may
                                      * also be used for parallel
                                      * programs, where each process
                                      * only generates output for its
                                      * share of the cells, even if
                                      * all processes can see all
                                      * cells.
                                      *
                                      * For this to work, the input
                                      * files for the present object
                                      * and the given argument need to
                                      * have the same number of output
                                      * vectors, and they need to use
                                      * the same number of
                                      * subdivisions per patch. The
                                      * output will probably look
                                      * rather funny if patches in
                                      * both objects overlap in space.
                                      *
                                      * If you call read() for this
                                      * object after merging in
                                      * patches, the previous state is
                                      * overwritten, and the merged-in
                                      * patches are lost.
				      *
                                      * This function will fail if
                                      * either this or the other
                                      * object did not yet set up any
                                      * patches.
				      *
				      * The use of this function is
				      * demonstrated in step-19.
                                      */
    void merge (const DataOutReader<dim,spacedim> &other);

                                     /**
                                      * Exception
                                      */
    DeclException0 (ExcNoPatches);
                                     /**
                                      * Exception
                                      */
    DeclException0 (ExcIncompatibleDatasetNames);
                                     /**
                                      * Exception
                                      */
    DeclException0 (ExcIncompatiblePatchLists);
				     /**
				      * Exception
				      */
    DeclException4 (ExcIncompatibleDimensions,
		    int, int, int, int,
		    << "Either the dimensions <" << arg1 << "> and <"
		    << arg2 << "> or the space dimensions <"
		    << arg3 << "> and <" << arg4
		    << "> do not match!");

  protected:
				     /**
				      * This is the function
				      * through which this class
				      * propagates preprocessed data in
				      * the form of Patch
				      * structures (declared in the
				      * base class DataOutBase) to
				      * the actual output
				      * function.
				      *
				      * It returns the patches as read
				      * the last time a stream was
				      * given to the read() function.
				      */
    virtual const std::vector<typename dealii::DataOutBase::Patch<dim,spacedim> > &
    get_patches () const;

				     /**
				      * Abstract virtual function
				      * through which the names of
				      * data sets are obtained by the
				      * output functions of the base
				      * class.
				      *
				      * Return the names of the
				      * variables as read the last
				      * time we read a file.
				      */
    virtual std::vector<std::string> get_dataset_names () const;

				     /**
				      * This functions returns
				      * information about how the
				      * individual components of
				      * output files that consist of
				      * more than one data set are to
				      * be interpreted.
				      *
				      * It returns a list of index
				      * pairs and corresponding name
				      * indicating which components of
				      * the output are to be
				      * considered vector-valued
				      * rather than just a collection
				      * of scalar data. The index
				      * pairs are inclusive; for
				      * example, if we have a Stokes
				      * problem in 2d with components
				      * (u,v,p), then the
				      * corresponding vector data
				      * range should be (0,1), and the
				      * returned list would consist of
				      * only a single element with a
				      * tuple such as (0,1,"velocity").
				      *
				      * Since some of the derived
				      * classes do not know about
				      * vector data, this function has
				      * a default implementation that
				      * simply returns an empty
				      * string, meaning that all data
				      * is to be considered a
				      * collection of scalar fields.
				      */
    virtual
    std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> >
    get_vector_data_ranges () const;

  private:
				     /**
				      * Arrays holding the set of
				      * patches as well as the names
				      * of output variables, all of
				      * which we read from an input
				      * stream.
				      */
    std::vector<typename dealii::DataOutBase::Patch<dim,spacedim> > patches;
    std::vector<std::string> dataset_names;

				     /**
				      * Information about whether
				      * certain components of the
				      * output field are to be
				      * considered vectors.
				      */
    std::vector<std_cxx1x::tuple<unsigned int, unsigned int, std::string> >
    vector_data_ranges;
};




/* -------------------- inline functions ------------------- */

inline
bool
DataOutBase::EpsFlags::RgbValues::is_grey () const
{
  return (red == green) && (red == blue);
}


/* -------------------- template functions ------------------- */

/**
 * Output operator for an object of type
 * <tt>DataOutBase::Patch</tt>. This operator dumps the intermediate
 * graphics format represented by the patch data structure. It may
 * later be converted into regular formats for a number of graphics
 * programs.
 *
 * @author Wolfgang Bangerth, 2005
 */
template <int dim, int spacedim>
std::ostream &
operator << (std::ostream                           &out,
	     const DataOutBase::Patch<dim,spacedim> &patch);



/**
 * Input operator for an object of type
 * <tt>DataOutBase::Patch</tt>. This operator reads the intermediate
 * graphics format represented by the patch data structure, using the
 * format in which it was written using the operator<<.
 *
 * @author Wolfgang Bangerth, 2005
 */
template <int dim, int spacedim>
std::istream &
operator >> (std::istream                     &in,
	     DataOutBase::Patch<dim,spacedim> &patch);



DEAL_II_NAMESPACE_CLOSE

#endif