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
|
/*=========================================================================
Program: Visualization Toolkit
Module: vtkMPASReader.h
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
/*=========================================================================
Copyright (c) 2002-2005 Los Alamos National Laboratory
This software and ancillary information known as vtk_ext (and herein
called "SOFTWARE") is made available under the terms described below.
The SOFTWARE has been approved for release with associated LA_CC
Number 99-44, granted by Los Alamos National Laboratory in July 1999.
Unless otherwise indicated, this SOFTWARE has been authored by an
employee or employees of the University of California, operator of the
Los Alamos National Laboratory under Contract No. W-7405-ENG-36 with
the United States Department of Energy.
The United States Government has rights to use, reproduce, and
distribute this SOFTWARE. The public may copy, distribute, prepare
derivative works and publicly display this SOFTWARE without charge,
provided that this Notice and any statement of authorship are
reproduced on all copies.
Neither the U. S. Government, the University of California, nor the
Advanced Computing Laboratory makes any warranty, either express or
implied, nor assumes any liability or responsibility for the use of
this SOFTWARE.
If SOFTWARE is modified to produce derivative works, such modified
SOFTWARE should be clearly marked, so as not to confuse it with the
version available from Los Alamos National Laboratory.
=========================================================================*/
// Christine Ahrens (cahrens@lanl.gov)
// Version 1.3
/*=========================================================================
NOTES
When using this reader, it is important that you remember to do the following:
1. When changing a selected variable, remember to select it also in the drop
down box to "color by". It doesn't color by that variable automatically.
2. When selecting multilayer sphere view, make layer thickness around
100,000.
3. When selecting multilayer lat/lon view, make layer thickness around 10.
4. Always click the -Z orientation after making a switch from lat/lon to
sphere, from single to multilayer or changing thickness.
5. Be conservative on the number of changes you make before hitting Apply,
since there may be bugs in this reader. Just make one change and then hit
Apply.
=========================================================================*/
#include "vtkMPASReader.h"
#include "vtkCallbackCommand.h"
#include "vtkCellData.h"
#include "vtkCellType.h"
#include "vtkCellArray.h"
#include "vtkDataArraySelection.h"
#include "vtkDataObject.h"
#include "vtkDoubleArray.h"
#include "vtkInformation.h"
#include "vtkInformationDoubleVectorKey.h"
#include "vtkInformationVector.h"
#include "vtkIntArray.h"
#include "vtkMath.h"
#include "vtkNew.h"
#include "vtkObjectFactory.h"
#include "vtkPointData.h"
#include "vtkSmartPointer.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include "vtkStringArray.h"
#include "vtkToolkits.h"
#include "vtkUnstructuredGrid.h"
#include "vtk_netcdfcpp.h"
#include <algorithm>
#include <cfloat>
#include <cmath>
#include <cstdarg>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <map>
#include <set>
#include <sstream>
#include <string>
#include <vector>
using namespace std;
// In VTK commit 64cb89e3e6ae08f440eb6d4cbfb308c41ab7d258, a lot of signatures
// in the netcdf library changed size values from 'long' to 'size_t'. However,
// upstream versions of netcdf are using long. The following typedef resolves
// this issue.
// We cannot just check for VTK_USE_SYSTEM_NETCDF because ParaView's superbuild
// installs a "system" netcdf with the same modifications...
#ifdef VTK_NETCDF_USE_SIZE_T
typedef size_t nc_size_t;
#else
typedef long nc_size_t;
#endif
// Restricted to the supported NcType-convertible types.
#define vtkNcTemplateMacro(call) \
vtkTemplateMacroCase(VTK_DOUBLE, double, call); \
vtkTemplateMacroCase(VTK_FLOAT, float, call); \
vtkTemplateMacroCase(VTK_INT, int, call); \
vtkTemplateMacroCase(VTK_SHORT, short, call); \
vtkTemplateMacroCase(VTK_CHAR, char, call); \
vtkTemplateMacroCase(VTK_SIGNED_CHAR, signed char, call) /* ncbyte */
#define vtkNcDispatch(type, call) \
switch (type) \
{ \
vtkNcTemplateMacro(call); \
default: \
vtkErrorMacro(<<"Unsupported data type: " << (type)); \
abort(); \
}
namespace {
std::string dimensionedArrayName(NcVar *var)
{
std::ostringstream out;
out << var->name() << "(";
for (int dim = 0; dim < var->num_dims(); ++dim)
{
if (dim != 0)
{
out << ", ";
}
out << var->get_dim(dim)->name();
}
out << ")";
return out.str();
}
struct DimMetaData
{
long curIdx;
long dimSize;
};
} // end anon namespace
//----------------------------------------------------------------------------
// Internal class to avoid name pollution
//----------------------------------------------------------------------------
class vtkMPASReader::Internal
{
public:
// variableIndex --> vtkDataArray
typedef std::map<int, vtkSmartPointer<vtkDataArray> > ArrayMap;
typedef std::map<std::string, DimMetaData> DimMetaDataMap;
Internal() : ncFile(NULL) {}
~Internal() { delete ncFile; }
NcFile* ncFile;
std::vector<NcVar*> pointVars;
std::vector<NcVar*> cellVars;
ArrayMap pointArrays;
ArrayMap cellArrays;
// Returns true if the dimension name is not nCells, nVertices, or Time.
bool isExtraDim(const std::string &name);
// Indices at which arbitrary trailing dimensions are fixed:
DimMetaDataMap dimMetaDataMap;
vtkTimeStamp dimMetaDataTime;
// Set of dimensions currently used by the selected arrays:
vtkNew<vtkStringArray> extraDims;
vtkTimeStamp extraDimTime;
};
bool vtkMPASReader::Internal::isExtraDim(const string &name)
{
return name != "nCells" && name != "nVertices" && name != "Time";
}
//----------------------------------------------------------------------------
// Macro to check malloc didn't return an error
//----------------------------------------------------------------------------
#define CHECK_MALLOC(ptr) \
if (ptr == NULL) \
{ \
vtkErrorMacro( << "malloc failed!" << endl); \
return(0); \
}
//----------------------------------------------------------------------------
// Macro to check if the named NetCDF dimension exists
//----------------------------------------------------------------------------
#define CHECK_DIM(ncFile, name) \
if (!isNcDim(ncFile, name)) \
{ \
vtkErrorMacro( << "Cannot find dimension: " << name << endl); \
return 0; \
}
//----------------------------------------------------------------------------
// Macro to check if the named NetCDF variable exists
//----------------------------------------------------------------------------
#define CHECK_VAR(ncFile, name) \
if (!isNcVar(ncFile, name)) \
{ \
vtkErrorMacro( << "Cannot find variable: " << name << endl); \
return 0; \
}
//----------------------------------------------------------------------------
// Function to check if there is a NetCDF variable by that name
//-----------------------------------------------------------------------------
static bool isNcVar(NcFile *ncFile, NcToken name)
{
int num_vars = ncFile->num_vars();
for (int i = 0; i < num_vars; i++)
{
NcVar* ncVar = ncFile->get_var(i);
if ((strcmp(ncVar->name(), name)) == 0)
{
// we have a match, so return
return true;
}
}
return false;
}
//----------------------------------------------------------------------------
// Check if there is a NetCDF dimension by that name
//----------------------------------------------------------------------------
static bool isNcDim(NcFile *ncFile, NcToken name)
{
int num_dims = ncFile->num_dims();
//cerr << "looking for: " << name << endl;
for (int i = 0; i < num_dims; i++)
{
NcDim* ncDim = ncFile->get_dim(i);
//cerr << "checking " << ncDim->name() << endl;
const char* ncDimName = ncDim->name();
if ((ncDimName && strcmp(ncDimName, name)) == 0)
{
// we have a match, so return
return true;
}
}
return false;
}
//----------------------------------------------------------------------------
// Check if there is a NetCDF attribute by that name
//----------------------------------------------------------------------------
static bool isNcAtt(NcFile *ncFile, NcToken name)
{
int num_atts = ncFile->num_atts();
for (int i = 0; i < num_atts; i++)
{
NcAtt *ncAtt = ncFile->get_att(i);
if ((strcmp(ncAtt->name(), name)) == 0)
{
return true;
}
}
return false;
}
//-----------------------------------------------------------------------------
// Function to convert cartesian coordinates to spherical, for use in
// computing points in different layers of multilayer spherical view
//----------------------------------------------------------------------------
static int CartesianToSpherical(double x, double y, double z, double* rho,
double* phi, double* theta)
{
double trho, ttheta, tphi;
trho = sqrt( (x*x) + (y*y) + (z*z));
ttheta = atan2(y, x);
tphi = acos(z/(trho));
if (vtkMath::IsNan(trho) || vtkMath::IsNan(ttheta) || vtkMath::IsNan(tphi))
{
return -1;
}
*rho = trho;
*theta = ttheta;
*phi = tphi;
return 0;
}
//----------------------------------------------------------------------------
// Function to convert spherical coordinates to cartesian, for use in
// computing points in different layers of multilayer spherical view
//----------------------------------------------------------------------------
static int SphericalToCartesian(double rho, double phi, double theta, double* x,
double* y, double* z)
{
double tx, ty, tz;
tx = rho* sin(phi) * cos(theta);
ty = rho* sin(phi) * sin(theta);
tz = rho* cos(phi);
if (vtkMath::IsNan(tx) || vtkMath::IsNan(ty) || vtkMath::IsNan(tz))
{
return -1;
}
*x = tx;
*y = ty;
*z = tz;
return 0;
}
vtkStandardNewMacro(vtkMPASReader);
//----------------------------------------------------------------------------
// Constructor for vtkMPASReader
//----------------------------------------------------------------------------
vtkMPASReader::vtkMPASReader()
{
this->Internals = new vtkMPASReader::Internal;
// Debugging
//this->DebugOn();
vtkDebugMacro(<< "Starting to create vtkMPASReader..." << endl);
this->SetNumberOfInputPorts(0);
this->SetNumberOfOutputPorts(1);
this->SetDefaults();
// Setup selection callback to modify this object when array selection changes
this->PointDataArraySelection = vtkDataArraySelection::New();
this->CellDataArraySelection = vtkDataArraySelection::New();
this->SelectionObserver = vtkCallbackCommand::New();
this->SelectionObserver->SetCallback(&vtkMPASReader::SelectionCallback);
this->SelectionObserver->SetClientData(this);
this->CellDataArraySelection->AddObserver(vtkCommand::ModifiedEvent,
this->SelectionObserver);
this->PointDataArraySelection->AddObserver(vtkCommand::ModifiedEvent,
this->SelectionObserver);
vtkDebugMacro(<< "Created vtkMPASReader" << endl);
}
//----------------------------------------------------------------------------
// Destroys data stored for variables, points, and cells, but
// doesn't destroy the list of variables or toplevel cell/pointVarDataArray.
//----------------------------------------------------------------------------
void vtkMPASReader::DestroyData()
{
vtkDebugMacro(<<"DestroyData...");
this->Internals->cellArrays.clear();
this->Internals->pointArrays.clear();
free(this->CellMap);
this->CellMap = NULL;
free(this->PointMap);
this->PointMap = NULL;
free(this->MaximumLevelPoint);
this->MaximumLevelPoint = NULL;
}
//----------------------------------------------------------------------------
// Destructor for MPAS Reader
//----------------------------------------------------------------------------
vtkMPASReader::~vtkMPASReader()
{
vtkDebugMacro(<< "Destructing vtkMPASReader..." << endl);
this->SetFileName(NULL);
delete this->Internals->ncFile;
this->Internals->ncFile = NULL;
this->DestroyData();
vtkDebugMacro(<< "Destructing other stuff..." << endl);
if (this->PointDataArraySelection)
{
this->PointDataArraySelection->Delete();
this->PointDataArraySelection = NULL;
}
if (this->CellDataArraySelection)
{
this->CellDataArraySelection->Delete();
this->CellDataArraySelection = NULL;
}
if (this->SelectionObserver)
{
this->SelectionObserver->Delete();
this->SelectionObserver = NULL;
}
delete this->Internals;
vtkDebugMacro(<< "Destructed vtkMPASReader" << endl);
}
//----------------------------------------------------------------------------
void vtkMPASReader::ReleaseNcData()
{
this->Internals->pointVars.clear();
this->Internals->pointArrays.clear();
this->Internals->cellVars.clear();
this->Internals->cellArrays.clear();
this->PointDataArraySelection->RemoveAllArrays();
this->CellDataArraySelection->RemoveAllArrays();
this->UpdateDimensions(true); // Reset extra dimension list.
free(this->PointX);
this->PointX = NULL;
free(this->PointY);
this->PointY = NULL;
free(this->PointZ);
this->PointZ = NULL;
free(this->OrigConnections);
this->OrigConnections = NULL;
free(this->ModConnections);
this->ModConnections = NULL;
free(this->CellMap);
this->CellMap = NULL;
free(this->PointMap);
this->PointMap = NULL;
free(this->MaximumLevelPoint);
this->MaximumLevelPoint = NULL;
delete this->Internals->ncFile;
this->Internals->ncFile = NULL;
}
//----------------------------------------------------------------------------
// Verify that the file exists, get dimension sizes and variables
//----------------------------------------------------------------------------
int vtkMPASReader::RequestInformation(
vtkInformation *reqInfo,
vtkInformationVector **inVector,
vtkInformationVector *outVector)
{
vtkDebugMacro(<< "In vtkMPASReader::RequestInformation" << endl);
this->ReleaseNcData();
if (!this->Superclass::RequestInformation(reqInfo, inVector, outVector))
{
return 0;
}
// Verify that file exists
if ( !this->FileName )
{
vtkErrorMacro("No filename specified");
return 0;
}
// Get ParaView information pointer
vtkInformation* outInfo = outVector->GetInformationObject(0);
this->Internals->ncFile = new NcFile(this->FileName);
if (!this->Internals->ncFile->is_valid())
{
vtkErrorMacro( << "Couldn't open file: " << this->FileName << endl);
this->ReleaseNcData();
return 0;
}
if (!this->GetNcDims())
{
this->ReleaseNcData();
return 0;
}
if (!this->GetNcAtts())
{
this->ReleaseNcData();
return 0;
}
if (!this->CheckParams())
{
this->ReleaseNcData();
return 0;
}
if (!this->BuildVarArrays())
{
this->ReleaseNcData();
return 0;
}
// Collect temporal information
// At this time, MPAS doesn't have fine-grained time value, just
// the number of the step, so that is what I store here for TimeSteps.
if (this->NumberOfTimeSteps > 0)
{
// Tell the pipeline what steps are available
std::vector<double> timeSteps;
for (int i = 0; i < this->NumberOfTimeSteps; ++i)
{
timeSteps.push_back(static_cast<double>(i));
}
outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),
&timeSteps[0], static_cast<int>(timeSteps.size()));
double tRange[2];
tRange[0] = 0.;
tRange[1] = static_cast<double>(this->NumberOfTimeSteps - 1);
outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(), tRange, 2);
}
else
{
outInfo->Remove(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
outInfo->Remove(vtkStreamingDemandDrivenPipeline::TIME_RANGE());
}
return 1;
}
//----------------------------------------------------------------------------
// Data is read into a vtkUnstructuredGrid
//----------------------------------------------------------------------------
int vtkMPASReader::RequestData(vtkInformation *vtkNotUsed(reqInfo),
vtkInformationVector **vtkNotUsed(inVector),
vtkInformationVector *outVector)
{
vtkDebugMacro(<< "In vtkMPASReader::RequestData" << endl);
// get the info object
vtkInformation *outInfo = outVector->GetInformationObject(0);
// Output will be an ImageData
vtkUnstructuredGrid *output = vtkUnstructuredGrid::SafeDownCast(
outInfo->Get(vtkDataObject::DATA_OBJECT()));
this->DestroyData();
if (!this->ReadAndOutputGrid())
{
this->DestroyData();
return 0;
}
// Collect the time step requested
this->DTime = 0;
if (outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
{
this->DTime =
outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
}
output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(), this->DTime);
// Examine each variable to see if it is selected
int numPointVars = static_cast<int>(this->Internals->pointVars.size());
for (int var = 0; var < numPointVars; var++)
{
// Is this variable requested
if (this->PointDataArraySelection->GetArraySetting(var))
{
vtkDataArray *array = this->LoadPointVarData(var);
if (!array)
{
vtkWarningMacro(<< "Error loading point variable '"
<< this->Internals->pointVars[var]->name() << "'.");
continue;
}
output->GetPointData()->AddArray(array);
}
}
int numCellVars = static_cast<int>(this->Internals->cellVars.size());
for (int var = 0; var < numCellVars; var++)
{
if (this->CellDataArraySelection->GetArraySetting(var))
{
vtkDataArray *array = this->LoadCellVarData(var);
if (!array)
{
vtkWarningMacro(<< "Error loading point variable '"
<< this->Internals->pointVars[var]->name() << "'.");
continue;
}
output->GetCellData()->AddArray(array);
}
}
this->LoadTimeFieldData(output);
vtkDebugMacro( << "Returning from RequestData" << endl);
return 1;
}
//----------------------------------------------------------------------------
// Set defaults for various parameters and initialize some variables
//----------------------------------------------------------------------------
void vtkMPASReader::SetDefaults() {
// put in defaults
this->VerticalDimension = "nVertLevels";
this->VerticalLevelRange[0] = 0;
this->VerticalLevelRange[1] = 1;
this->LayerThicknessRange[0] = 0;
this->LayerThicknessRange[1] = 200000;
this->LayerThickness = 10000;
vtkDebugMacro
( << "SetDefaults: LayerThickness set to " << LayerThickness << endl);
this->CenterLonRange[0] = 0;
this->CenterLonRange[1] = 360;
this->CenterLon = 180;
this->Geometry = Spherical;
this->IsAtmosphere = false;
this->ProjectLatLon = false;
this->OnASphere = false;
this->ShowMultilayerView = false;
this->IsZeroCentered = false;
this->IncludeTopography = false;
this->DoBugFix = false;
this->CenterRad = CenterLon * vtkMath::Pi() / 180.0;
this->UseDimensionedArrayNames = false;
this->PointX = NULL;
this->PointY = NULL;
this->PointZ = NULL;
this->OrigConnections = NULL;
this->ModConnections = NULL;
this->CellMap = NULL;
this->PointMap = NULL;
this->MaximumLevelPoint = NULL;
this->FileName = NULL;
this->DTime = 0;
this->MaximumPoints = 0;
this->MaximumCells = 0;
}
//----------------------------------------------------------------------------
// Get dimensions of key NetCDF variables
//----------------------------------------------------------------------------
int vtkMPASReader::GetNcDims()
{
NcFile *pnf = this->Internals->ncFile;
CHECK_DIM(pnf, "nCells");
NcDim* nCells = pnf->get_dim("nCells");
this->NumberOfPoints = nCells->size();
this->PointOffset = 1;
CHECK_DIM(pnf, "nVertices");
NcDim* nVertices = pnf->get_dim("nVertices");
this->NumberOfCells = nVertices->size();
this->CellOffset = 0;
CHECK_DIM(pnf, "vertexDegree");
NcDim* vertexDegree = pnf->get_dim("vertexDegree");
this->PointsPerCell = vertexDegree->size();
CHECK_DIM(pnf, "Time");
NcDim* Time = pnf->get_dim("Time");
this->NumberOfTimeSteps = Time->size();
if (isNcDim(pnf, this->VerticalDimension.c_str()))
{
NcDim* nVertLevels = pnf->get_dim(this->VerticalDimension.c_str());
this->MaximumNVertLevels = nVertLevels->size();
}
else
{
this->MaximumNVertLevels = 0;
}
return 1;
}
//----------------------------------------------------------------------------
int vtkMPASReader::GetNcAtts()
{
NcFile *pnf = this->Internals->ncFile;
if (!isNcAtt(pnf, "on_a_sphere"))
{
vtkWarningMacro("Attribute 'on_a_sphere' missing in file " << this->FileName
<< ". Assuming \"YES\".");
this->OnASphere = true;
}
else
{
NcAtt *onASphere = pnf->get_att("on_a_sphere");
char *onASphereString = onASphere->as_string(0);
this->OnASphere = (strcmp(onASphereString, "YES") == 0);
delete [] onASphereString;
onASphereString = NULL;
}
return 1;
}
//----------------------------------------------------------------------------
// Check parameters are valid
//----------------------------------------------------------------------------
int vtkMPASReader::CheckParams()
{
if ((this->PointsPerCell != 3) && (this->PointsPerCell != 4))
{
vtkErrorMacro
("This code is only for hexagonal or quad primal grids" << endl);
return(0);
}
this->VerticalLevelRange[0] = 0;
this->VerticalLevelRange[1] = this->MaximumNVertLevels-1;
if (this->OnASphere)
{
if (this->ProjectLatLon)
{
this->Geometry = Projected;
}
else
{
this->Geometry = Spherical;
}
}
else
{
this->Geometry = Planar;
if (this->ProjectLatLon)
{
vtkWarningMacro("Ignoring ProjectLatLong -- Data is not on_a_sphere.");
}
}
return(1);
}
//----------------------------------------------------------------------------
// Get the NetCDF variables on cell or vertex
//----------------------------------------------------------------------------
int vtkMPASReader::GetNcVars (const char* cellDimName, const char* pointDimName)
{
this->Internals->pointArrays.clear();
this->Internals->pointVars.clear();
this->Internals->cellArrays.clear();
this->Internals->cellVars.clear();
NcFile* ncFile = this->Internals->ncFile;
int numVars = ncFile->num_vars();
for (int i = 0; i < numVars; i++)
{
NcVar* var = ncFile->get_var(i);
// Variables must have the following dimension specification:
// [Time, ] (nCells | nVertices), [arbitraryDim1, [arbitraryDim2, [...]]]
bool isPointData = false;
bool isCellData = false;
int numDims = var->num_dims();
std::vector<std::string> dimNames;
for (int dim = 0; dim < std::min(numDims, 2); ++dim)
{
dimNames.push_back(var->get_dim(dim)->name());
}
if (numDims < 1)
{
vtkWarningMacro(<<"Variable '" << var->name()
<< "' has invalid number of dimensions: " << numDims);
continue;
}
if (dimNames[0] == "Time" && dimNames.size() >= 2)
{
if (dimNames[1] == pointDimName)
{
isPointData = true;
}
else if (dimNames[1] == cellDimName)
{
isCellData = true;
}
}
else if (dimNames[0] == pointDimName)
{
isPointData = true;
}
else if (dimNames[0] == cellDimName)
{
isCellData = true;
}
// Add to cell or point var array
if (isCellData)
{
this->Internals->cellVars.push_back(var);
}
else if (isPointData)
{
this->Internals->pointVars.push_back(var);
}
}
return(1);
}
//----------------------------------------------------------------------------
// Build the selection Arrays for points and cells in the GUI.
//----------------------------------------------------------------------------
int vtkMPASReader::BuildVarArrays()
{
// figure out what variables to visualize -
if (!GetNcVars("nVertices", "nCells"))
{
return 0;
}
for (size_t v = 0; v < this->Internals->pointVars.size(); v++)
{
NcVar *var = this->Internals->pointVars[v];
string name = this->UseDimensionedArrayNames ? dimensionedArrayName(var)
: var->name();
this->PointDataArraySelection->EnableArray(name.c_str());
// Register the dimensions:
for (int d = 0; d < var->num_dims(); ++d)
{
this->InitializeDimension(var->get_dim(d));
}
vtkDebugMacro(<<"Adding point var: " << name);
}
for (size_t v = 0; v < this->Internals->cellVars.size(); v++)
{
NcVar *var = this->Internals->cellVars[v];
string name = this->UseDimensionedArrayNames ? dimensionedArrayName(var)
: var->name();
this->CellDataArraySelection->EnableArray(name.c_str());
// Register the dimensions:
for (int d = 0; d < var->num_dims(); ++d)
{
this->InitializeDimension(var->get_dim(d));
}
vtkDebugMacro(<<"Adding cell var: " << name);
}
return 1;
}
//----------------------------------------------------------------------------
// Read the data from the ncfile, allocate the geometry and create the
// vtk data structures for points and cells.
//----------------------------------------------------------------------------
int vtkMPASReader::ReadAndOutputGrid()
{
switch (this->Geometry)
{
case vtkMPASReader::Spherical:
if (!this->AllocSphericalGeometry())
{
return 0;
}
this->FixPoints();
break;
case vtkMPASReader::Projected:
if (!this->AllocProjectedGeometry())
{
return 0;
}
this->ShiftLonData();
this->FixPoints();
if (!this->EliminateXWrap())
{
return 0;
}
break;
case vtkMPASReader::Planar:
if (!this->AllocPlanarGeometry())
{
return 0;
}
this->FixPoints();
if (!this->EliminateXWrap())
{
return 0;
}
break;
default:
vtkErrorMacro("Invalid geometry type (" << this->Geometry << ").");
return 0;
}
this->OutputPoints();
this->OutputCells();
return 1;
}
//----------------------------------------------------------------------------
// Allocate into sphere view of dual geometry
//----------------------------------------------------------------------------
int vtkMPASReader::AllocSphericalGeometry()
{
NcFile* ncFile = this->Internals->ncFile;
CHECK_VAR(ncFile, "xCell");
this->PointX = (double*)malloc((this->NumberOfPoints+this->PointOffset) *
sizeof(double));
CHECK_MALLOC(this->PointX);
NcVar* xCellVar = ncFile->get_var("xCell");
if (!this->ValidateDimensions(xCellVar, false, 1, "nCells"))
{
return 0;
}
xCellVar->get(this->PointX + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointX[0] = 0.0;
CHECK_VAR(ncFile, "yCell");
this->PointY = (double*)malloc((this->NumberOfPoints+this->PointOffset) *
sizeof(double));
CHECK_MALLOC(this->PointY);
NcVar* yCellVar = ncFile->get_var("yCell");
if (!this->ValidateDimensions(yCellVar, false, 1, "nCells"))
{
return 0;
}
yCellVar->get(this->PointY + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointY[0] = 0.0;
CHECK_VAR(ncFile, "zCell");
this->PointZ = (double*)malloc((this->NumberOfPoints+this->PointOffset) *
sizeof(double));
CHECK_MALLOC(this->PointZ);
NcVar* zCellVar = ncFile->get_var("zCell");
if (!this->ValidateDimensions(zCellVar, false, 1, "nCells"))
{
return 0;
}
zCellVar->get(this->PointZ + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointZ[0] = 0.0;
CHECK_VAR(ncFile, "cellsOnVertex");
this->OrigConnections = (int *) malloc(this->NumberOfCells *
this->PointsPerCell * sizeof(int));
CHECK_MALLOC(this->OrigConnections);
NcVar *connectionsVar = ncFile->get_var("cellsOnVertex");
// TODO Spec says dims should be '3', 'nVertices', but my example files
// use nVertices, vertexDegree...
if (!this->ValidateDimensions(connectionsVar, false, 2,
"nVertices", "vertexDegree"))
{
return 0;
}
connectionsVar->get(this->OrigConnections, this->NumberOfCells,
this->PointsPerCell);
if (isNcVar(ncFile, "maxLevelCell"))
{
this->IncludeTopography = true;
this->MaximumLevelPoint = (int*)malloc((this->NumberOfPoints +
this->PointOffset) * sizeof(int));
CHECK_MALLOC(this->MaximumLevelPoint);
NcVar *maxLevelPointVar = ncFile->get_var("maxLevelCell");
if (!this->ValidateDimensions(maxLevelPointVar, false, 1, "nCells"))
{
return 0;
}
maxLevelPointVar->get(this->MaximumLevelPoint + this->PointOffset,
this->NumberOfPoints);
}
this->CurrentExtraPoint = this->NumberOfPoints + this->PointOffset;
this->CurrentExtraCell = this->NumberOfCells + this->CellOffset;
if (this->ShowMultilayerView)
{
this->MaximumCells = this->CurrentExtraCell*this->MaximumNVertLevels;
vtkDebugMacro
(<< "alloc sphere: multilayer: setting MaximumCells to " << this->MaximumCells);
this->MaximumPoints = this->CurrentExtraPoint*(this->MaximumNVertLevels+1);
vtkDebugMacro
(<< "alloc sphere: multilayer: setting MaximumPoints to " << this->MaximumPoints);
}
else
{
this->MaximumCells = this->CurrentExtraCell;
this->MaximumPoints = this->CurrentExtraPoint;
vtkDebugMacro
(<< "alloc sphere: singlelayer: setting MaximumPoints to " << this->MaximumPoints);
}
return 1;
}
//----------------------------------------------------------------------------
// Allocate the lat/lon projection of dual geometry.
//----------------------------------------------------------------------------
int vtkMPASReader::AllocProjectedGeometry()
{
NcFile* ncFile = this->Internals->ncFile;
const float BLOATFACTOR = .5;
this->ModNumPoints = (int)floor(this->NumberOfPoints*(1.0 + BLOATFACTOR));
this->ModNumCells = (int)floor(this->NumberOfCells*(1.0 + BLOATFACTOR))+1;
CHECK_VAR(ncFile, "lonCell");
this->PointX = (double*)malloc(this->ModNumPoints * sizeof(double));
CHECK_MALLOC(this->PointX);
NcVar* xCellVar = ncFile->get_var("lonCell");
if (!this->ValidateDimensions(xCellVar, false, 1, "nCells"))
{
return 0;
}
xCellVar->get(this->PointX + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointX[0] = 0.0;
CHECK_VAR(ncFile, "latCell");
this->PointY = (double*)malloc(this->ModNumPoints * sizeof(double));
CHECK_MALLOC(this->PointY);
NcVar* yCellVar = ncFile->get_var("latCell");
if (!this->ValidateDimensions(yCellVar, false, 1, "nCells"))
{
return 0;
}
yCellVar->get(this->PointY+this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointY[0] = 0.0;
CHECK_VAR(ncFile, "cellsOnVertex");
this->OrigConnections = (int *) malloc(this->NumberOfCells * this->PointsPerCell *
sizeof(int));
CHECK_MALLOC(this->OrigConnections);
NcVar *connectionsVar = ncFile->get_var("cellsOnVertex");
// TODO Spec says dims should be '3', 'nVertices', but my example files
// use nVertices, vertexDegree...
if (!this->ValidateDimensions(connectionsVar, false, 2,
"nVertices", "vertexDegree"))
{
return 0;
}
connectionsVar->get(this->OrigConnections, this->NumberOfCells,
this->PointsPerCell);
// create my own list to include modified origConnections (due to
// eliminating wraparound in the lat/lon projection) plus additional
// cells added when mirroring cells that had previously wrapped around
this->ModConnections = (int *) malloc(this->ModNumCells * this->PointsPerCell
* sizeof(int));
CHECK_MALLOC(this->ModConnections);
// allocate an array to map the extra points and cells to the original
// so that when obtaining data, we know where to get it
this->PointMap = (int*)malloc((int)floor(this->NumberOfPoints*BLOATFACTOR)
* sizeof(int));
CHECK_MALLOC(this->PointMap);
this->CellMap = (int*)malloc((int)floor(this->NumberOfCells*BLOATFACTOR)
* sizeof(int));
CHECK_MALLOC(this->CellMap);
if (isNcVar(ncFile, "maxLevelCell"))
{
this->IncludeTopography = true;
this->MaximumLevelPoint = (int*)malloc((this->NumberOfPoints + this->NumberOfPoints) * sizeof(int));
CHECK_MALLOC(this->MaximumLevelPoint);
NcVar *maxLevelPointVar = ncFile->get_var("maxLevelCell");
if (!this->ValidateDimensions(maxLevelPointVar, false, 1, "nCells"))
{
return 0;
}
maxLevelPointVar->get(this->MaximumLevelPoint + this->PointOffset,
this->NumberOfPoints);
}
this->CurrentExtraPoint = this->NumberOfPoints + this->PointOffset;
this->CurrentExtraCell = this->NumberOfCells + this->CellOffset;
if (ShowMultilayerView)
{
this->MaximumCells = this->CurrentExtraCell*this->MaximumNVertLevels;
this->MaximumPoints = this->CurrentExtraPoint*(this->MaximumNVertLevels+1);
vtkDebugMacro
(<< "alloc latlon: multilayer: setting this->MaximumPoints to " << this->MaximumPoints
<< endl);
}
else
{
this->MaximumCells = this->CurrentExtraCell;
this->MaximumPoints = this->CurrentExtraPoint;
vtkDebugMacro
(<< "alloc latlon: singlelayer: setting this->MaximumPoints to " << this->MaximumPoints
<< endl);
}
return 1;
}
int vtkMPASReader::AllocPlanarGeometry()
{
NcFile* ncFile = this->Internals->ncFile;
const float BLOATFACTOR = .5f;
this->ModNumPoints = (int)floor(this->NumberOfPoints*(1.0 + BLOATFACTOR));
this->ModNumCells = (int)floor(this->NumberOfCells*(1.0 + BLOATFACTOR))+1;
CHECK_VAR(ncFile, "xCell");
this->PointX = (double*)malloc(this->ModNumPoints * sizeof(double));
CHECK_MALLOC(this->PointX);
NcVar* xCellVar = ncFile->get_var("xCell");
if (!this->ValidateDimensions(xCellVar, false, 1, "nCells"))
{
return 0;
}
xCellVar->get(this->PointX + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointX[0] = 0.0;
CHECK_VAR(ncFile, "yCell");
this->PointY = (double*)malloc(this->ModNumPoints * sizeof(double));
CHECK_MALLOC(this->PointY);
NcVar* yCellVar = ncFile->get_var("yCell");
if (!this->ValidateDimensions(yCellVar, false, 1, "nCells"))
{
return 0;
}
yCellVar->get(this->PointY + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointY[0] = 0.0;
CHECK_VAR(ncFile, "zCell");
this->PointZ = (double*)malloc(this->ModNumPoints * sizeof(double));
CHECK_MALLOC(this->PointZ);
NcVar *zCellVar = ncFile->get_var("zCell");
if (!this->ValidateDimensions(zCellVar, false, 1, "nCells"))
{
return 0;
}
zCellVar->get(this->PointZ + this->PointOffset, this->NumberOfPoints);
// point 0 is 0.0
this->PointZ[0] = 0.0;
CHECK_VAR(ncFile, "cellsOnVertex");
this->OrigConnections = (int *) malloc(this->NumberOfCells * this->PointsPerCell *
sizeof(int));
CHECK_MALLOC(this->OrigConnections);
NcVar *connectionsVar = ncFile->get_var("cellsOnVertex");
// TODO Spec says dims should be '3', 'nVertices', but my example files
// use nVertices, vertexDegree...
if (!this->ValidateDimensions(connectionsVar, false, 2,
"nVertices", "vertexDegree"))
{
return 0;
}
connectionsVar->get(this->OrigConnections, this->NumberOfCells,
this->PointsPerCell);
// create list to include modified origConnections (due to eliminating
// wraparound) plus additional cells added when mirroring cells that had
// previously wrapped around
this->ModConnections = (int *) malloc(this->ModNumCells * this->PointsPerCell
* sizeof(int));
CHECK_MALLOC(this->ModConnections);
// allocate an array to map the extra points and cells to the original
// so that when obtaining data, we know where to get it
this->PointMap = (int*)malloc((int)floor(this->NumberOfPoints*BLOATFACTOR)
* sizeof(int));
CHECK_MALLOC(this->PointMap);
this->CellMap = (int*)malloc((int)floor(this->NumberOfCells*BLOATFACTOR)
* sizeof(int));
CHECK_MALLOC(this->CellMap);
if (isNcVar(ncFile, "maxLevelCell"))
{
this->IncludeTopography = true;
this->MaximumLevelPoint = (int*)malloc(2 * this->NumberOfPoints *
sizeof(int));
CHECK_MALLOC(this->MaximumLevelPoint);
NcVar *maxLevelPointVar = ncFile->get_var("maxLevelCell");
if (!this->ValidateDimensions(maxLevelPointVar, false, 1, "nCells"))
{
return 0;
}
maxLevelPointVar->get(this->MaximumLevelPoint + this->PointOffset,
this->NumberOfPoints);
}
this->CurrentExtraPoint = this->NumberOfPoints + this->PointOffset;
this->CurrentExtraCell = this->NumberOfCells + this->CellOffset;
if (this->ShowMultilayerView)
{
this->MaximumCells = this->CurrentExtraCell * this->MaximumNVertLevels;
this->MaximumPoints = this->CurrentExtraPoint *
(this->MaximumNVertLevels + 1);
}
else
{
this->MaximumCells = this->CurrentExtraCell;
this->MaximumPoints = this->CurrentExtraPoint;
}
return 1;
}
//----------------------------------------------------------------------------
// Shift data if center longitude needs to change.
//----------------------------------------------------------------------------
void vtkMPASReader::ShiftLonData()
{
vtkDebugMacro(<< "In ShiftLonData..." << endl);
// if atmospheric data, or zero centered, set center to 180 instead of 0
if (IsAtmosphere || IsZeroCentered)
{
for (int j = this->PointOffset; j < this->NumberOfPoints + this->PointOffset; j++)
{
// need to shift over the point so center is at PI
if (this->PointX[j] < 0)
{
this->PointX[j] += 2*vtkMath::Pi();
}
}
}
if (CenterLon != 180)
{
for (int j = this->PointOffset; j < this->NumberOfPoints + this->PointOffset; j++)
{
// need to shift over the point if centerLon dictates
if (this->CenterRad < vtkMath::Pi())
{
if (this->PointX[j] > (this->CenterRad + vtkMath::Pi()))
{
this->PointX[j] = -((2*vtkMath::Pi()) - this->PointX[j]);
}
}
else if (this->CenterRad > vtkMath::Pi())
{
if (this->PointX[j] < (this->CenterRad - vtkMath::Pi()))
{
this->PointX[j] += 2*vtkMath::Pi();
}
}
}
}
vtkDebugMacro(<< "Leaving ShiftLonData..." << endl);
}
//----------------------------------------------------------------------------
// Add a "mirror point" -- a point on the opposite side of the lat/lon
// projection.
//----------------------------------------------------------------------------
int vtkMPASReader::AddMirrorPoint(int index, double dividerX, double offset)
{
double X = this->PointX[index];
double Y = this->PointY[index];
// add on east
if (X < dividerX)
{
X += offset;
}
else
{
// add on west
X -= offset;
}
assert(this->CurrentExtraPoint < this->ModNumPoints);
this->PointX[this->CurrentExtraPoint] = X;
this->PointY[this->CurrentExtraPoint] = Y;
int mirrorPoint = this->CurrentExtraPoint;
// record mapping
*(this->PointMap + (this->CurrentExtraPoint - this->NumberOfPoints - this->PointOffset)) = index;
this->CurrentExtraPoint++;
return mirrorPoint;
}
//----------------------------------------------------------------------------
// Check for out-of-range values and do bugfix
//----------------------------------------------------------------------------
void vtkMPASReader::FixPoints()
{
vtkDebugMacro(<< "In FixPoints..." << endl);
for (int j = this->CellOffset; j < this->NumberOfCells + this->CellOffset; j++ )
{
int *conns = this->OrigConnections + (j * this->PointsPerCell);
// go through and make sure none of the referenced points are
// out of range
// if so, set all to point 0
for (int k = 0; k < this->PointsPerCell; k++)
{
if ((conns[k] <= 0) || (conns[k] > this->NumberOfPoints))
{
for (int m = 0; m < this->PointsPerCell; m++)
{
conns[m] = 0;
}
break;
}
}
if (this->DoBugFix)
{
//BUG FIX for problem where cells are stretching to a faraway point
int lastk = this->PointsPerCell-1;
const double thresh = .06981317007977; // 4 degrees
for (int k = 0; k < this->PointsPerCell; k++)
{
double ydiff = abs(this->PointY[conns[k]]
- this->PointY[conns[lastk]]);
// Don't look at cells at map border
if (ydiff > thresh)
{
for (int m = 0; m < this->PointsPerCell; m++)
{
conns[m] = 0;
}
break;
}
}
}
}
vtkDebugMacro(<< "Leaving FixPoints..." << endl);
}
//----------------------------------------------------------------------------
// Eliminate wraparound at east/west edges of lat/lon projection
//----------------------------------------------------------------------------
int vtkMPASReader::EliminateXWrap()
{
if (this->NumberOfPoints == 0)
{
return 1;
}
double xLength;
double xCenter;
switch (this->Geometry)
{
case vtkMPASReader::Spherical:
vtkErrorMacro("EliminateXWrap called for spherical geometry.");
return 0;
case vtkMPASReader::Projected:
xLength = 2 * vtkMath::Pi();
xCenter = this->CenterRad;
break;
case vtkMPASReader::Planar:
{
// Determine the bounds in the x-dimension
double xRange[2] = { this->PointX[this->PointOffset],
this->PointX[this->PointOffset] };
for (int i = 1; i < this->NumberOfPoints; ++i)
{
double x = this->PointX[this->PointOffset + i];
xRange[0] = std::min(xRange[0], x);
xRange[1] = std::max(xRange[1], x);
}
xLength = xRange[1] - xRange[0];
xCenter = (xRange[0] + xRange[1]) * 0.5;
}
break;
default:
vtkErrorMacro("Unrecognized geometry type (" << this->Geometry << ").");
return 0;
}
// Assume a point is wrapped if it is more than 3/4 across the dataset:
double tolerance = xLength * 0.75;
// For each cell, examine vertices
// Add new points and cells where needed to account for wraparound.
for (int j = this->CellOffset; j < this->NumberOfCells + this->CellOffset; j++ )
{
int *conns = this->OrigConnections + (j * this->PointsPerCell);
int *modConns = this->ModConnections + (j * this->PointsPerCell);
// Determine if we are wrapping in X direction
int lastk = this->PointsPerCell-1;
bool xWrap = false;
for (int k = 0; k < this->PointsPerCell; k++)
{
if (abs(this->PointX[conns[k]] - this->PointX[conns[lastk]]) > tolerance)
{
xWrap = true;
break;
}
lastk = k;
}
// If we wrapped in X direction, modify cell and add mirror cell
if (xWrap)
{
// first point is anchor it doesn't move
double anchorX = this->PointX[conns[0]];
modConns[0] = conns[0];
// modify existing cell, so it doesn't wrap
// move points to one side
for (int k = 1; k < this->PointsPerCell; k++)
{
int neigh = conns[k];
// add a new point, figure out east or west
if (abs(this->PointX[neigh] - anchorX) > tolerance)
{
modConns[k] = this->AddMirrorPoint(neigh, anchorX, xLength);
}
else
{
// use existing kth point
modConns[k] = neigh;
}
}
// move addedConns to this->ModConnections extra cells area
int* addedConns = this->ModConnections
+ (this->CurrentExtraCell * this->PointsPerCell);
// add a mirroring cell to other side
// add mirrored anchor first
addedConns[0] = this->AddMirrorPoint(conns[0], xCenter, xLength);
anchorX = this->PointX[addedConns[0]];
// add mirror cell points if needed
for (int k = 1; k < this->PointsPerCell; k++)
{
int neigh = conns[k];
// add a new point for neighbor, figure out east or west
if (abs(this->PointX[neigh] - anchorX) > tolerance)
{
addedConns[k] = this->AddMirrorPoint(neigh, anchorX, xLength);
}
else
{
// use existing kth point
addedConns[k] = neigh;
}
}
*(this->CellMap + (this->CurrentExtraCell - this->NumberOfCells - this->CellOffset)) = j;
this->CurrentExtraCell++;
}
else
{
// just add cell "as is" to this->ModConnections
for (int k=0; k< this->PointsPerCell; k++)
{
modConns[k] = conns[k];
}
}
if (this->CurrentExtraCell > this->ModNumCells)
{
vtkErrorMacro( << "Exceeded storage for extra cells!" << endl);
return(0);
}
if (this->CurrentExtraPoint > this->ModNumPoints)
{
vtkErrorMacro( << "Exceeded storage for extra points!" << endl);
return(0);
}
}
if (!ShowMultilayerView)
{
this->MaximumCells = this->CurrentExtraCell;
this->MaximumPoints = this->CurrentExtraPoint;
vtkDebugMacro
(<< "elim xwrap: singlelayer: setting this->MaximumPoints to " << this->MaximumPoints
<< endl);
}
else
{
this->MaximumCells = this->CurrentExtraCell*this->MaximumNVertLevels;
this->MaximumPoints = this->CurrentExtraPoint*(this->MaximumNVertLevels+1);
vtkDebugMacro
(<< "elim xwrap: multilayer: setting this->MaximumPoints to " <<
this->MaximumPoints << endl);
}
return 1;
}
//----------------------------------------------------------------------------
// Add points to vtk data structures
//----------------------------------------------------------------------------
void vtkMPASReader::OutputPoints()
{
vtkUnstructuredGrid* output = this->GetOutput();
double adjustedLayerThickness = this->IsAtmosphere
? static_cast<double>(-this->LayerThickness)
: static_cast<double>(this->LayerThickness);
vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
points->Allocate(this->MaximumPoints);
output->SetPoints(points);
for (int j = 0; j < this->CurrentExtraPoint; j++ )
{
double x, y, z;
switch (this->Geometry)
{
case vtkMPASReader::Planar:
case vtkMPASReader::Spherical:
x = this->PointX[j];
y = this->PointY[j];
z = this->PointZ[j];
break;
case vtkMPASReader::Projected:
x = this->PointX[j] * 180.0 / vtkMath::Pi();
y = this->PointY[j] * 180.0 / vtkMath::Pi();
z = 0.0;
break;
default:
vtkErrorMacro("Unrecognized geometry type (" << this->Geometry << ").");
return;
}
if (!this->ShowMultilayerView)
{
points->InsertNextPoint(x, y, z);
}
else
{
double rho=0.0, rholevel=0.0, theta=0.0, phi=0.0;
int retval = -1;
if (this->Geometry == Spherical)
{
if ((x != 0.0) || (y != 0.0) || (z != 0.0))
{
retval = CartesianToSpherical(x, y, z, &rho, &phi, &theta);
if (retval)
{
vtkWarningMacro("Can't create point for layered view.");
}
}
}
for (int levelNum = 0; levelNum < this->MaximumNVertLevels+1; levelNum++)
{
if (this->Geometry == Spherical)
{
if (!retval && ((x != 0.0) || (y != 0.0) || (z != 0.0)))
{
rholevel = rho - (adjustedLayerThickness * levelNum);
retval = SphericalToCartesian(rholevel, phi, theta, &x, &y, &z);
if (retval)
{
vtkWarningMacro("Can't create point for layered view.");
}
}
}
else
{
z = -levelNum * adjustedLayerThickness;
}
points->InsertNextPoint(x, y, z);
}
}
}
if (this->PointX)
{
free(this->PointX);
this->PointX = NULL;
}
if (this->PointY)
{
free(this->PointY);
this->PointY = NULL;
}
if (this->PointZ)
{
free(this->PointZ);
this->PointZ = NULL;
}
}
//----------------------------------------------------------------------------
// Determine if cell is one of VTK_TRIANGLE, VTK_WEDGE, VTK_QUAD or
// VTK_HEXAHEDRON
//----------------------------------------------------------------------------
unsigned char vtkMPASReader::GetCellType()
{
// write cell types
unsigned char cellType = VTK_TRIANGLE;
switch (this->PointsPerCell) {
case 3:
if (!ShowMultilayerView)
{
cellType = VTK_TRIANGLE;
}
else
{
cellType = VTK_WEDGE;
}
break;
case 4:
if (!ShowMultilayerView)
{
cellType = VTK_QUAD;
}
else
{
cellType = VTK_HEXAHEDRON;
}
break;
default:
break;
}
return cellType;
}
//------------------------------------------------------------------------------
bool vtkMPASReader::ValidateDimensions(NcVar *var, bool silent, int ndims, ...)
{
if (var->num_dims() != ndims)
{
if (!silent)
{
vtkWarningMacro(<< "Expected variable '" << var->name() << "' to have "
<< ndims << "dimension(s), but it has " << var->num_dims()
<< ".");
}
return false;
}
va_list args;
va_start(args, ndims);
for (int i = 0; i < ndims; ++i)
{
NcDim *dim = var->get_dim(i);
std::string dimName(va_arg(args, const char *));
if (dimName != dim->name())
{
if (!silent)
{
vtkWarningMacro(<< "Expected variable '" << var->name() << "' to have '"
<< dimName << "' at dimension index " << i << ", not '"
<< dim->name() << "'.");
}
va_end(args);
return false;
}
}
va_end(args);
return true;
}
//------------------------------------------------------------------------------
long vtkMPASReader::GetCursorForDimension(const NcDim *dim)
{
std::string dimName = dim->name();
if (dimName == "nCells" || dimName == "nVertices")
{
return 0;
}
else if (dimName == "Time")
{
return std::min(static_cast<long>(std::floor(this->DTime)),
static_cast<long>(this->NumberOfTimeSteps-1));
}
else if (this->ShowMultilayerView &&
dimName == this->VerticalDimension)
{
return 0;
}
else
{
return this->InitializeDimension(dim);
}
}
//------------------------------------------------------------------------------
size_t vtkMPASReader::GetCountForDimension(const NcDim *dim)
{
std::string dimName = dim->name();
if (dimName == "nCells")
{
return this->NumberOfPoints;
}
else if (dimName == "nVertices")
{
return this->NumberOfCells;
}
else if (this->ShowMultilayerView && dimName == this->VerticalDimension)
{
return this->MaximumNVertLevels;
}
else
{
return 1;
}
}
//------------------------------------------------------------------------------
long vtkMPASReader::InitializeDimension(const NcDim *dim)
{
Internal::DimMetaDataMap::const_iterator match =
this->Internals->dimMetaDataMap.find(dim->name());
long result = 0;
if (match == this->Internals->dimMetaDataMap.end())
{
DimMetaData metaData;
metaData.curIdx = result;
metaData.dimSize = dim->size();
this->Internals->dimMetaDataMap.insert(
std::make_pair(std::string(dim->name()), metaData));
this->Internals->dimMetaDataTime.Modified();
}
else
{
result = match->second.curIdx;
}
return result;
}
//----------------------------------------------------------------------------
// Add cells to vtk data structures
//----------------------------------------------------------------------------
void vtkMPASReader::OutputCells()
{
vtkDebugMacro(<< "In OutputCells..." << endl);
vtkUnstructuredGrid* output = GetOutput();
output->Allocate(this->MaximumCells, this->MaximumCells);
int cellType = GetCellType();
int val;
int pointsPerPolygon;
if (this->ShowMultilayerView)
{
pointsPerPolygon = 2 * this->PointsPerCell;
}
else
{
pointsPerPolygon = this->PointsPerCell;
}
vtkDebugMacro
(<< "OutputCells: this->MaximumCells: " << this->MaximumCells
<< " cellType: " << cellType << " this->MaximumNVertLevels: " << this->MaximumNVertLevels
<< " LayerThickness: " << LayerThickness << " ProjectLatLon: "
<< ProjectLatLon << " ShowMultilayerView: " << ShowMultilayerView);
std::vector<vtkIdType> polygon(pointsPerPolygon);
for (int j = 0; j < this->CurrentExtraCell ; j++)
{
int* conns;
if (this->Geometry == Spherical)
{
conns = this->OrigConnections + (j * this->PointsPerCell);
}
else
{
conns = this->ModConnections + (j * this->PointsPerCell);
}
int minLevel= 0;
if (this->IncludeTopography)
{
int* connections;
//check if it is a mirror cell, if so, get original
if (j >= this->NumberOfCells + this->CellOffset)
{
int origCellNum = *(this->CellMap + (j - this->NumberOfCells - this->CellOffset));
connections = this->OrigConnections + (origCellNum*this->PointsPerCell);
}
else
{
connections = this->OrigConnections + (j * this->PointsPerCell);
}
minLevel = this->MaximumLevelPoint[connections[0]];
// Take the min of the this->MaximumLevelPoint of each point
for (int k = 1; k < this->PointsPerCell; k++)
{
minLevel = min(minLevel, this->MaximumLevelPoint[connections[k]]);
}
}
// singlelayer
if (!this->ShowMultilayerView)
{
// If that min is greater than or equal to this output level,
// include the cell, otherwise set all points to zero.
if (this->IncludeTopography && ((minLevel-1) < this->GetVerticalLevel()))
{
//cerr << "Setting all points to zero" << endl;
val = 0;
for (int k = 0; k < this->PointsPerCell; k++)
{
polygon[k] = val;
}
}
else
{
for (int k = 0; k < this->PointsPerCell; k++)
{
polygon[k] = conns[k];
}
}
output->InsertNextCell(cellType, pointsPerPolygon, &polygon[0]);
}
else
{ // multilayer
// for each level, write the cell
for (int levelNum = 0; levelNum < this->MaximumNVertLevels; levelNum++)
{
if (this->IncludeTopography && ((minLevel-1) < levelNum))
{
// setting all points to zero
val = 0;
for (int k = 0; k < pointsPerPolygon; k++)
{
polygon[k] = val;
}
}
else
{
for (int k = 0; k < this->PointsPerCell; k++)
{
val = (conns[k]*(this->MaximumNVertLevels+1)) + levelNum;
polygon[k] = val;
}
for (int k = 0; k < this->PointsPerCell; k++)
{
val = (conns[k]*(this->MaximumNVertLevels+1)) + levelNum + 1;
polygon[k+this->PointsPerCell] = val;
}
}
//vtkDebugMacro
//("InsertingCell j: " << j << " level: " << levelNum << endl);
output->InsertNextCell(cellType, pointsPerPolygon, &polygon[0]);
}
}
}
free(this->ModConnections); this->ModConnections = NULL;
free(this->OrigConnections); this->OrigConnections = NULL;
vtkDebugMacro(<< "Leaving OutputCells..." << endl);
}
//------------------------------------------------------------------------------
vtkDataArray *vtkMPASReader::CreateDataArray(int typeNc)
{
int typeVtk = this->NcTypeToVtkType(static_cast<NcType>(typeNc));
return vtkDataArray::CreateDataArray(typeVtk);
}
//------------------------------------------------------------------------------
vtkIdType vtkMPASReader::ComputeNumberOfTuples(NcVar *ncVar)
{
int numDims = ncVar->num_dims();
vtkIdType size = 0;
for (int dim = 0; dim < numDims; ++dim)
{
vtkIdType count = static_cast<vtkIdType>(
this->GetCountForDimension(ncVar->get_dim(dim)));
if (size == 0)
{
size = count;
}
else
{
size *= count;
}
}
return size;
}
//------------------------------------------------------------------------------
template <typename ValueType>
bool vtkMPASReader::LoadDataArray(NcVar *ncVar, vtkDataArray *array,
bool resize)
{
if (array->GetDataType() != this->NcTypeToVtkType(ncVar->type()))
{
vtkWarningMacro("Invalid array type.");
return false;
}
int numDims = ncVar->num_dims();
std::vector<long> cursor;
std::vector<nc_size_t> counts;
vtkIdType size = 0;
for (int dim = 0; dim < numDims; ++dim)
{
cursor.push_back(this->GetCursorForDimension(ncVar->get_dim(dim)));
counts.push_back(static_cast<nc_size_t>(
this->GetCountForDimension(ncVar->get_dim(dim))));
if (size == 0)
{
size = counts.back();
}
else
{
size *= counts.back();
}
}
if (resize)
{
array->SetNumberOfComponents(1);
array->SetNumberOfTuples(size);
}
else
{
if (array->GetNumberOfComponents() != 1)
{
vtkWarningMacro("Invalid number of components: "
<< array->GetNumberOfComponents() << ".");
return false;
}
else if (array->GetNumberOfTuples() < size)
{
vtkWarningMacro("Array only has " << array->GetNumberOfTuples()
<< " allocated, but we need " << size << ".");
return false;
}
}
ValueType *dataBlock = static_cast<ValueType*>(array->GetVoidPointer(0));
if (!dataBlock)
{
vtkWarningMacro("GetVoidPointer returned NULL.");
return false;
}
if (!ncVar->set_cur(&cursor[0]))
{
vtkWarningMacro("Setting cursor failed.");
return false;
}
if (!ncVar->get(dataBlock, &counts[0]))
{
vtkWarningMacro("Reading " << size << " elements failed.");
return false;
}
return true;
}
//------------------------------------------------------------------------------
template <typename ValueType>
int vtkMPASReader::LoadPointVarDataImpl(NcVar *ncVar, vtkDataArray *array)
{
// Don't resize, we've preallocated extra room for multilayer (if needed):
if (!this->LoadDataArray<ValueType>(ncVar, array, /*resize=*/false))
{
return 0;
}
// Check if this variable contains the vertical dimension:
bool hasVerticalDimension = false;
int numDims = ncVar->num_dims();
if (this->ShowMultilayerView)
{
for (int d = 0; d < numDims; ++d)
{
if (this->VerticalDimension == ncVar->get_dim(d)->name())
{
hasVerticalDimension = true;
break;
}
}
}
vtkIdType varSize = this->ComputeNumberOfTuples(ncVar);
ValueType *dataBlock = static_cast<ValueType*>(array->GetVoidPointer(0));
std::vector<ValueType> tempData; // Used for Multilayer
// singlelayer
if (!this->ShowMultilayerView)
{
// Account for point offset:
if (this->PointOffset != 0)
{
assert(this->NumberOfPoints <= array->GetNumberOfTuples() &&
"Source array too small.");
assert(this->PointOffset + this->NumberOfPoints <=
array->GetNumberOfTuples() && "Destination array too small.");
if (this->PointOffset < this->NumberOfPoints)
{
std::copy_backward(dataBlock, dataBlock + this->NumberOfPoints,
dataBlock + this->PointOffset +
this->NumberOfPoints);
}
else
{
std::copy(dataBlock, dataBlock + this->NumberOfPoints,
dataBlock + this->PointOffset);
}
}
dataBlock[0] = dataBlock[1];
// data is all in place, don't need to do next step
}
else
{ // multilayer
if (this->MaximumPoints == 0)
{
return 0; // No points
}
tempData.resize(this->MaximumPoints);
size_t vertPointOffset = this->MaximumNVertLevels * this->PointOffset;
ValueType *dataPtr = &tempData[0] + vertPointOffset;
assert(varSize < array->GetNumberOfTuples());
assert(varSize < static_cast<vtkIdType>(this->MaximumPoints -
vertPointOffset));
std::copy(dataBlock, dataBlock + varSize, dataPtr);
if (!hasVerticalDimension)
{
// need to replicate data over all vertical layers
// layout in memory needs to be:
// pt1, pt1, ..., (VertLevels times), pt2, pt2, ..., (VertLevels times),
// need to go backwards through the points in order to not overwrite
// anything.
for(int i = this->NumberOfPoints; i > 0; i--)
{
// point to copy
ValueType pt = *(dataPtr + i - 1);
// where to start copying
ValueType *copyPtr = dataPtr + (i-1)*this->MaximumNVertLevels;
std::fill(copyPtr, copyPtr + this->MaximumNVertLevels, pt);
}
}
}
vtkDebugMacro(<<"Got point data.");
int i = 0;
int k = 0;
if (this->ShowMultilayerView)
{
// put in dummy points
assert(this->MaximumNVertLevels * 2 <= this->MaximumPoints);
assert(this->MaximumNVertLevels <= array->GetNumberOfTuples());
std::copy(tempData.begin() + this->MaximumNVertLevels,
tempData.begin() + (2*this->MaximumNVertLevels),
dataBlock);
// write highest level dummy point (duplicate of last level)
assert(this->MaximumNVertLevels < array->GetNumberOfTuples());
assert(2*this->MaximumNVertLevels - 1 < this->MaximumPoints);
dataBlock[this->MaximumNVertLevels] =
tempData[2*this->MaximumNVertLevels - 1];
vtkDebugMacro(<<"Wrote dummy point data.");
// put in other points
for (int j = this->PointOffset;
j < this->NumberOfPoints + this->PointOffset;
j++)
{
i = j*(this->MaximumNVertLevels+1);
k = j*(this->MaximumNVertLevels);
// write data for one point -- lowest level to highest
assert(k + this->MaximumNVertLevels <= this->MaximumPoints);
assert(i + this->MaximumNVertLevels <= array->GetNumberOfTuples());
std::copy(tempData.begin() + k,
tempData.begin() + k + this->MaximumNVertLevels,
dataBlock + i);
// for last layer of points, repeat last level's values
// Need Mark's input on this one
dataBlock[i++] = tempData[--k];
//vtkDebugMacro (<< "Wrote j:" << j << endl);
}
}
vtkDebugMacro(<<"Wrote next points.");
vtkDebugMacro(<<"this->NumberOfPoints: " << this->NumberOfPoints << " "
<<"this->CurrentExtraPoint: " << this->CurrentExtraPoint);
// put out data for extra points
for (int j = this->PointOffset + this->NumberOfPoints;
j < this->CurrentExtraPoint;
j++)
{
// use map to find out what point data we are using
if (!this->ShowMultilayerView)
{
k = this->PointMap[j - this->NumberOfPoints - this->PointOffset];
assert(j < array->GetNumberOfTuples());
assert(k < array->GetNumberOfTuples());
dataBlock[j] = dataBlock[k];
}
else
{
k = this->PointMap[j - this->NumberOfPoints - this->PointOffset] *
this->MaximumNVertLevels;
// write data for one point -- lowest level to highest
assert(k + this->MaximumNVertLevels <= this->MaximumPoints);
assert(i + this->MaximumNVertLevels <= array->GetNumberOfTuples());
std::copy(tempData.begin() + k,
tempData.begin() + k + this->MaximumNVertLevels,
dataBlock + i);
// for last layer of points, repeat last level's values
// Need Mark's input on this one
dataBlock[i++] = tempData[--k];
}
}
vtkDebugMacro(<<"wrote extra point data.");
return 1;
}
//----------------------------------------------------------------------------
// Load the data for a point variable
//----------------------------------------------------------------------------
vtkDataArray *vtkMPASReader::LoadPointVarData(int variableIndex)
{
NcVar* ncVar = this->Internals->pointVars[variableIndex];
if (ncVar == NULL)
{
vtkErrorMacro(<<"No NetCDF data for pointVar @ index " << variableIndex);
return 0;
}
vtkDebugMacro(<<"Loading point data array named: " << ncVar->name());
// Get data type:
NcType typeNc = ncVar->type();
int typeVtk = this->NcTypeToVtkType(typeNc);
// Allocate data array pointer for this variable:
vtkSmartPointer<vtkDataArray> array =
this->LookupPointDataArray(variableIndex);
if (array == NULL)
{
vtkDebugMacro(<<"Allocating data array.");
array = vtkSmartPointer<vtkDataArray>::Take(
vtkDataArray::CreateDataArray(typeVtk));
}
array->SetName(ncVar->name());
array->SetNumberOfComponents(1);
array->SetNumberOfTuples(this->MaximumPoints);
int success = false;
vtkNcDispatch(typeVtk,
success = this->LoadPointVarDataImpl<VTK_TT>(ncVar, array););
if (success)
{
this->Internals->pointArrays[variableIndex] = array;
return array;
}
return NULL;
}
//------------------------------------------------------------------------------
template <typename ValueType>
int vtkMPASReader::LoadCellVarDataImpl(NcVar *ncVar, vtkDataArray *array)
{
// Don't resize, we've preallocated extra room for multilayer (if needed):
if (!this->LoadDataArray<ValueType>(ncVar, array, /*resize=*/false))
{
return 0;
}
ValueType *dataBlock = static_cast<ValueType*>(array->GetVoidPointer(0));
// put out data for extra cells
for (int j = this->CellOffset + this->NumberOfCells;
j < this->CurrentExtraCell;
j++)
{
// use map to find out what cell data we are using
if (!this->ShowMultilayerView)
{
int k = this->CellMap[j - this->NumberOfCells - this->CellOffset];
assert(j < array->GetNumberOfTuples());
assert(k < array->GetNumberOfTuples());
dataBlock[j] = dataBlock[k];
}
else
{
int i = j*this->MaximumNVertLevels;
int k = this->CellMap[j - this->NumberOfCells - this->CellOffset] *
this->MaximumNVertLevels;
// write data for one cell -- lowest level to highest
assert(i < array->GetNumberOfTuples());
assert(k + this->MaximumNVertLevels <= array->GetNumberOfTuples());
std::copy(dataBlock + k, dataBlock + k + this->MaximumNVertLevels,
dataBlock + i);
}
}
vtkDebugMacro(<<"Stored data.");
return 1;
}
//----------------------------------------------------------------------------
// Load the data for a cell variable
//----------------------------------------------------------------------------
vtkDataArray* vtkMPASReader::LoadCellVarData(int variableIndex)
{
NcVar* ncVar = this->Internals->cellVars[variableIndex];
if (ncVar == NULL)
{
vtkErrorMacro(<<"No NetCDF data for cellVar @ index " << variableIndex);
return 0;
}
vtkDebugMacro(<<"Loading cell data array named: " << ncVar->name());
// Get data type:
NcType typeNc = ncVar->type();
int typeVtk = this->NcTypeToVtkType(typeNc);
// Allocate data array pointer for this variable:
vtkSmartPointer<vtkDataArray> array =
this->LookupCellDataArray(variableIndex);
if (array == NULL)
{
vtkDebugMacro(<<"Allocating data array.");
array = vtkSmartPointer<vtkDataArray>::Take(
vtkDataArray::CreateDataArray(typeVtk));
}
array->SetName(ncVar->name());
array->SetNumberOfComponents(1);
array->SetNumberOfTuples(this->MaximumCells);
int success = false;
vtkNcDispatch(typeVtk,
success = this->LoadCellVarDataImpl<VTK_TT>(ncVar, array););
if (success)
{
this->Internals->cellArrays[variableIndex] = array;
return array;
}
return NULL;
}
//------------------------------------------------------------------------------
vtkDataArray *vtkMPASReader::LookupPointDataArray(int varIdx)
{
Internal::ArrayMap::iterator it = this->Internals->pointArrays.find(varIdx);
return it != this->Internals->pointArrays.end() ? it->second : NULL;
}
//------------------------------------------------------------------------------
vtkDataArray *vtkMPASReader::LookupCellDataArray(int varIdx)
{
Internal::ArrayMap::iterator it = this->Internals->cellArrays.find(varIdx);
return it != this->Internals->cellArrays.end() ? it->second : NULL;
}
//------------------------------------------------------------------------------
void vtkMPASReader::LoadTimeFieldData(vtkUnstructuredGrid *dataset)
{
vtkStringArray *array = NULL;
vtkFieldData *fd = dataset->GetFieldData();
if (!fd)
{
fd = vtkFieldData::New();
dataset->SetFieldData(fd);
fd->Delete();
}
if (vtkDataArray *da = fd->GetArray("Time"))
{
if (!(array = vtkArrayDownCast<vtkStringArray>(da)))
{
vtkWarningMacro("Not creating \"Time\" field data array: a data array "
"with this name already exists.");
return;
}
}
if (!array)
{
array = vtkStringArray::New();
array->SetName("Time");
fd->AddArray(array);
array->Delete();
}
// If the xtime variable exists, use its value at the current timestep:
std::string time;
if (isNcVar(this->Internals->ncFile, "xtime"))
{
NcVar *var = this->Internals->ncFile->get_var("xtime");
if (var && this->ValidateDimensions(var, false, 2, "Time", "StrLen"))
{
NcDim *strLenDim = this->Internals->ncFile->get_dim("StrLen");
assert(strLenDim);
long strLen = strLenDim->size();
if (strLen > 0)
{
time.resize(strLen);
var->set_cur(this->GetCursorForDimension(var->get_dim(0)), 0);
if (var->get(&time[0], 1, strLen))
{
// Trim off trailing whitespace:
size_t realLength = time.find_last_not_of(" ");
if (realLength != vtkStdString::npos)
{
time.resize(realLength + 1);
}
}
else
{
vtkWarningMacro("Error reading xtime variable from file.");
time.clear();
}
}
}
}
// If no string time is available or the read fails, just insert the timestep:
if (time.empty())
{
std::ostringstream timeStr;
timeStr << "Timestep " << std::floor(this->DTime)
<< "/" << this->NumberOfTimeSteps;
time = timeStr.str();
}
assert(array != NULL);
array->SetNumberOfComponents(1);
array->SetNumberOfTuples(1);
array->SetValue(0, time);
}
//------------------------------------------------------------------------------
inline int vtkMPASReader::NcTypeToVtkType(int ncType)
{
switch (static_cast<NcType>(ncType))
{
case ncByte:
return VTK_SIGNED_CHAR;
case ncChar:
return VTK_CHAR;
case ncShort:
return VTK_SHORT;
case ncInt:
return VTK_INT;
case ncFloat:
return VTK_FLOAT;
case ncDouble:
return VTK_DOUBLE;
case ncNoType:
default: // Shouldn't happen...
vtkGenericWarningMacro(<<"Invalid NcType: " << ncType);
return VTK_VOID;
}
}
//----------------------------------------------------------------------------
// Callback if the user selects a variable.
//----------------------------------------------------------------------------
void vtkMPASReader::SelectionCallback(vtkObject*,
unsigned long vtkNotUsed(eventid),
void* clientdata,
void* vtkNotUsed(calldata))
{
static_cast<vtkMPASReader*>(clientdata)->Modified();
}
//----------------------------------------------------------------------------
void vtkMPASReader::UpdateDimensions(bool force)
{
if (!force &&
this->Internals->dimMetaDataTime < this->Internals->extraDimTime)
{
return;
}
this->Internals->extraDims->Reset();
if (!this->Internals->ncFile)
{
this->Internals->extraDimTime.Modified();
return;
}
std::set<std::string> dimSet;
typedef Internal::DimMetaDataMap::const_iterator Iter;
const Internal::DimMetaDataMap &map = this->Internals->dimMetaDataMap;
for (Iter it = map.begin(), itEnd = map.end(); it != itEnd; ++it)
{
if (this->Internals->isExtraDim(it->first))
{
dimSet.insert(it->first);
}
}
typedef std::set<std::string>::const_iterator SetIter;
this->Internals->extraDims->Allocate(dimSet.size());
for (SetIter it = dimSet.begin(), itEnd = dimSet.end(); it != itEnd; ++it)
{
this->Internals->extraDims->InsertNextValue(it->c_str());
}
this->Internals->extraDimTime.Modified();
}
//----------------------------------------------------------------------------
// Return the output.
//----------------------------------------------------------------------------
vtkUnstructuredGrid* vtkMPASReader::GetOutput()
{
return this->GetOutput(0);
}
//----------------------------------------------------------------------------
// Returns the output given an id.
//----------------------------------------------------------------------------
vtkUnstructuredGrid* vtkMPASReader::GetOutput(int idx)
{
if (idx)
{
return NULL;
}
else
{
return vtkUnstructuredGrid::SafeDownCast( this->GetOutputDataObject(idx) );
}
}
//----------------------------------------------------------------------------
// Get number of point arrays.
//----------------------------------------------------------------------------
int vtkMPASReader::GetNumberOfPointArrays()
{
return this->PointDataArraySelection->GetNumberOfArrays();
}
//----------------------------------------------------------------------------
// Get number of cell arrays.
//----------------------------------------------------------------------------
int vtkMPASReader::GetNumberOfCellArrays()
{
return this->CellDataArraySelection->GetNumberOfArrays();
}
//----------------------------------------------------------------------------
// Make all point selections available.
//----------------------------------------------------------------------------
void vtkMPASReader::EnableAllPointArrays()
{
this->PointDataArraySelection->EnableAllArrays();
}
//----------------------------------------------------------------------------
// Make all point selections unavailable.
//----------------------------------------------------------------------------
void vtkMPASReader::DisableAllPointArrays()
{
this->PointDataArraySelection->DisableAllArrays();
}
//----------------------------------------------------------------------------
// Make all cell selections available.
//----------------------------------------------------------------------------
void vtkMPASReader::EnableAllCellArrays()
{
this->CellDataArraySelection->EnableAllArrays();
}
//----------------------------------------------------------------------------
// Make all cell selections unavailable.
//----------------------------------------------------------------------------
void vtkMPASReader::DisableAllCellArrays()
{
this->CellDataArraySelection->DisableAllArrays();
}
//----------------------------------------------------------------------------
// Get name of indexed point variable
//----------------------------------------------------------------------------
const char* vtkMPASReader::GetPointArrayName(int index)
{
return this->PointDataArraySelection->GetArrayName(index);
}
//----------------------------------------------------------------------------
// Get status of named point variable selection
//----------------------------------------------------------------------------
int vtkMPASReader::GetPointArrayStatus(const char* name)
{
return this->PointDataArraySelection->ArrayIsEnabled(name);
}
//----------------------------------------------------------------------------
// Set status of named point variable selection.
//----------------------------------------------------------------------------
void vtkMPASReader::SetPointArrayStatus(const char* name, int status)
{
if (status)
{
this->PointDataArraySelection->EnableArray(name);
}
else
{
this->PointDataArraySelection->DisableArray(name);
}
}
//----------------------------------------------------------------------------
// Get name of indexed cell variable
//----------------------------------------------------------------------------
const char* vtkMPASReader::GetCellArrayName(int index)
{
return this->CellDataArraySelection->GetArrayName(index);
}
//----------------------------------------------------------------------------
// Get status of named cell variable selection.
//----------------------------------------------------------------------------
int vtkMPASReader::GetCellArrayStatus(const char* name)
{
return this->CellDataArraySelection->ArrayIsEnabled(name);
}
//----------------------------------------------------------------------------
// Set status of named cell variable selection.
//----------------------------------------------------------------------------
void vtkMPASReader::SetCellArrayStatus(const char* name, int status)
{
if (status)
{
this->CellDataArraySelection->EnableArray(name);
}
else
{
this->CellDataArraySelection->DisableArray(name);
}
}
//----------------------------------------------------------------------------
int vtkMPASReader::GetNumberOfDimensions()
{
this->UpdateDimensions();
return this->Internals->extraDims->GetNumberOfTuples();
}
//----------------------------------------------------------------------------
std::string vtkMPASReader::GetDimensionName(int idx)
{
this->UpdateDimensions();
return this->Internals->extraDims->GetValue(idx);
}
//----------------------------------------------------------------------------
vtkStringArray *vtkMPASReader::GetAllDimensions()
{
this->UpdateDimensions();
return this->Internals->extraDims.GetPointer();
}
//----------------------------------------------------------------------------
int vtkMPASReader::GetDimensionCurrentIndex(const std::string &dim)
{
this->UpdateDimensions();
typedef Internal::DimMetaDataMap::const_iterator Iter;
Iter it = this->Internals->dimMetaDataMap.find(dim);
if (it == this->Internals->dimMetaDataMap.end())
{
return -1;
}
return it->second.curIdx;
}
//----------------------------------------------------------------------------
void vtkMPASReader::SetDimensionCurrentIndex(const string &dim, int idx)
{
this->UpdateDimensions();
typedef Internal::DimMetaDataMap::iterator Iter;
Iter it = this->Internals->dimMetaDataMap.find(dim);
if (it != this->Internals->dimMetaDataMap.end() &&
idx < it->second.dimSize)
{
it->second.curIdx = idx;
this->Modified();
}
}
//----------------------------------------------------------------------------
int vtkMPASReader::GetDimensionSize(const string &dim)
{
this->UpdateDimensions();
typedef Internal::DimMetaDataMap::const_iterator Iter;
Iter it = this->Internals->dimMetaDataMap.find(dim);
if (it == this->Internals->dimMetaDataMap.end())
{
return -1;
}
return it->second.dimSize;
}
//----------------------------------------------------------------------------
// Set vertical level to be viewed.
//----------------------------------------------------------------------------
void vtkMPASReader::SetVerticalLevel(int level)
{
this->SetDimensionCurrentIndex(this->VerticalDimension, level);
}
//------------------------------------------------------------------------------
int vtkMPASReader::GetVerticalLevel()
{
return this->GetDimensionCurrentIndex(this->VerticalDimension);
}
//----------------------------------------------------------------------------
// Set center longitude for lat/lon projection
//----------------------------------------------------------------------------
void vtkMPASReader::SetCenterLon(int val)
{
vtkDebugMacro(<< "SetCenterLon: is " << this->CenterLon << endl);
if (CenterLon != val)
{
this->CenterLon = val;
this->CenterRad = val * vtkMath::Pi() / 180.0;
this->Modified();
vtkDebugMacro(<< "SetCenterLon: set to " << this->CenterLon << endl);
vtkDebugMacro(<< "CenterRad set to " << this->CenterRad << endl);
}
}
//----------------------------------------------------------------------------
// Determine if this reader can read the given file (if it is an MPAS format)
// NetCDF file
//----------------------------------------------------------------------------
int vtkMPASReader::CanReadFile(const char *filename)
{
NcFile ncFile(filename);
if (!ncFile.is_valid())
{
return 0;
}
bool ret = true;
ret &= isNcDim(&ncFile, "nCells");
ret &= isNcDim(&ncFile, "nVertices");
ret &= isNcDim(&ncFile, "vertexDegree");
ret &= isNcDim(&ncFile, "Time");
return ret;
}
//------------------------------------------------------------------------------
vtkMTimeType vtkMPASReader::GetMTime()
{
vtkMTimeType result = this->Superclass::GetMTime();
result = std::max(result, this->CellDataArraySelection->GetMTime());
result = std::max(result, this->PointDataArraySelection->GetMTime());
// Excluded, as this just manages a cache:
// result = std::max(result, this->Internals->extraDimTime.GetMTime());
result = std::max(result, this->Internals->dimMetaDataTime.GetMTime());
return result;
}
//----------------------------------------------------------------------------
// Print self.
//----------------------------------------------------------------------------
void vtkMPASReader::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "FileName: "
<< (this->FileName?this->FileName:"NULL") << "\n";
os << indent << "VerticalLevelRange: "
<< this->VerticalLevelRange[0] << ","
<< this->VerticalLevelRange[1] << "\n";
os << indent << "this->MaximumPoints: " << this->MaximumPoints << "\n";
os << indent << "this->MaximumCells: " << this->MaximumCells << "\n";
os << indent << "ProjectLatLon: "
<< (this->ProjectLatLon?"ON":"OFF") << endl;
os << indent << "OnASphere: "
<< (this->OnASphere?"ON":"OFF") << endl;
os << indent << "ShowMultilayerView: "
<< (this->ShowMultilayerView?"ON":"OFF") << endl;
os << indent << "CenterLonRange: "
<< this->CenterLonRange[0] << "," << this->CenterLonRange[1] << endl;
os << indent << "IsAtmosphere: "
<< (this->IsAtmosphere?"ON":"OFF") << endl;
os << indent << "IsZeroCentered: "
<< (this->IsZeroCentered?"ON":"OFF") << endl;
os << indent << "LayerThicknessRange: "
<< this->LayerThicknessRange[0] << "," << this->LayerThicknessRange[1] << endl;
}
//------------------------------------------------------------------------------
int vtkMPASReader::GetNumberOfCellVars()
{
return static_cast<int>(this->Internals->cellVars.size());
}
//------------------------------------------------------------------------------
int vtkMPASReader::GetNumberOfPointVars()
{
return static_cast<int>(this->Internals->pointVars.size());
}
|