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 2896 2897 2898 2899 2900 2901 2902 2903 2904 2905 2906 2907 2908 2909 2910 2911 2912 2913 2914 2915 2916 2917 2918 2919 2920 2921 2922 2923 2924 2925 2926 2927 2928 2929 2930 2931 2932 2933 2934 2935 2936 2937 2938 2939 2940 2941 2942 2943 2944 2945 2946 2947 2948 2949 2950 2951 2952 2953 2954 2955 2956 2957 2958 2959 2960 2961 2962 2963 2964 2965 2966 2967 2968 2969 2970 2971 2972 2973 2974 2975 2976 2977 2978 2979 2980 2981 2982 2983 2984 2985 2986 2987 2988 2989 2990 2991 2992 2993 2994 2995 2996 2997 2998 2999 3000 3001 3002 3003 3004 3005 3006 3007 3008 3009 3010 3011 3012 3013 3014 3015 3016 3017 3018 3019 3020 3021 3022 3023 3024 3025 3026 3027 3028 3029 3030 3031 3032 3033 3034 3035 3036 3037 3038 3039 3040 3041 3042 3043 3044 3045 3046 3047 3048 3049 3050 3051 3052 3053 3054 3055 3056 3057 3058 3059 3060 3061 3062 3063 3064 3065 3066 3067 3068 3069 3070 3071 3072 3073 3074 3075 3076 3077 3078 3079 3080 3081 3082 3083 3084 3085 3086 3087 3088 3089 3090 3091 3092 3093 3094 3095 3096 3097 3098 3099 3100 3101 3102 3103 3104 3105 3106 3107 3108 3109 3110 3111 3112 3113 3114 3115 3116 3117 3118 3119 3120 3121 3122 3123 3124 3125 3126 3127 3128 3129 3130 3131 3132 3133 3134 3135 3136 3137 3138 3139 3140 3141 3142 3143 3144 3145 3146 3147 3148 3149 3150 3151 3152 3153 3154 3155 3156 3157 3158 3159 3160 3161 3162 3163 3164 3165 3166 3167 3168 3169 3170 3171 3172 3173 3174 3175 3176 3177 3178 3179 3180 3181 3182 3183 3184 3185 3186 3187 3188 3189 3190 3191 3192 3193 3194 3195 3196 3197 3198 3199 3200 3201 3202 3203 3204 3205 3206 3207 3208 3209 3210 3211 3212 3213 3214 3215 3216 3217 3218 3219 3220 3221 3222 3223 3224 3225 3226 3227 3228 3229 3230 3231 3232 3233 3234 3235 3236 3237 3238 3239 3240 3241 3242 3243 3244 3245 3246 3247 3248 3249 3250 3251 3252 3253 3254 3255 3256 3257 3258 3259 3260 3261 3262 3263 3264 3265 3266 3267 3268 3269 3270 3271 3272 3273 3274 3275 3276 3277 3278 3279 3280 3281 3282 3283 3284 3285 3286 3287 3288 3289 3290 3291 3292 3293
|
<!doctype html>
<html lang="en">
<head>
<meta charset="utf-8">
<meta name="viewport" content="width=device-width, initial-scale=1, minimum-scale=1">
<meta name="generator" content="pdoc3 0.11.1">
<title>netCDF4 API documentation</title>
<meta name="description" content="Version 1.7.2
…">
<link rel="stylesheet" href="https://cdnjs.cloudflare.com/ajax/libs/10up-sanitize.css/13.0.0/sanitize.min.css" integrity="sha512-y1dtMcuvtTMJc1yPgEqF0ZjQbhnc/bFhyvIyVNb9Zk5mIGtqVaAB1Ttl28su8AvFMOY0EwRbAe+HCLqj6W7/KA==" crossorigin>
<link rel="stylesheet" href="https://cdnjs.cloudflare.com/ajax/libs/10up-sanitize.css/13.0.0/typography.min.css" integrity="sha512-Y1DYSb995BAfxobCkKepB1BqJJTPrOp3zPL74AWFugHHmmdcvO+C48WLrUOlhGMc0QG7AE3f7gmvvcrmX2fDoA==" crossorigin>
<link rel="stylesheet" href="https://cdnjs.cloudflare.com/ajax/libs/highlight.js/11.9.0/styles/default.min.css" crossorigin>
<style>:root{--highlight-color:#fe9}.flex{display:flex !important}body{line-height:1.5em}#content{padding:20px}#sidebar{padding:1.5em;overflow:hidden}#sidebar > *:last-child{margin-bottom:2cm}.http-server-breadcrumbs{font-size:130%;margin:0 0 15px 0}#footer{font-size:.75em;padding:5px 30px;border-top:1px solid #ddd;text-align:right}#footer p{margin:0 0 0 1em;display:inline-block}#footer p:last-child{margin-right:30px}h1,h2,h3,h4,h5{font-weight:300}h1{font-size:2.5em;line-height:1.1em}h2{font-size:1.75em;margin:2em 0 .50em 0}h3{font-size:1.4em;margin:1.6em 0 .7em 0}h4{margin:0;font-size:105%}h1:target,h2:target,h3:target,h4:target,h5:target,h6:target{background:var(--highlight-color);padding:.2em 0}a{color:#058;text-decoration:none;transition:color .2s ease-in-out}a:visited{color:#503}a:hover{color:#b62}.title code{font-weight:bold}h2[id^="header-"]{margin-top:2em}.ident{color:#900;font-weight:bold}pre code{font-size:.8em;line-height:1.4em;padding:1em;display:block}code{background:#f3f3f3;font-family:"DejaVu Sans Mono",monospace;padding:1px 4px;overflow-wrap:break-word}h1 code{background:transparent}pre{border-top:1px solid #ccc;border-bottom:1px solid #ccc;margin:1em 0}#http-server-module-list{display:flex;flex-flow:column}#http-server-module-list div{display:flex}#http-server-module-list dt{min-width:10%}#http-server-module-list p{margin-top:0}.toc ul,#index{list-style-type:none;margin:0;padding:0}#index code{background:transparent}#index h3{border-bottom:1px solid #ddd}#index ul{padding:0}#index h4{margin-top:.6em;font-weight:bold}@media (min-width:200ex){#index .two-column{column-count:2}}@media (min-width:300ex){#index .two-column{column-count:3}}dl{margin-bottom:2em}dl dl:last-child{margin-bottom:4em}dd{margin:0 0 1em 3em}#header-classes + dl > dd{margin-bottom:3em}dd dd{margin-left:2em}dd p{margin:10px 0}.name{background:#eee;font-size:.85em;padding:5px 10px;display:inline-block;min-width:40%}.name:hover{background:#e0e0e0}dt:target .name{background:var(--highlight-color)}.name > span:first-child{white-space:nowrap}.name.class > span:nth-child(2){margin-left:.4em}.inherited{color:#999;border-left:5px solid #eee;padding-left:1em}.inheritance em{font-style:normal;font-weight:bold}.desc h2{font-weight:400;font-size:1.25em}.desc h3{font-size:1em}.desc dt code{background:inherit}.source summary,.git-link-div{color:#666;text-align:right;font-weight:400;font-size:.8em;text-transform:uppercase}.source summary > *{white-space:nowrap;cursor:pointer}.git-link{color:inherit;margin-left:1em}.source pre{max-height:500px;overflow:auto;margin:0}.source pre code{font-size:12px;overflow:visible}.hlist{list-style:none}.hlist li{display:inline}.hlist li:after{content:',\2002'}.hlist li:last-child:after{content:none}.hlist .hlist{display:inline;padding-left:1em}img{max-width:100%}td{padding:0 .5em}.admonition{padding:.1em 1em;margin-bottom:1em}.admonition-title{font-weight:bold}.admonition.note,.admonition.info,.admonition.important{background:#aef}.admonition.todo,.admonition.versionadded,.admonition.tip,.admonition.hint{background:#dfd}.admonition.warning,.admonition.versionchanged,.admonition.deprecated{background:#fd4}.admonition.error,.admonition.danger,.admonition.caution{background:lightpink}</style>
<style media="screen and (min-width: 700px)">@media screen and (min-width:700px){#sidebar{width:30%;height:100vh;overflow:auto;position:sticky;top:0}#content{width:70%;max-width:100ch;padding:3em 4em;border-left:1px solid #ddd}pre code{font-size:1em}.name{font-size:1em}main{display:flex;flex-direction:row-reverse;justify-content:flex-end}.toc ul ul,#index ul ul{padding-left:1em}.toc > ul > li{margin-top:.5em}}</style>
<style media="print">@media print{#sidebar h1{page-break-before:always}.source{display:none}}@media print{*{background:transparent !important;color:#000 !important;box-shadow:none !important;text-shadow:none !important}a[href]:after{content:" (" attr(href) ")";font-size:90%}a[href][title]:after{content:none}abbr[title]:after{content:" (" attr(title) ")"}.ir a:after,a[href^="javascript:"]:after,a[href^="#"]:after{content:""}pre,blockquote{border:1px solid #999;page-break-inside:avoid}thead{display:table-header-group}tr,img{page-break-inside:avoid}img{max-width:100% !important}@page{margin:0.5cm}p,h2,h3{orphans:3;widows:3}h1,h2,h3,h4,h5,h6{page-break-after:avoid}}</style>
<script defer src="https://cdnjs.cloudflare.com/ajax/libs/highlight.js/11.9.0/highlight.min.js" integrity="sha512-D9gUyxqja7hBtkWpPWGt9wfbfaMGVt9gnyCvYa+jojwwPHLCzUm5i8rpk7vD7wNee9bA35eYIjobYPaQuKS1MQ==" crossorigin></script>
<script>window.addEventListener('DOMContentLoaded', () => {
hljs.configure({languages: ['bash', 'css', 'diff', 'graphql', 'ini', 'javascript', 'json', 'plaintext', 'python', 'python-repl', 'rust', 'shell', 'sql', 'typescript', 'xml', 'yaml']});
hljs.highlightAll();
})</script>
</head>
<body>
<main>
<article id="content">
<header>
<h1 class="title">Package <code>netCDF4</code></h1>
</header>
<section id="section-intro">
<h2 id="version-172">Version 1.7.2</h2>
<h1 id="introduction">Introduction</h1>
<p>netcdf4-python is a Python interface to the netCDF C library.</p>
<p><a href="http://www.unidata.ucar.edu/software/netcdf/">netCDF</a> version 4 has many features
not found in earlier versions of the library and is implemented on top of
<a href="http://www.hdfgroup.org/HDF5">HDF5</a>. This module can read and write
files in both the new netCDF 4 and the old netCDF 3 format, and can create
files that are readable by HDF5 clients. The API modelled after
<a href="http://dirac.cnrs-orleans.fr/ScientificPython/">Scientific.IO.NetCDF</a>,
and should be familiar to users of that module.</p>
<p>Most new features of netCDF 4 are implemented, such as multiple
unlimited dimensions, groups and data compression.
All the new
numeric data types (such as 64 bit and unsigned integer types) are
implemented. Compound (struct), variable length (vlen) and
enumerated (enum) data types are supported, but not the opaque data type.
Mixtures of compound, vlen and enum data types (such as
compound types containing enums, or vlens containing compound
types) are not supported.</p>
<h2 id="quick-install">Quick Install</h2>
<ul>
<li>the easiest way to get going is to install via <code>pip install <a title="netCDF4" href="#netCDF4">netCDF4</a></code>.
(or if you use the <a href="http://conda.io">conda</a> package manager <code>conda install -c conda-forge netCDF4</code>).</li>
</ul>
<h2 id="developer-install">Developer Install</h2>
<ul>
<li>Clone the <a href="http://github.com/Unidata/netcdf4-python">github repository</a>. Make
sure you either clone recursively, or run <code>git submodule update --init</code> to
ensure all the submodules are also checked out.</li>
<li>Make sure the dependencies are satisfied (Python 3.8 or later,
<a href="http://numpy.scipy.org">numpy</a>,
<a href="http://cython.org">Cython</a>,
<a href="https://github.com/Unidata/cftime">cftime</a>,
<a href="https://pypi.python.org/pypi/setuptools">setuptools</a>,
the <a href="https://www.hdfgroup.org/solutions/hdf5/">HDF5 C library</a>,
and the <a href="https://www.unidata.ucar.edu/software/netcdf/">netCDF C library</a>).
For MPI parallel IO support, an MPI-enabled versions of the netcdf library
is required, as is <a href="http://mpi4py.scipy.org">mpi4py</a>.
Parallel IO further depends on the existence of MPI-enabled HDF5 or the
<a href="https://parallel-netcdf.github.io/">PnetCDF</a> library.</li>
<li>By default, the utility <code>nc-config</code> (installed with netcdf-c)
will be run used to determine where all the dependencies live.</li>
<li>If <code>nc-config</code> is not in your default <code>PATH</code>, you can set the <code>NETCDF4_DIR</code>
environment variable and <code>setup.py</code> will look in <code>$NETCDF4_DIR/bin</code>.
You can also use the file <code>setup.cfg</code> to set the path to <code>nc-config</code>, or
enter the paths to the libraries and include files manually. Just
edit the <code>setup.cfg</code> file
in a text editor and follow the instructions in the comments.
To disable the use of <code>nc-config</code>, set the env var <code>USE_NCCONFIG</code> to 0.
To disable the use of <code>setup.cfg</code>, set <code>USE_SETUPCFG</code> to 0.
As a last resort, the library and include paths can be set via environment variables.
If you go this route, set <code>USE_NCCONFIG</code> and <code>USE_SETUPCFG</code> to 0, and specify
<code>NETCDF4_LIBDIR</code>, <code>NETCDF4_INCDIR</code>, <code>HDF5_LIBDIR</code> and <code>HDF5_INCDIR</code>.
If the dependencies are not found
in any of the paths specified by environment variables, then standard locations
(such as <code>/usr</code> and <code>/usr/local</code>) are searched.</li>
<li>if the env var <code>NETCDF_PLUGIN_DIR</code> is set to point to the location of the netcdf-c compression
plugins built by netcdf >= 4.9.0, they will be installed inside the package.
In this
case <code>HDF5_PLUGIN_PATH</code> will be set to the package installation path on import,
so the extra compression algorithms available in netcdf-c >= 4.9.0 will automatically
be available.
Otherwise, the user will have to set <code>HDF5_PLUGIN_PATH</code> explicitly
to have access to the extra compression plugins.</li>
<li>run <code>pip install -v .</code> (as root if necessary)</li>
<li>run the tests in the 'test' directory by running <code>python run_all.py</code>.</li>
</ul>
<h1 id="tutorial">Tutorial</h1>
<ul>
<li><a href="#creatingopeningclosing-a-netcdf-file">Creating/Opening/Closing a netCDF file</a></li>
<li><a href="#groups-in-a-netcdf-file">Groups in a netCDF file</a></li>
<li><a href="#dimensions-in-a-netcdf-file">Dimensions in a netCDF file</a></li>
<li><a href="#variables-in-a-netcdf-file">Variables in a netCDF file</a></li>
<li><a href="#attributes-in-a-netcdf-file">Attributes in a netCDF file</a></li>
<li><a href="#dealing-with-time-coordinates">Dealing with time coordinates</a></li>
<li><a href="#writing-data-to-and-retrieving-data-from-a-netcdf-variable">Writing data to and retrieving data from a netCDF variable</a></li>
<li><a href="#reading-data-from-a-multi-file-netcdf-dataset">Reading data from a multi-file netCDF dataset</a></li>
<li><a href="#efficient-compression-of-netcdf-variables">Efficient compression of netCDF variables</a></li>
<li><a href="#beyond-homogeneous-arrays-of-a-fixed-type-compound-data-types">Beyond homogeneous arrays of a fixed type - compound data types</a></li>
<li><a href="#variable-length-vlen-data-types">Variable-length (vlen) data types</a></li>
<li><a href="#enum-data-type">Enum data type</a></li>
<li><a href="#parallel-io">Parallel IO</a></li>
<li><a href="#dealing-with-strings">Dealing with strings</a></li>
<li><a href="#in-memory-diskless-datasets">In-memory (diskless) Datasets</a></li>
</ul>
<p>All of the code in this tutorial is available in <code>examples/tutorial.py</code>, except
the parallel IO example, which is in <code>examples/mpi_example.py</code>.
Unit tests are in the <code>test</code> directory.</p>
<h2 id="creatingopeningclosing-a-netcdf-file">Creating/Opening/Closing a netCDF file</h2>
<p>To create a netCDF file from python, you simply call the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>
constructor. This is also the method used to open an existing netCDF
file.
If the file is open for write access (<code>mode='w', 'r+'</code> or <code>'a'</code>), you may
write any type of data including new dimensions, groups, variables and
attributes.
netCDF files come in five flavors (<code>NETCDF3_CLASSIC</code>,
<code>NETCDF3_64BIT_OFFSET</code>, <code>NETCDF3_64BIT_DATA</code>, <code>NETCDF4_CLASSIC</code>, and <code>NETCDF4</code>).
<code>NETCDF3_CLASSIC</code> was the original netcdf binary format, and was limited
to file sizes less than 2 Gb. <code>NETCDF3_64BIT_OFFSET</code> was introduced
in version 3.6.0 of the library, and extended the original binary format
to allow for file sizes greater than 2 Gb.
<code>NETCDF3_64BIT_DATA</code> is a new format that requires version 4.4.0 of
the C library - it extends the <code>NETCDF3_64BIT_OFFSET</code> binary format to
allow for unsigned/64 bit integer data types and 64-bit dimension sizes.
<code>NETCDF3_64BIT</code> is an alias for <code>NETCDF3_64BIT_OFFSET</code>.
<code>NETCDF4_CLASSIC</code> files use the version 4 disk format (HDF5), but omits features
not found in the version 3 API. They can be read by netCDF 3 clients
only if they have been relinked against the netCDF 4 library. They can
also be read by HDF5 clients. <code>NETCDF4</code> files use the version 4 disk
format (HDF5) and use the new features of the version 4 API.
The
netCDF4 module can read and write files in any of these formats. When
creating a new file, the format may be specified using the <code>format</code>
keyword in the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> constructor.
The default format is
<code>NETCDF4</code>. To see how a given file is formatted, you can examine the
<code>data_model</code> attribute.
Closing the netCDF file is
accomplished via the <code><a title="netCDF4.Dataset.close" href="#netCDF4.Dataset.close">Dataset.close()</a></code> method of the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>
instance.</p>
<p>Here's an example:</p>
<pre><code class="language-python">>>> from netCDF4 import Dataset
>>> rootgrp = Dataset("test.nc", "w", format="NETCDF4")
>>> print(rootgrp.data_model)
NETCDF4
>>> rootgrp.close()
</code></pre>
<p>Remote <a href="http://opendap.org">OPeNDAP</a>-hosted datasets can be accessed for
reading over http if a URL is provided to the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> constructor instead of a
filename.
However, this requires that the netCDF library be built with
OPenDAP support, via the <code>--enable-dap</code> configure option (added in
version 4.0.1).</p>
<h2 id="groups-in-a-netcdf-file">Groups in a netCDF file</h2>
<p>netCDF version 4 added support for organizing data in hierarchical
groups, which are analogous to directories in a filesystem. Groups serve
as containers for variables, dimensions and attributes, as well as other
groups.
A <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> creates a special group, called
the 'root group', which is similar to the root directory in a unix
filesystem.
To create <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instances, use the
<code><a title="netCDF4.Dataset.createGroup" href="#netCDF4.Dataset.createGroup">Dataset.createGroup()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>
instance. <code><a title="netCDF4.Dataset.createGroup" href="#netCDF4.Dataset.createGroup">Dataset.createGroup()</a></code> takes a single argument, a
python string containing the name of the new group. The new <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>
instances contained within the root group can be accessed by name using
the <code>groups</code> dictionary attribute of the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance.
Only
<code>NETCDF4</code> formatted files support Groups, if you try to create a Group
in a netCDF 3 file you will get an error message.</p>
<pre><code class="language-python">>>> rootgrp = Dataset("test.nc", "a")
>>> fcstgrp = rootgrp.createGroup("forecasts")
>>> analgrp = rootgrp.createGroup("analyses")
>>> print(rootgrp.groups)
{'forecasts': <class 'netCDF4._netCDF4.Group'>
group /forecasts:
dimensions(sizes):
variables(dimensions):
groups: , 'analyses': <class 'netCDF4._netCDF4.Group'>
group /analyses:
dimensions(sizes):
variables(dimensions):
groups: }
>>>
</code></pre>
<p>Groups can exist within groups in a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>, just as directories
exist within directories in a unix filesystem. Each <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance
has a <code>groups</code> attribute dictionary containing all of the group
instances contained within that group. Each <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance also has a
<code>path</code> attribute that contains a simulated unix directory path to
that group.
To simplify the creation of nested groups, you can
use a unix-like path as an argument to <code><a title="netCDF4.Dataset.createGroup" href="#netCDF4.Dataset.createGroup">Dataset.createGroup()</a></code>.</p>
<pre><code class="language-python">>>> fcstgrp1 = rootgrp.createGroup("/forecasts/model1")
>>> fcstgrp2 = rootgrp.createGroup("/forecasts/model2")
</code></pre>
<p>If any of the intermediate elements of the path do not exist, they are created,
just as with the unix command <code>'mkdir -p'</code>. If you try to create a group
that already exists, no error will be raised, and the existing group will be
returned.</p>
<p>Here's an example that shows how to navigate all the groups in a
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>. The function <code>walktree</code> is a Python generator that is used
to walk the directory tree. Note that printing the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>
object yields summary information about it's contents.</p>
<pre><code class="language-python">>>> def walktree(top):
... yield top.groups.values()
... for value in top.groups.values():
... yield from walktree(value)
>>> print(rootgrp)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes):
variables(dimensions):
groups: forecasts, analyses
>>> for children in walktree(rootgrp):
... for child in children:
... print(child)
<class 'netCDF4._netCDF4.Group'>
group /forecasts:
dimensions(sizes):
variables(dimensions):
groups: model1, model2
<class 'netCDF4._netCDF4.Group'>
group /analyses:
dimensions(sizes):
variables(dimensions):
groups:
<class 'netCDF4._netCDF4.Group'>
group /forecasts/model1:
dimensions(sizes):
variables(dimensions):
groups:
<class 'netCDF4._netCDF4.Group'>
group /forecasts/model2:
dimensions(sizes):
variables(dimensions):
groups:
</code></pre>
<h2 id="dimensions-in-a-netcdf-file">Dimensions in a netCDF file</h2>
<p>netCDF defines the sizes of all variables in terms of dimensions, so
before any variables can be created the dimensions they use must be
created first. A special case, not often used in practice, is that of a
scalar variable, which has no dimensions. A dimension is created using
the <code><a title="netCDF4.Dataset.createDimension" href="#netCDF4.Dataset.createDimension">Dataset.createDimension()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>
or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance. A Python string is used to set the name of the
dimension, and an integer value is used to set the size. To create an
unlimited dimension (a dimension that can be appended to), the size
value is set to <code>None</code> or 0. In this example, there both the <code>time</code> and
<code>level</code> dimensions are unlimited.
Having more than one unlimited
dimension is a new netCDF 4 feature, in netCDF 3 files there may be only
one, and it must be the first (leftmost) dimension of the variable.</p>
<pre><code class="language-python">>>> level = rootgrp.createDimension("level", None)
>>> time = rootgrp.createDimension("time", None)
>>> lat = rootgrp.createDimension("lat", 73)
>>> lon = rootgrp.createDimension("lon", 144)
</code></pre>
<p>All of the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instances are stored in a python dictionary.</p>
<pre><code class="language-python">>>> print(rootgrp.dimensions)
{'level': <class 'netCDF4._netCDF4.Dimension'> (unlimited): name = 'level', size = 0, 'time': <class 'netCDF4._netCDF4.Dimension'> (unlimited): name = 'time', size = 0, 'lat': <class 'netCDF4._netCDF4.Dimension'>: name = 'lat', size = 73, 'lon': <class 'netCDF4._netCDF4.Dimension'>: name = 'lon', size = 144}
</code></pre>
<p>Using the python <code>len</code> function with a <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance returns
current size of that dimension.
<code><a title="netCDF4.Dimension.isunlimited" href="#netCDF4.Dimension.isunlimited">Dimension.isunlimited()</a></code> method of a <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance
be used to determine if the dimensions is unlimited, or appendable.</p>
<pre><code class="language-python">>>> print(len(lon))
144
>>> print(lon.isunlimited())
False
>>> print(time.isunlimited())
True
</code></pre>
<p>Printing the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> object
provides useful summary info, including the name and length of the dimension,
and whether it is unlimited.</p>
<pre><code class="language-python">>>> for dimobj in rootgrp.dimensions.values():
... print(dimobj)
<class 'netCDF4._netCDF4.Dimension'> (unlimited): name = 'level', size = 0
<class 'netCDF4._netCDF4.Dimension'> (unlimited): name = 'time', size = 0
<class 'netCDF4._netCDF4.Dimension'>: name = 'lat', size = 73
<class 'netCDF4._netCDF4.Dimension'>: name = 'lon', size = 144
</code></pre>
<p><code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> names can be changed using the
<code><a title="netCDF4.Dataset.renameDimension" href="#netCDF4.Dataset.renameDimension">Dataset.renameDimension()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance.</p>
<h2 id="variables-in-a-netcdf-file">Variables in a netCDF file</h2>
<p>netCDF variables behave much like python multidimensional array objects
supplied by the <a href="http://numpy.scipy.org">numpy module</a>. However,
unlike numpy arrays, netCDF4 variables can be appended to along one or
more 'unlimited' dimensions. To create a netCDF variable, use the
<code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance. The <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code>j method
has two mandatory arguments, the variable name (a Python string), and
the variable datatype. The variable's dimensions are given by a tuple
containing the dimension names (defined previously with
<code><a title="netCDF4.Dataset.createDimension" href="#netCDF4.Dataset.createDimension">Dataset.createDimension()</a></code>). To create a scalar
variable, simply leave out the dimensions keyword. The variable
primitive datatypes correspond to the dtype attribute of a numpy array.
You can specify the datatype as a numpy dtype object, or anything that
can be converted to a numpy dtype object.
Valid datatype specifiers
include: <code>'f4'</code> (32-bit floating point), <code>'f8'</code> (64-bit floating
point), <code>'i4'</code> (32-bit signed integer), <code>'i2'</code> (16-bit signed
integer), <code>'i8'</code> (64-bit signed integer), <code>'i1'</code> (8-bit signed
integer), <code>'u1'</code> (8-bit unsigned integer), <code>'u2'</code> (16-bit unsigned
integer), <code>'u4'</code> (32-bit unsigned integer), <code>'u8'</code> (64-bit unsigned
integer), or <code>'S1'</code> (single-character string).
The old Numeric
single-character typecodes (<code>'f'</code>,<code>'d'</code>,<code>'h'</code>,
<code>'s'</code>,<code>'b'</code>,<code>'B'</code>,<code>'c'</code>,<code>'i'</code>,<code>'l'</code>), corresponding to
(<code>'f4'</code>,<code>'f8'</code>,<code>'i2'</code>,<code>'i2'</code>,<code>'i1'</code>,<code>'i1'</code>,<code>'S1'</code>,<code>'i4'</code>,<code>'i4'</code>),
will also work. The unsigned integer types and the 64-bit integer type
can only be used if the file format is <code>NETCDF4</code>.</p>
<p>The dimensions themselves are usually also defined as variables, called
coordinate variables. The <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code>
method returns an instance of the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> class whose methods can be
used later to access and set variable data and attributes.</p>
<pre><code class="language-python">>>> times = rootgrp.createVariable("time","f8",("time",))
>>> levels = rootgrp.createVariable("level","i4",("level",))
>>> latitudes = rootgrp.createVariable("lat","f4",("lat",))
>>> longitudes = rootgrp.createVariable("lon","f4",("lon",))
>>> # two dimensions unlimited
>>> temp = rootgrp.createVariable("temp","f4",("time","level","lat","lon",))
>>> temp.units = "K"
</code></pre>
<p>To get summary info on a <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance in an interactive session,
just print it.</p>
<pre><code class="language-python">>>> print(temp)
<class 'netCDF4._netCDF4.Variable'>
float32 temp(time, level, lat, lon)
units: K
unlimited dimensions: time, level
current shape = (0, 0, 73, 144)
filling on, default _FillValue of 9.969209968386869e+36 used
</code></pre>
<p>You can use a path to create a Variable inside a hierarchy of groups.</p>
<pre><code class="language-python">>>> ftemp = rootgrp.createVariable("/forecasts/model1/temp","f4",("time","level","lat","lon",))
</code></pre>
<p>If the intermediate groups do not yet exist, they will be created.</p>
<p>You can also query a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance directly to obtain <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instances using paths.</p>
<pre><code class="language-python">>>> print(rootgrp["/forecasts/model1"]) # a Group instance
<class 'netCDF4._netCDF4.Group'>
group /forecasts/model1:
dimensions(sizes):
variables(dimensions): float32 temp(time,level,lat,lon)
groups:
>>> print(rootgrp["/forecasts/model1/temp"]) # a Variable instance
<class 'netCDF4._netCDF4.Variable'>
float32 temp(time, level, lat, lon)
path = /forecasts/model1
unlimited dimensions: time, level
current shape = (0, 0, 73, 144)
filling on, default _FillValue of 9.969209968386869e+36 used
</code></pre>
<p>All of the variables in the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> are stored in a
Python dictionary, in the same way as the dimensions:</p>
<pre><code class="language-python">>>> print(rootgrp.variables)
{'time': <class 'netCDF4._netCDF4.Variable'>
float64 time(time)
unlimited dimensions: time
current shape = (0,)
filling on, default _FillValue of 9.969209968386869e+36 used, 'level': <class 'netCDF4._netCDF4.Variable'>
int32 level(level)
unlimited dimensions: level
current shape = (0,)
filling on, default _FillValue of -2147483647 used, 'lat': <class 'netCDF4._netCDF4.Variable'>
float32 lat(lat)
unlimited dimensions:
current shape = (73,)
filling on, default _FillValue of 9.969209968386869e+36 used, 'lon': <class 'netCDF4._netCDF4.Variable'>
float32 lon(lon)
unlimited dimensions:
current shape = (144,)
filling on, default _FillValue of 9.969209968386869e+36 used, 'temp': <class 'netCDF4._netCDF4.Variable'>
float32 temp(time, level, lat, lon)
units: K
unlimited dimensions: time, level
current shape = (0, 0, 73, 144)
filling on, default _FillValue of 9.969209968386869e+36 used}
</code></pre>
<p><code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> names can be changed using the
<code><a title="netCDF4.Dataset.renameVariable" href="#netCDF4.Dataset.renameVariable">Dataset.renameVariable()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>
instance.</p>
<p>Variables can be sliced similar to numpy arrays, but there are some differences.
See
<a href="#writing-data-to-and-retrieving-data-from-a-netcdf-variable">Writing data to and retrieving data from a netCDF variable</a> for more details.</p>
<h2 id="attributes-in-a-netcdf-file">Attributes in a netCDF file</h2>
<p>There are two types of attributes in a netCDF file, global and variable.
Global attributes provide information about a group, or the entire
dataset, as a whole. <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> attributes provide information about
one of the variables in a group. Global attributes are set by assigning
values to <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance variables. <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>
attributes are set by assigning values to <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instances
variables. Attributes can be strings, numbers or sequences. Returning to
our example,</p>
<pre><code class="language-python">>>> import time
>>> rootgrp.description = "bogus example script"
>>> rootgrp.history = "Created " + time.ctime(time.time())
>>> rootgrp.source = "netCDF4 python module tutorial"
>>> latitudes.units = "degrees north"
>>> longitudes.units = "degrees east"
>>> levels.units = "hPa"
>>> temp.units = "K"
>>> times.units = "hours since 0001-01-01 00:00:00.0"
>>> times.calendar = "gregorian"
</code></pre>
<p>The <code><a title="netCDF4.Dataset.ncattrs" href="#netCDF4.Dataset.ncattrs">Dataset.ncattrs()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>, <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance can be used to retrieve the names of all the netCDF
attributes. This method is provided as a convenience, since using the
built-in <code>dir</code> Python function will return a bunch of private methods
and attributes that cannot (or should not) be modified by the user.</p>
<pre><code class="language-python">>>> for name in rootgrp.ncattrs():
... print("Global attr {} = {}".format(name, getattr(rootgrp, name)))
Global attr description = bogus example script
Global attr history = Created Mon Jul 8 14:19:41 2019
Global attr source = netCDF4 python module tutorial
</code></pre>
<p>The <code>__dict__</code> attribute of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>, <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>
instance provides all the netCDF attribute name/value pairs in a python
dictionary:</p>
<pre><code class="language-python">>>> print(rootgrp.__dict__)
{'description': 'bogus example script', 'history': 'Created Mon Jul 8 14:19:41 2019', 'source': 'netCDF4 python module tutorial'}
</code></pre>
<p>Attributes can be deleted from a netCDF <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>, <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> using the python <code>del</code> statement (i.e. <code>del grp.foo</code>
removes the attribute <code>foo</code> the the group <code>grp</code>).</p>
<h2 id="writing-data-to-and-retrieving-data-from-a-netcdf-variable">Writing data to and retrieving data from a netCDF variable</h2>
<p>Now that you have a netCDF <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance, how do you put data
into it? You can just treat it like an array and assign data to a slice.</p>
<pre><code class="language-python">>>> import numpy as np
>>> lats = np.arange(-90,91,2.5)
>>> lons = np.arange(-180,180,2.5)
>>> latitudes[:] = lats
>>> longitudes[:] = lons
>>> print("latitudes =\n{}".format(latitudes[:]))
latitudes =
[-90. -87.5 -85. -82.5 -80. -77.5 -75. -72.5 -70. -67.5 -65. -62.5
-60. -57.5 -55. -52.5 -50. -47.5 -45. -42.5 -40. -37.5 -35. -32.5
-30. -27.5 -25. -22.5 -20. -17.5 -15. -12.5 -10. -7.5 -5. -2.5
0. 2.5 5. 7.5 10. 12.5 15. 17.5 20. 22.5 25. 27.5
30. 32.5 35. 37.5 40. 42.5 45. 47.5 50. 52.5 55. 57.5
60. 62.5 65. 67.5 70. 72.5 75. 77.5 80. 82.5 85. 87.5
90. ]
</code></pre>
<p>Unlike NumPy's array objects, netCDF <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>
objects with unlimited dimensions will grow along those dimensions if you
assign data outside the currently defined range of indices.</p>
<pre><code class="language-python">>>> # append along two unlimited dimensions by assigning to slice.
>>> nlats = len(rootgrp.dimensions["lat"])
>>> nlons = len(rootgrp.dimensions["lon"])
>>> print("temp shape before adding data = {}".format(temp.shape))
temp shape before adding data = (0, 0, 73, 144)
>>>
>>> from numpy.random import uniform
>>> temp[0:5, 0:10, :, :] = uniform(size=(5, 10, nlats, nlons))
>>> print("temp shape after adding data = {}".format(temp.shape))
temp shape after adding data = (5, 10, 73, 144)
>>>
>>> # levels have grown, but no values yet assigned.
>>> print("levels shape after adding pressure data = {}".format(levels.shape))
levels shape after adding pressure data = (10,)
</code></pre>
<p>Note that the size of the levels variable grows when data is appended
along the <code>level</code> dimension of the variable <code>temp</code>, even though no
data has yet been assigned to levels.</p>
<pre><code class="language-python">>>> # now, assign data to levels dimension variable.
>>> levels[:] = [1000.,850.,700.,500.,300.,250.,200.,150.,100.,50.]
</code></pre>
<p>However, that there are some differences between NumPy and netCDF
variable slicing rules. Slices behave as usual, being specified as a
<code>start:stop:step</code> triplet. Using a scalar integer index <code>i</code> takes the ith
element and reduces the rank of the output array by one. Boolean array and
integer sequence indexing behaves differently for netCDF variables
than for numpy arrays.
Only 1-d boolean arrays and integer sequences are
allowed, and these indices work independently along each dimension (similar
to the way vector subscripts work in fortran).
This means that</p>
<pre><code class="language-python">>>> temp[0, 0, [0,1,2,3], [0,1,2,3]].shape
(4, 4)
</code></pre>
<p>returns an array of shape (4,4) when slicing a netCDF variable, but for a
numpy array it returns an array of shape (4,).
Similarly, a netCDF variable of shape <code>(2,3,4,5)</code> indexed
with <code>[0, array([True, False, True]), array([False, True, True, True]), :]</code>
would return a <code>(2, 3, 5)</code> array. In NumPy, this would raise an error since
it would be equivalent to <code>[0, [0,1], [1,2,3], :]</code>. When slicing with integer
sequences, the indices <strong><em>need not be sorted</em></strong> and <strong><em>may contain
duplicates</em></strong> (both of these are new features in version 1.2.1).
While this behaviour may cause some confusion for those used to NumPy's 'fancy indexing' rules,
it provides a very powerful way to extract data from multidimensional netCDF
variables by using logical operations on the dimension arrays to create slices.</p>
<p>For example,</p>
<pre><code class="language-python">>>> tempdat = temp[::2, [1,3,6], lats>0, lons>0]
</code></pre>
<p>will extract time indices 0,2 and 4, pressure levels
850, 500 and 200 hPa, all Northern Hemisphere latitudes and Eastern
Hemisphere longitudes, resulting in a numpy array of shape
(3, 3, 36, 71).</p>
<pre><code class="language-python">>>> print("shape of fancy temp slice = {}".format(tempdat.shape))
shape of fancy temp slice = (3, 3, 36, 71)
</code></pre>
<p><strong><em>Special note for scalar variables</em></strong>: To extract data from a scalar variable
<code>v</code> with no associated dimensions, use <code>numpy.asarray(v)</code> or <code>v[…]</code>.
The result will be a numpy scalar array.</p>
<p>By default, netcdf4-python returns numpy masked arrays with values equal to the
<code>missing_value</code> or <code>_FillValue</code> variable attributes masked for primitive and
enum data types.
The <code><a title="netCDF4.Dataset.set_auto_mask" href="#netCDF4.Dataset.set_auto_mask">Dataset.set_auto_mask()</a></code> <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> and <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> methods
can be used to disable this feature so that
numpy arrays are always returned, with the missing values included. Prior to
version 1.4.0 the default behavior was to only return masked arrays when the
requested slice contained missing values.
This behavior can be recovered
using the <code><a title="netCDF4.Dataset.set_always_mask" href="#netCDF4.Dataset.set_always_mask">Dataset.set_always_mask()</a></code> method. If a masked array is
written to a netCDF variable, the masked elements are filled with the
value specified by the <code>missing_value</code> attribute.
If the variable has
no <code>missing_value</code>, the <code>_FillValue</code> is used instead.</p>
<h2 id="dealing-with-time-coordinates">Dealing with time coordinates</h2>
<p>Time coordinate values pose a special challenge to netCDF users.
Most
metadata standards (such as CF) specify that time should be
measure relative to a fixed date using a certain calendar, with units
specified like <code>hours since YY-MM-DD hh:mm:ss</code>.
These units can be
awkward to deal with, without a utility to convert the values to and
from calendar dates.
The functions <a href="https://unidata.github.io/cftime/api.html">num2date</a>
and <a href="https://unidata.github.io/cftime/api.html">date2num</a> are
provided by <a href="https://unidata.github.io/cftime">cftime</a> to do just that.
Here's an example of how they can be used:</p>
<pre><code class="language-python">>>> # fill in times.
>>> from datetime import datetime, timedelta
>>> from cftime import num2date, date2num
>>> dates = [datetime(2001,3,1)+n*timedelta(hours=12) for n in range(temp.shape[0])]
>>> times[:] = date2num(dates,units=times.units,calendar=times.calendar)
>>> print("time values (in units {}):\n{}".format(times.units, times[:]))
time values (in units hours since 0001-01-01 00:00:00.0):
[17533104. 17533116. 17533128. 17533140. 17533152.]
>>> dates = num2date(times[:],units=times.units,calendar=times.calendar)
>>> print("dates corresponding to time values:\n{}".format(dates))
[cftime.DatetimeGregorian(2001, 3, 1, 0, 0, 0, 0, has_year_zero=False)
cftime.DatetimeGregorian(2001, 3, 1, 12, 0, 0, 0, has_year_zero=False)
cftime.DatetimeGregorian(2001, 3, 2, 0, 0, 0, 0, has_year_zero=False)
cftime.DatetimeGregorian(2001, 3, 2, 12, 0, 0, 0, has_year_zero=False)
cftime.DatetimeGregorian(2001, 3, 3, 0, 0, 0, 0, has_year_zero=False)]
</code></pre>
<p><code><a title="netCDF4.num2date" href="#netCDF4.num2date">num2date()</a></code> converts numeric values of time in the specified <code>units</code>
and <code>calendar</code> to datetime objects, and <code><a title="netCDF4.date2num" href="#netCDF4.date2num">date2num()</a></code> does the reverse.
All the calendars currently defined in the
<a href="http://cfconventions.org">CF metadata convention</a> are supported.
A function called <code><a title="netCDF4.date2index" href="#netCDF4.date2index">date2index()</a></code> is also provided which returns the indices
of a netCDF time variable corresponding to a sequence of datetime instances.</p>
<h2 id="reading-data-from-a-multi-file-netcdf-dataset">Reading data from a multi-file netCDF dataset</h2>
<p>If you want to read data from a variable that spans multiple netCDF files,
you can use the <code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code> class to read the data as if it were
contained in a single file. Instead of using a single filename to create
a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance, create a <code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code> instance with either a list
of filenames, or a string with a wildcard (which is then converted to
a sorted list of files using the python glob module).
Variables in the list of files that share the same unlimited
dimension are aggregated together, and can be sliced across multiple
files.
To illustrate this, let's first create a bunch of netCDF files with
the same variable (with the same unlimited dimension).
The files
must in be in <code>NETCDF3_64BIT_OFFSET</code>, <code>NETCDF3_64BIT_DATA</code>, <code>NETCDF3_CLASSIC</code> or
<code>NETCDF4_CLASSIC</code> format (<code>NETCDF4</code> formatted multi-file
datasets are not supported).</p>
<pre><code class="language-python">>>> for nf in range(10):
... with Dataset("mftest%s.nc" % nf, "w", format="NETCDF4_CLASSIC") as f:
... _ = f.createDimension("x",None)
... x = f.createVariable("x","i",("x",))
... x[0:10] = np.arange(nf*10,10*(nf+1))
</code></pre>
<p>Now read all the files back in at once with <code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code></p>
<pre><code class="language-python">>>> from netCDF4 import MFDataset
>>> f = MFDataset("mftest*nc")
>>> print(f.variables["x"][:])
[ 0 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]
</code></pre>
<p>Note that <code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code> can only be used to read, not write, multi-file
datasets.</p>
<h2 id="efficient-compression-of-netcdf-variables">Efficient compression of netCDF variables</h2>
<p>Data stored in netCDF <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> objects can be compressed and
decompressed on the fly. The compression algorithm used is determined
by the <code>compression</code> keyword argument to the <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method.
<code>zlib</code> compression is always available, <code>szip</code> is available if the linked HDF5
library supports it, and <code>zstd</code>, <code>bzip2</code>, <code>blosc_lz</code>,<code>blosc_lz4</code>,<code>blosc_lz4hc</code>,
<code>blosc_zlib</code> and <code>blosc_zstd</code> are available via optional external plugins.
The <code>complevel</code> keyword regulates the
speed and efficiency of the compression for <code>zlib</code>, <code>bzip</code> and <code>zstd</code> (1 being fastest, but lowest
compression ratio, 9 being slowest but best compression ratio). The
default value of <code>complevel</code> is 4. Setting <code>shuffle=False</code> will turn
off the HDF5 shuffle filter, which de-interlaces a block of data before
<code>zlib</code> compression by reordering the bytes.
The shuffle filter can
significantly improve compression ratios, and is on by default if <code>compression=zlib</code>.
Setting
<code>fletcher32</code> keyword argument to
<code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> to <code>True</code> (it's <code>False</code> by
default) enables the Fletcher32 checksum algorithm for error detection.
It's also possible to set the HDF5 chunking parameters and endian-ness
of the binary data stored in the HDF5 file with the <code>chunksizes</code>
and <code>endian</code> keyword arguments to
<code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code>.
These keyword arguments only
are relevant for <code>NETCDF4</code> and <code>NETCDF4_CLASSIC</code> files (where the
underlying file format is HDF5) and are silently ignored if the file
format is <code>NETCDF3_CLASSIC</code>, <code>NETCDF3_64BIT_OFFSET</code> or <code>NETCDF3_64BIT_DATA</code>.
If the HDF5 library is built with szip support, compression=<code>szip</code> can also
be used (in conjunction with the <code>szip_coding</code> and <code>szip_pixels_per_block</code> keyword
arguments).</p>
<p>If your data only has a certain number of digits of precision (say for
example, it is temperature data that was measured with a precision of
0.1 degrees), you can dramatically improve compression by
quantizing (or truncating) the data. There are two methods supplied for
doing this.
You can use the <code>least_significant_digit</code>
keyword argument to <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> to specify
the power of ten of the smallest decimal place in
the data that is a reliable value. For example if the data has a
precision of 0.1, then setting <code>least_significant_digit=1</code> will cause
data the data to be quantized using <code>numpy.around(scale*data)/scale</code>, where
scale = 2**bits, and bits is determined so that a precision of 0.1 is
retained (in this case bits=4).
This is done at the python level and is
not a part of the underlying C library.
Starting with netcdf-c version 4.9.0,
a quantization capability is provided in the library.
This can be
used via the <code>significant_digits</code> <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> kwarg (new in
version 1.6.0).
The interpretation of <code>significant_digits</code> is different than <code>least_signficant_digit</code>
in that it specifies the absolute number of significant digits independent
of the magnitude of the variable (the floating point exponent).
Either of these approaches makes the compression
'lossy' instead of 'lossless', that is some precision in the data is
sacrificed for the sake of disk space.</p>
<p>In our example, try replacing the line</p>
<pre><code class="language-python">>>> temp = rootgrp.createVariable("temp","f4",("time","level","lat","lon",))
</code></pre>
<p>with</p>
<pre><code class="language-python">>>> temp = rootgrp.createVariable("temp","f4",("time","level","lat","lon",),compression='zlib')
</code></pre>
<p>and then</p>
<pre><code class="language-python">>>> temp = rootgrp.createVariable("temp","f4",("time","level","lat","lon",),compression='zlib',least_significant_digit=3)
</code></pre>
<p>or with netcdf-c >= 4.9.0</p>
<pre><code class="language-python">>>> temp = rootgrp.createVariable("temp","f4",("time","level","lat","lon",),compression='zlib',significant_digits=4)
</code></pre>
<p>and see how much smaller the resulting files are.</p>
<h2 id="beyond-homogeneous-arrays-of-a-fixed-type-compound-data-types">Beyond homogeneous arrays of a fixed type - compound data types</h2>
<p>Compound data types map directly to numpy structured (a.k.a 'record')
arrays.
Structured arrays are akin to C structs, or derived types
in Fortran. They allow for the construction of table-like structures
composed of combinations of other data types, including other
compound types. Compound types might be useful for representing multiple
parameter values at each point on a grid, or at each time and space
location for scattered (point) data. You can then access all the
information for a point by reading one variable, instead of reading
different parameters from different variables.
Compound data types
are created from the corresponding numpy data type using the
<code><a title="netCDF4.Dataset.createCompoundType" href="#netCDF4.Dataset.createCompoundType">Dataset.createCompoundType()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance.
Since there is no native complex data type in netcdf (but see
<a href="#support-for-complex-numbers">Support for complex numbers</a>), compound
types are handy for storing numpy complex arrays. Here's an example:</p>
<pre><code class="language-python">>>> f = Dataset("complex.nc","w")
>>> size = 3 # length of 1-d complex array
>>> # create sample complex data.
>>> datac = np.exp(1j*(1.+np.linspace(0, np.pi, size)))
>>> # create complex128 compound data type.
>>> complex128 = np.dtype([("real",np.float64),("imag",np.float64)])
>>> complex128_t = f.createCompoundType(complex128,"complex128")
>>> # create a variable with this data type, write some data to it.
>>> x_dim = f.createDimension("x_dim",None)
>>> v = f.createVariable("cmplx_var",complex128_t,"x_dim")
>>> data = np.empty(size,complex128) # numpy structured array
>>> data["real"] = datac.real; data["imag"] = datac.imag
>>> v[:] = data # write numpy structured array to netcdf compound var
>>> # close and reopen the file, check the contents.
>>> f.close(); f = Dataset("complex.nc")
>>> v = f.variables["cmplx_var"]
>>> datain = v[:] # read in all the data into a numpy structured array
>>> # create an empty numpy complex array
>>> datac2 = np.empty(datain.shape,np.complex128)
>>> # .. fill it with contents of structured array.
>>> datac2.real = datain["real"]; datac2.imag = datain["imag"]
>>> print('{}: {}'.format(datac.dtype, datac)) # original data
complex128: [ 0.54030231+0.84147098j -0.84147098+0.54030231j -0.54030231-0.84147098j]
>>>
>>> print('{}: {}'.format(datac2.dtype, datac2)) # data from file
complex128: [ 0.54030231+0.84147098j -0.84147098+0.54030231j -0.54030231-0.84147098j]
</code></pre>
<p>Compound types can be nested, but you must create the 'inner'
ones first. All possible numpy structured arrays cannot be
represented as Compound variables - an error message will be
raise if you try to create one that is not supported.
All of the compound types defined for a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> are stored
in a Python dictionary, just like variables and dimensions. As always, printing
objects gives useful summary information in an interactive session:</p>
<pre><code class="language-python">>>> print(f)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes): x_dim(3)
variables(dimensions): {'names':['real','imag'], 'formats':['<f8','<f8'], 'offsets':[0,8], 'itemsize':16, 'aligned':True} cmplx_var(x_dim)
groups:
>>> print(f.variables["cmplx_var"])
<class 'netCDF4._netCDF4.Variable'>
compound cmplx_var(x_dim)
compound data type: {'names':['real','imag'], 'formats':['<f8','<f8'], 'offsets':[0,8], 'itemsize':16, 'aligned':True}
unlimited dimensions: x_dim
current shape = (3,)
>>> print(f.cmptypes)
{'complex128': <class 'netCDF4._netCDF4.CompoundType'>: name = 'complex128', numpy dtype = {'names':['real','imag'], 'formats':['<f8','<f8'], 'offsets':[0,8], 'itemsize':16, 'aligned':True}}
>>> print(f.cmptypes["complex128"])
<class 'netCDF4._netCDF4.CompoundType'>: name = 'complex128', numpy dtype = {'names':['real','imag'], 'formats':['<f8','<f8'], 'offsets':[0,8], 'itemsize':16, 'aligned':True}
</code></pre>
<h2 id="variable-length-vlen-data-types">Variable-length (vlen) data types</h2>
<p>NetCDF 4 has support for variable-length or "ragged" arrays.
These are arrays
of variable length sequences having the same type. To create a variable-length
data type, use the <code><a title="netCDF4.Dataset.createVLType" href="#netCDF4.Dataset.createVLType">Dataset.createVLType()</a></code> method
method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance.</p>
<pre><code class="language-python">>>> f = Dataset("tst_vlen.nc","w")
>>> vlen_t = f.createVLType(np.int32, "phony_vlen")
</code></pre>
<p>The numpy datatype of the variable-length sequences and the name of the
new datatype must be specified. Any of the primitive datatypes can be
used (signed and unsigned integers, 32 and 64 bit floats, and characters),
but compound data types cannot.
A new variable can then be created using this datatype.</p>
<pre><code class="language-python">>>> x = f.createDimension("x",3)
>>> y = f.createDimension("y",4)
>>> vlvar = f.createVariable("phony_vlen_var", vlen_t, ("y","x"))
</code></pre>
<p>Since there is no native vlen datatype in numpy, vlen arrays are represented
in python as object arrays (arrays of dtype <code>object</code>). These are arrays whose
elements are Python object pointers, and can contain any type of python object.
For this application, they must contain 1-D numpy arrays all of the same type
but of varying length.
In this case, they contain 1-D numpy <code>int32</code> arrays of random length between
1 and 10.</p>
<pre><code class="language-python">>>> import random
>>> random.seed(54321)
>>> data = np.empty(len(y)*len(x),object)
>>> for n in range(len(y)*len(x)):
... data[n] = np.arange(random.randint(1,10),dtype="int32")+1
>>> data = np.reshape(data,(len(y),len(x)))
>>> vlvar[:] = data
>>> print("vlen variable =\n{}".format(vlvar[:]))
vlen variable =
[[array([1, 2, 3, 4, 5, 6, 7, 8], dtype=int32) array([1, 2], dtype=int32)
array([1, 2, 3, 4], dtype=int32)]
[array([1, 2, 3], dtype=int32)
array([1, 2, 3, 4, 5, 6, 7, 8, 9], dtype=int32)
array([1, 2, 3, 4, 5, 6, 7, 8, 9], dtype=int32)]
[array([1, 2, 3, 4, 5, 6, 7], dtype=int32) array([1, 2, 3], dtype=int32)
array([1, 2, 3, 4, 5, 6], dtype=int32)]
[array([1, 2, 3, 4, 5, 6, 7, 8, 9], dtype=int32)
array([1, 2, 3, 4, 5], dtype=int32) array([1, 2], dtype=int32)]]
>>> print(f)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes): x(3), y(4)
variables(dimensions): int32 phony_vlen_var(y,x)
groups:
>>> print(f.variables["phony_vlen_var"])
<class 'netCDF4._netCDF4.Variable'>
vlen phony_vlen_var(y, x)
vlen data type: int32
unlimited dimensions:
current shape = (4, 3)
>>> print(f.vltypes["phony_vlen"])
<class 'netCDF4._netCDF4.VLType'>: name = 'phony_vlen', numpy dtype = int32
</code></pre>
<p>Numpy object arrays containing python strings can also be written as vlen
variables,
For vlen strings, you don't need to create a vlen data type.
Instead, simply use the python <code>str</code> builtin (or a numpy string datatype
with fixed length greater than 1) when calling the
<code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method.</p>
<pre><code class="language-python">>>> z = f.createDimension("z",10)
>>> strvar = f.createVariable("strvar", str, "z")
</code></pre>
<p>In this example, an object array is filled with random python strings with
random lengths between 2 and 12 characters, and the data in the object
array is assigned to the vlen string variable.</p>
<pre><code class="language-python">>>> chars = "1234567890aabcdefghijklmnopqrstuvwxyzABCDEFGHIJKLMNOPQRSTUVWXYZ"
>>> data = np.empty(10,"O")
>>> for n in range(10):
... stringlen = random.randint(2,12)
... data[n] = "".join([random.choice(chars) for i in range(stringlen)])
>>> strvar[:] = data
>>> print("variable-length string variable:\n{}".format(strvar[:]))
variable-length string variable:
['Lh' '25F8wBbMI' '53rmM' 'vvjnb3t63ao' 'qjRBQk6w' 'aJh' 'QF'
'jtIJbJACaQk4' '3Z5' 'bftIIq']
>>> print(f)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes): x(3), y(4), z(10)
variables(dimensions): int32 phony_vlen_var(y,x), <class 'str'> strvar(z)
groups:
>>> print(f.variables["strvar"])
<class 'netCDF4._netCDF4.Variable'>
vlen strvar(z)
vlen data type: <class 'str'>
unlimited dimensions:
current shape = (10,)
</code></pre>
<p>It is also possible to set contents of vlen string variables with numpy arrays
of any string or unicode data type. Note, however, that accessing the contents
of such variables will always return numpy arrays with dtype <code>object</code>.</p>
<h2 id="enum-data-type">Enum data type</h2>
<p>netCDF4 has an enumerated data type, which is an integer datatype that is
restricted to certain named values. Since Enums don't map directly to
a numpy data type, they are read and written as integer arrays.</p>
<p>Here's an example of using an Enum type to hold cloud type data.
The base integer data type and a python dictionary describing the allowed
values and their names are used to define an Enum data type using
<code><a title="netCDF4.Dataset.createEnumType" href="#netCDF4.Dataset.createEnumType">Dataset.createEnumType()</a></code>.</p>
<pre><code class="language-python">>>> nc = Dataset('clouds.nc','w')
>>> # python dict with allowed values and their names.
>>> enum_dict = {'Altocumulus': 7, 'Missing': 255,
... 'Stratus': 2, 'Clear': 0,
... 'Nimbostratus': 6, 'Cumulus': 4, 'Altostratus': 5,
... 'Cumulonimbus': 1, 'Stratocumulus': 3}
>>> # create the Enum type called 'cloud_t'.
>>> cloud_type = nc.createEnumType(np.uint8,'cloud_t',enum_dict)
>>> print(cloud_type)
<class 'netCDF4._netCDF4.EnumType'>: name = 'cloud_t', numpy dtype = uint8, fields/values ={'Altocumulus': 7, 'Missing': 255, 'Stratus': 2, 'Clear': 0, 'Nimbostratus': 6, 'Cumulus': 4, 'Altostratus': 5, 'Cumulonimbus': 1, 'Stratocumulus': 3}
</code></pre>
<p>A new variable can be created in the usual way using this data type.
Integer data is written to the variable that represents the named
cloud types in enum_dict. A <code>ValueError</code> will be raised if an attempt
is made to write an integer value not associated with one of the
specified names.</p>
<pre><code class="language-python">>>> time = nc.createDimension('time',None)
>>> # create a 1d variable of type 'cloud_type'.
>>> # The fill_value is set to the 'Missing' named value.
>>> cloud_var = nc.createVariable('primary_cloud',cloud_type,'time',
... fill_value=enum_dict['Missing'])
>>> # write some data to the variable.
>>> cloud_var[:] = [enum_dict[k] for k in ['Clear', 'Stratus', 'Cumulus',
... 'Missing', 'Cumulonimbus']]
>>> nc.close()
>>> # reopen the file, read the data.
>>> nc = Dataset('clouds.nc')
>>> cloud_var = nc.variables['primary_cloud']
>>> print(cloud_var)
<class 'netCDF4._netCDF4.Variable'>
enum primary_cloud(time)
_FillValue: 255
enum data type: uint8
unlimited dimensions: time
current shape = (5,)
>>> print(cloud_var.datatype.enum_dict)
{'Altocumulus': 7, 'Missing': 255, 'Stratus': 2, 'Clear': 0, 'Nimbostratus': 6, 'Cumulus': 4, 'Altostratus': 5, 'Cumulonimbus': 1, 'Stratocumulus': 3}
>>> print(cloud_var[:])
[0 2 4 -- 1]
>>> nc.close()
</code></pre>
<h2 id="parallel-io">Parallel IO</h2>
<p>If MPI parallel enabled versions of netcdf and hdf5 or pnetcdf are detected,
and <a href="https://mpi4py.scipy.org">mpi4py</a> is installed, netcdf4-python will
be built with parallel IO capabilities enabled. Parallel IO of NETCDF4 or
NETCDF4_CLASSIC formatted files is only available if the MPI parallel HDF5
library is available. Parallel IO of classic netcdf-3 file formats is only
available if the <a href="https://parallel-netcdf.github.io/">PnetCDF</a> library is
available. To use parallel IO, your program must be running in an MPI
environment using <a href="https://mpi4py.scipy.org">mpi4py</a>.</p>
<pre><code class="language-python">>>> from mpi4py import MPI
>>> import numpy as np
>>> from netCDF4 import Dataset
>>> rank = MPI.COMM_WORLD.rank # The process ID (integer 0-3 for 4-process run)
</code></pre>
<p>To run an MPI-based parallel program like this, you must use <code>mpiexec</code> to launch several
parallel instances of Python (for example, using <code>mpiexec -np 4 python mpi_example.py</code>).
The parallel features of netcdf4-python are mostly transparent -
when a new dataset is created or an existing dataset is opened,
use the <code>parallel</code> keyword to enable parallel access.</p>
<pre><code class="language-python">>>> nc = Dataset('parallel_test.nc','w',parallel=True)
</code></pre>
<p>The optional <code>comm</code> keyword may be used to specify a particular
MPI communicator (<code>MPI_COMM_WORLD</code> is used by default).
Each process (or rank)
can now write to the file independently.
In this example the process rank is
written to a different variable index on each task</p>
<pre><code class="language-python">>>> d = nc.createDimension('dim',4)
>>> v = nc.createVariable('var', np.int64, 'dim')
>>> v[rank] = rank
>>> nc.close()
% ncdump parallel_test.nc
netcdf parallel_test {
dimensions:
dim = 4 ;
variables:
int64 var(dim) ;
data:
var = 0, 1, 2, 3 ;
}
</code></pre>
<p>There are two types of parallel IO, independent (the default) and collective.
Independent IO means that each process can do IO independently. It should not
depend on or be affected by other processes. Collective IO is a way of doing
IO defined in the MPI-IO standard; unlike independent IO, all processes must
participate in doing IO. To toggle back and forth between
the two types of IO, use the <code><a title="netCDF4.Variable.set_collective" href="#netCDF4.Variable.set_collective">Variable.set_collective()</a></code>
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> method. All metadata
operations (such as creation of groups, types, variables, dimensions, or attributes)
are collective.
There are a couple of important limitations of parallel IO:</p>
<ul>
<li>parallel IO for NETCDF4 or NETCDF4_CLASSIC formatted files is only available
if the netcdf library was compiled with MPI enabled HDF5.</li>
<li>parallel IO for all classic netcdf-3 file formats is only available if the
netcdf library was compiled with <a href="https://parallel-netcdf.github.io">PnetCDF</a>.</li>
<li>If a variable has an unlimited dimension, appending data must be done in collective mode.
If the write is done in independent mode, the operation will fail with a
a generic "HDF Error".</li>
<li>You can write compressed data in parallel only with netcdf-c >= 4.7.4
and hdf5 >= 1.10.3 (although you can read in parallel with earlier versions). To write
compressed data in parallel, the variable must be in 'collective IO mode'.
This is done
automatically on variable creation if compression is turned on, but if you are appending
to a variable in an existing file, you must use <code><a title="netCDF4.Variable.set_collective" href="#netCDF4.Variable.set_collective">Variable.set_collective()</a>(True)</code> before attempting
to write to it.</li>
<li>You cannot use variable-length (VLEN) data types.</li>
</ul>
<h2 id="dealing-with-strings">Dealing with strings</h2>
<p>The most flexible way to store arrays of strings is with the
<a href="#variable-length-vlen-data-type">Variable-length (vlen) string data type</a>. However, this requires
the use of the NETCDF4 data model, and the vlen type does not map very well
numpy arrays (you have to use numpy arrays of dtype=<code>object</code>, which are arrays of
arbitrary python objects). numpy does have a fixed-width string array
data type, but unfortunately the netCDF data model does not.
Instead fixed-width byte strings are typically stored as <a href="https://www.unidata.ucar.edu/software/netcdf/docs/BestPractices.html#bp_Strings-and-Variables-of-type-char">arrays of 8-bit
characters</a>.
To perform the conversion to and from character arrays to fixed-width numpy string arrays, the
following convention is followed by the python interface.
If the <code>_Encoding</code> special attribute is set for a character array
(dtype <code>S1</code>) variable, the <code><a title="netCDF4.chartostring" href="#netCDF4.chartostring">chartostring()</a></code> utility function is used to convert the array of
characters to an array of strings with one less dimension (the last dimension is
interpreted as the length of each string) when reading the data. The character
set (usually ascii) is specified by the <code>_Encoding</code> attribute. If <code>_Encoding</code>
is 'none' or 'bytes', then the character array is converted to a numpy
fixed-width byte string array (dtype <code>S#</code>), otherwise a numpy unicode (dtype
<code>U#</code>) array is created.
When writing the data,
<code><a title="netCDF4.stringtochar" href="#netCDF4.stringtochar">stringtochar()</a></code> is used to convert the numpy string array to an array of
characters with one more dimension. For example,</p>
<pre><code class="language-python">>>> from netCDF4 import stringtochar
>>> nc = Dataset('stringtest.nc','w',format='NETCDF4_CLASSIC')
>>> _ = nc.createDimension('nchars',3)
>>> _ = nc.createDimension('nstrings',None)
>>> v = nc.createVariable('strings','S1',('nstrings','nchars'))
>>> datain = np.array(['foo','bar'],dtype='S3')
>>> v[:] = stringtochar(datain) # manual conversion to char array
>>> print(v[:]) # data returned as char array
[[b'f' b'o' b'o']
[b'b' b'a' b'r']]
>>> v._Encoding = 'ascii' # this enables automatic conversion
>>> v[:] = datain # conversion to char array done internally
>>> print(v[:]) # data returned in numpy string array
['foo' 'bar']
>>> nc.close()
</code></pre>
<p>Even if the <code>_Encoding</code> attribute is set, the automatic conversion of char
arrays to/from string arrays can be disabled with
<code><a title="netCDF4.Variable.set_auto_chartostring" href="#netCDF4.Variable.set_auto_chartostring">Variable.set_auto_chartostring()</a></code>.</p>
<p>A similar situation is often encountered with numpy structured arrays with subdtypes
containing fixed-wdith byte strings (dtype=<code>S#</code>). Since there is no native fixed-length string
netCDF datatype, these numpy structure arrays are mapped onto netCDF compound
types with character array elements.
In this case the string <-> char array
conversion is handled automatically (without the need to set the <code>_Encoding</code>
attribute) using <a href="https://docs.scipy.org/doc/numpy/reference/generated/numpy.ndarray.view.html">numpy
views</a>.
The structured array dtype (including the string elements) can even be used to
define the compound data type - the string dtype will be converted to
character array dtype under the hood when creating the netcdf compound type.
Here's an example:</p>
<pre><code class="language-python">>>> nc = Dataset('compoundstring_example.nc','w')
>>> dtype = np.dtype([('observation', 'f4'),
... ('station_name','S10')])
>>> station_data_t = nc.createCompoundType(dtype,'station_data')
>>> _ = nc.createDimension('station',None)
>>> statdat = nc.createVariable('station_obs', station_data_t, ('station',))
>>> data = np.empty(2,dtype)
>>> data['observation'][:] = (123.,3.14)
>>> data['station_name'][:] = ('Boulder','New York')
>>> print(statdat.dtype) # strings actually stored as character arrays
{'names':['observation','station_name'], 'formats':['<f4',('S1', (10,))], 'offsets':[0,4], 'itemsize':16, 'aligned':True}
>>> statdat[:] = data # strings converted to character arrays internally
>>> print(statdat[:]) # character arrays converted back to strings
[(123. , b'Boulder') ( 3.14, b'New York')]
>>> print(statdat[:].dtype)
{'names':['observation','station_name'], 'formats':['<f4','S10'], 'offsets':[0,4], 'itemsize':16, 'aligned':True}
>>> statdat.set_auto_chartostring(False) # turn off auto-conversion
>>> statdat[:] = data.view(dtype=[('observation', 'f4'),('station_name','S1',10)])
>>> print(statdat[:]) # now structured array with char array subtype is returned
[(123. , [b'B', b'o', b'u', b'l', b'd', b'e', b'r', b'', b'', b''])
( 3.14, [b'N', b'e', b'w', b' ', b'Y', b'o', b'r', b'k', b'', b''])]
>>> nc.close()
</code></pre>
<p>Note that there is currently no support for mapping numpy structured arrays with
unicode elements (dtype <code>U#</code>) onto netCDF compound types, nor is there support
for netCDF compound types with vlen string components.</p>
<h2 id="in-memory-diskless-datasets">In-memory (diskless) Datasets</h2>
<p>You can create netCDF Datasets whose content is held in memory
instead of in a disk file.
There are two ways to do this.
If you
don't need access to the memory buffer containing the Dataset from
within python, the best way is to use the <code>diskless=True</code> keyword
argument when creating the Dataset.
If you want to save the Dataset
to disk when you close it, also set <code>persist=True</code>.
If you want to
create a new read-only Dataset from an existing python memory buffer, use the
<code>memory</code> keyword argument to pass the memory buffer when creating the Dataset.
If you want to create a new in-memory Dataset, and then access the memory buffer
directly from Python, use the <code>memory</code> keyword argument to specify the
estimated size of the Dataset in bytes when creating the Dataset with
<code>mode='w'</code>.
Then, the <code><a title="netCDF4.Dataset.close" href="#netCDF4.Dataset.close">Dataset.close()</a></code> method will return a python memoryview
object representing the Dataset. Below are examples illustrating both
approaches.</p>
<pre><code class="language-python">>>> # create a diskless (in-memory) Dataset,
>>> # and persist the file to disk when it is closed.
>>> nc = Dataset('diskless_example.nc','w',diskless=True,persist=True)
>>> d = nc.createDimension('x',None)
>>> v = nc.createVariable('v',np.int32,'x')
>>> v[0:5] = np.arange(5)
>>> print(nc)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes): x(5)
variables(dimensions): int32 v(x)
groups:
>>> print(nc['v'][:])
[0 1 2 3 4]
>>> nc.close() # file saved to disk
>>> # create an in-memory dataset from an existing python
>>> # python memory buffer.
>>> # read the newly created netcdf file into a python
>>> # bytes object.
>>> with open('diskless_example.nc', 'rb') as f:
... nc_bytes = f.read()
>>> # create a netCDF in-memory dataset from the bytes object.
>>> nc = Dataset('inmemory.nc', memory=nc_bytes)
>>> print(nc)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes): x(5)
variables(dimensions): int32 v(x)
groups:
>>> print(nc['v'][:])
[0 1 2 3 4]
>>> nc.close()
>>> # create an in-memory Dataset and retrieve memory buffer
>>> # estimated size is 1028 bytes - this is actually only
>>> # used if format is NETCDF3
>>> # (ignored for NETCDF4/HDF5 files).
>>> nc = Dataset('inmemory.nc', mode='w',memory=1028)
>>> d = nc.createDimension('x',None)
>>> v = nc.createVariable('v',np.int32,'x')
>>> v[0:5] = np.arange(5)
>>> nc_buf = nc.close() # close returns memoryview
>>> print(type(nc_buf))
<class 'memoryview'>
>>> # save nc_buf to disk, read it back in and check.
>>> with open('inmemory.nc', 'wb') as f:
... f.write(nc_buf)
>>> nc = Dataset('inmemory.nc')
>>> print(nc)
<class 'netCDF4._netCDF4.Dataset'>
root group (NETCDF4 data model, file format HDF5):
dimensions(sizes): x(5)
variables(dimensions): int32 v(x)
groups:
>>> print(nc['v'][:])
[0 1 2 3 4]
>>> nc.close()
</code></pre>
<h2 id="support-for-complex-numbers">Support for complex numbers</h2>
<p>Although there is no native support for complex numbers in netCDF, there are
some common conventions for storing them. Two of the most common are to either
use a compound datatype for the real and imaginary components, or a separate
dimension. <code><a title="netCDF4" href="#netCDF4">netCDF4</a></code> supports reading several of these conventions, as well as
writing using one of two conventions (depending on file format). This support
for complex numbers is enabled by setting <code>auto_complex=True</code> when opening a
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>:</p>
<pre><code class="language-python">>>> complex_array = np.array([0 + 0j, 1 + 0j, 0 + 1j, 1 + 1j, 0.25 + 0.75j])
>>> with netCDF4.Dataset("complex.nc", "w", auto_complex=True) as nc:
... nc.createDimension("x", size=len(complex_array))
... var = nc.createVariable("data", "c16", ("x",))
... var[:] = complex_array
... print(var)
<class 'netCDF4._netCDF4.Variable'>
compound data(x)
compound data type: complex128
unlimited dimensions:
current shape = (5,)
</code></pre>
<p>When reading files using <code>auto_complex=True</code>, <code><a title="netCDF4" href="#netCDF4">netCDF4</a></code> will interpret variables
stored using the following conventions as complex numbers:</p>
<ul>
<li>compound datatypes with two <code>float</code> or <code>double</code> members who names begin with
<code>r</code> and <code>i</code> (case insensitive)</li>
<li>a dimension of length 2 named <code>complex</code> or <code>ri</code></li>
</ul>
<p>When writing files using <code>auto_complex=True</code>, <code><a title="netCDF4" href="#netCDF4">netCDF4</a></code> will use:</p>
<ul>
<li>a compound datatype named <code>_PFNC_DOUBLE_COMPLEX_TYPE</code> (or <code>*FLOAT*</code> as
appropriate) with members <code>r</code> and <code>i</code> for netCDF4 formats;</li>
<li>or a dimension of length 2 named <code>_pfnc_complex</code> for netCDF3 or classic
formats.</li>
</ul>
<p>Support for complex numbers is handled via the
<a href="https://github.com/PlasmaFAIR/nc-complex"><code>nc-complex</code></a> library. See there for
further details.</p>
<p><strong>contact</strong>: Jeffrey Whitaker <a href="mailto:jeffrey.s.whitaker@noaa.gov">jeffrey.s.whitaker@noaa.gov</a></p>
<p><strong>copyright</strong>: 2008 by Jeffrey Whitaker.</p>
<p><strong>license</strong>: Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:</p>
<p>The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.</p>
<p>THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.</p>
</section>
<section>
</section>
<section>
</section>
<section>
<h2 class="section-title" id="header-functions">Functions</h2>
<dl>
<dt id="netCDF4.chartostring"><code class="name flex">
<span>def <span class="ident">chartostring</span></span>(<span>b, encoding='utf-8')</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>chartostring(b,encoding='utf-8')</code></strong></p>
<p>convert a character array to a string array with one less dimension.</p>
<p><strong><code>b</code></strong>:
Input character array (numpy datatype <code>'S1'</code> or <code>'U1'</code>).
Will be converted to a array of strings, where each string has a fixed
length of <code>b.shape[-1]</code> characters.</p>
<p>optional kwarg <code>encoding</code> can be used to specify character encoding (default
<code>utf-8</code>). If <code>encoding</code> is 'none' or 'bytes', a <code>numpy.string_</code> byte array is
returned.</p>
<p>returns a numpy string array with datatype <code>'UN'</code> (or <code>'SN'</code>) and shape
<code>b.shape[:-1]</code> where where <code>N=b.shape[-1]</code>.</p></div>
</dd>
<dt id="netCDF4.date2index"><code class="name flex">
<span>def <span class="ident">date2index</span></span>(<span>dates, nctime, calendar=None, select='exact', has_year_zero=None)</span>
</code></dt>
<dd>
<div class="desc"><p>date2index(dates, nctime, calendar=None, select=u'exact', has_year_zero=None)</p>
<p>Return indices of a netCDF time variable corresponding to the given dates.</p>
<p><strong>dates</strong>: A datetime object or a sequence of datetime objects.
The datetime objects should not include a time-zone offset.</p>
<p><strong>nctime</strong>: A netCDF time variable object. The nctime object must have a
<strong>units</strong> attribute.</p>
<p><strong>calendar</strong>: describes the calendar to be used in the time calculations.
All the values currently defined in the
<code>CF metadata convention <http://cfconventions.org/cf-conventions/cf-conventions#calendar></code>__ are supported.
Valid calendars <strong>'standard', 'gregorian', 'proleptic_gregorian'
'noleap', '365_day', '360_day', 'julian', 'all_leap', '366_day'</strong>.
Default is <code>None</code> which means the calendar associated with the first
input datetime instance will be used.</p>
<p><strong>select</strong>: <strong>'exact', 'before', 'after', 'nearest'</strong>
The index selection method. <strong>exact</strong> will return the indices perfectly
matching the dates given. <strong>before</strong> and <strong>after</strong> will return the indices
corresponding to the dates just before or just after the given dates if
an exact match cannot be found. <strong>nearest</strong> will return the indices that
correspond to the closest dates.</p>
<p><strong>has_year_zero</strong>: if set to True, astronomical year numbering
is used and the year zero exists.
If set to False for real-world
calendars, then historical year numbering is used and the year 1 is
preceded by year -1 and no year zero exists.
The defaults are set to conform with
CF version 1.9 conventions (False for 'julian', 'gregorian'/'standard', True
for 'proleptic_gregorian' (ISO 8601) and True for the idealized
calendars 'noleap'/'365_day', '360_day', 366_day'/'all_leap')
The defaults can only be over-ridden for the real-world calendars,
for the the idealized calendars the year zero
always exists and the has_year_zero kwarg is ignored.
This kwarg is not needed to define calendar systems allowed by CF
(the calendar-specific defaults do this).</p>
<p>returns an index (indices) of the netCDF time variable corresponding
to the given datetime object(s).</p></div>
</dd>
<dt id="netCDF4.date2num"><code class="name flex">
<span>def <span class="ident">date2num</span></span>(<span>dates, units, calendar=None, has_year_zero=None, longdouble=False)</span>
</code></dt>
<dd>
<div class="desc"><p>date2num(dates, units, calendar=None, has_year_zero=None, longdouble=False)</p>
<p>Return numeric time values given datetime objects. The units
of the numeric time values are described by the <strong>units</strong> argument
and the <strong>calendar</strong> keyword. The datetime objects must
be in UTC with no time-zone offset.
If there is a
time-zone offset in <strong>units</strong>, it will be applied to the
returned numeric values.</p>
<p><strong>dates</strong>: A datetime object or a sequence of datetime objects.
The datetime objects should not include a time-zone offset. They
can be either native python datetime instances (which use
the proleptic gregorian calendar) or cftime.datetime instances.</p>
<p><strong>units</strong>: a string of the form <strong><time units> since <reference time></strong>
describing the time units. <strong><time units></strong> can be days, hours, minutes,
seconds, milliseconds or microseconds. <strong><reference time></strong> is the time
origin. <strong>months since</strong> is allowed <em>only</em> for the <strong>360_day</strong> calendar
and <strong>common_years since</strong> is allowed <em>only</em> for the <strong>365_day</strong> calendar.</p>
<p><strong>calendar</strong>: describes the calendar to be used in the time calculations.
All the values currently defined in the
<code>CF metadata convention <http://cfconventions.org/cf-conventions/cf-conventions#calendar></code>__ are supported.
Valid calendars <strong>'standard', 'gregorian', 'proleptic_gregorian'
'noleap', '365_day', '360_day', 'julian', 'all_leap', '366_day'</strong>.
Default is <code>None</code> which means the calendar associated with the first
input datetime instance will be used.</p>
<p><strong>has_year_zero</strong>: If set to True, astronomical year numbering
is used and the year zero exists.
If set to False for real-world
calendars, then historical year numbering is used and the year 1 is
preceded by year -1 and no year zero exists.
The defaults are set to conform with
CF version 1.9 conventions (False for 'julian', 'gregorian'/'standard', True
for 'proleptic_gregorian' (ISO 8601) and True for the idealized
calendars 'noleap'/'365_day', '360_day', 366_day'/'all_leap')
Note that CF v1.9 does not specifically mention whether year zero
is allowed in the proleptic_gregorian calendar, but ISO-8601 has
a year zero so we have adopted this as the default.
The defaults can only be over-ridden for the real-world calendars,
for the the idealized calendars the year zero
always exists and the has_year_zero kwarg is ignored.
This kwarg is not needed to define calendar systems allowed by CF
(the calendar-specific defaults do this).</p>
<p><strong>longdouble</strong>: If set True, output is in the long double float type
(numpy.float128) instead of float (numpy.float64), allowing microsecond
accuracy when converting a time value to a date and back again. Otherwise
this is only possible if the discretization of the time variable is an
integer multiple of the units.</p>
<p>returns a numeric time value, or an array of numeric time values
with approximately 1 microsecond accuracy.</p></div>
</dd>
<dt id="netCDF4.get_alignment"><code class="name flex">
<span>def <span class="ident">get_alignment</span></span>(<span>)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>get_alignment()</code></strong></p>
<p>return current netCDF alignment within HDF5 files in a tuple
(threshold,alignment). See netcdf C library documentation for
<code>nc_get_alignment</code> for details. Values can be reset with
<code><a title="netCDF4.set_alignment" href="#netCDF4.set_alignment">set_alignment()</a></code>.</p>
<p>This function was added in netcdf 4.9.0.</p></div>
</dd>
<dt id="netCDF4.get_chunk_cache"><code class="name flex">
<span>def <span class="ident">get_chunk_cache</span></span>(<span>)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>get_chunk_cache()</code></strong></p>
<p>return current netCDF chunk cache information in a tuple (size,nelems,preemption).
See netcdf C library documentation for <code>nc_get_chunk_cache</code> for
details. Values can be reset with <code><a title="netCDF4.set_chunk_cache" href="#netCDF4.set_chunk_cache">set_chunk_cache()</a></code>.</p></div>
</dd>
<dt id="netCDF4.getlibversion"><code class="name flex">
<span>def <span class="ident">getlibversion</span></span>(<span>)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>getlibversion()</code></strong></p>
<p>returns a string describing the version of the netcdf library
used to build the module, and when it was built.</p></div>
</dd>
<dt id="netCDF4.num2date"><code class="name flex">
<span>def <span class="ident">num2date</span></span>(<span>times, units, calendar='standard', only_use_cftime_datetimes=True, only_use_python_datetimes=False, has_year_zero=None)</span>
</code></dt>
<dd>
<div class="desc"><p>num2date(times, units, calendar=u'standard', only_use_cftime_datetimes=True, only_use_python_datetimes=False, has_year_zero=None)</p>
<p>Return datetime objects given numeric time values. The units
of the numeric time values are described by the <strong>units</strong> argument
and the <strong>calendar</strong> keyword. The returned datetime objects represent
UTC with no time-zone offset, even if the specified
<strong>units</strong> contain a time-zone offset.</p>
<p><strong>times</strong>: numeric time values.</p>
<p><strong>units</strong>: a string of the form <strong><time units> since <reference time></strong>
describing the time units. <strong><time units></strong> can be days, hours, minutes,
seconds, milliseconds or microseconds. <strong><reference time></strong> is the time
origin. <strong>months since</strong> is allowed <em>only</em> for the <strong>360_day</strong> calendar
and <strong>common_years since</strong> is allowed <em>only</em> for the <strong>365_day</strong> calendar.</p>
<p><strong>calendar</strong>: describes the calendar used in the time calculations.
All the values currently defined in the
<code>CF metadata convention <http://cfconventions.org/cf-conventions/cf-conventions#calendar></code>__ are supported.
Valid calendars <strong>'standard', 'gregorian', 'proleptic_gregorian'
'noleap', '365_day', '360_day', 'julian', 'all_leap', '366_day'</strong>.
Default is <strong>'standard'</strong>, which is a mixed Julian/Gregorian calendar.</p>
<p><strong>only_use_cftime_datetimes</strong>: if False, python datetime.datetime
objects are returned from num2date where possible; if True dates which
subclass cftime.datetime are returned for all calendars. Default <strong>True</strong>.</p>
<p><strong>only_use_python_datetimes</strong>: always return python datetime.datetime
objects and raise an error if this is not possible. Ignored unless
<strong>only_use_cftime_datetimes=False</strong>. Default <strong>False</strong>.</p>
<p><strong>has_year_zero</strong>: if set to True, astronomical year numbering
is used and the year zero exists.
If set to False for real-world
calendars, then historical year numbering is used and the year 1 is
preceded by year -1 and no year zero exists.
The defaults are set to conform with
CF version 1.9 conventions (False for 'julian', 'gregorian'/'standard', True
for 'proleptic_gregorian' (ISO 8601) and True for the idealized
calendars 'noleap'/'365_day', '360_day', 366_day'/'all_leap')
The defaults can only be over-ridden for the real-world calendars,
for the the idealized calendars the year zero
always exists and the has_year_zero kwarg is ignored.
This kwarg is not needed to define calendar systems allowed by CF
(the calendar-specific defaults do this).</p>
<p>returns a datetime instance, or an array of datetime instances with
microsecond accuracy, if possible.</p>
<p><strong><em>Note</em></strong>: If only_use_cftime_datetimes=False and
use_only_python_datetimes=False, the datetime instances
returned are 'real' python datetime
objects if <strong>calendar='proleptic_gregorian'</strong>, or
<strong>calendar='standard'</strong> or <strong>'gregorian'</strong>
and the date is after the breakpoint between the Julian and
Gregorian calendars (1582-10-15). Otherwise, they are ctime.datetime
objects which support some but not all the methods of native python
datetime objects. The datetime instances
do not contain a time-zone offset, even if the specified <strong>units</strong>
contains one.</p></div>
</dd>
<dt id="netCDF4.rc_get"><code class="name flex">
<span>def <span class="ident">rc_get</span></span>(<span>key)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>rc_get(key)</code></strong></p>
<p>Returns the internal netcdf-c rc table value corresponding to key.
See <a href="https://docs.unidata.ucar.edu/netcdf-c/current/md_auth.html">https://docs.unidata.ucar.edu/netcdf-c/current/md_auth.html</a>
for more information on rc files and values.</p></div>
</dd>
<dt id="netCDF4.rc_set"><code class="name flex">
<span>def <span class="ident">rc_set</span></span>(<span>key, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>rc_set(key, value)</code></strong></p>
<p>Sets the internal netcdf-c rc table value corresponding to key.
See <a href="https://docs.unidata.ucar.edu/netcdf-c/current/md_auth.html">https://docs.unidata.ucar.edu/netcdf-c/current/md_auth.html</a>
for more information on rc files and values.</p></div>
</dd>
<dt id="netCDF4.set_alignment"><code class="name flex">
<span>def <span class="ident">set_alignment</span></span>(<span>threshold, alignment)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_alignment(threshold,alignment)</code></strong></p>
<p>Change the HDF5 file alignment.
See netcdf C library documentation for <code>nc_set_alignment</code> for
details.</p>
<p>This function was added in netcdf 4.9.0.</p></div>
</dd>
<dt id="netCDF4.set_chunk_cache"><code class="name flex">
<span>def <span class="ident">set_chunk_cache</span></span>(<span>size=None, nelems=None, preemption=None)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_chunk_cache(size=None,nelems=None,preemption=None)</code></strong></p>
<p>change netCDF4 chunk cache settings.
See netcdf C library documentation for <code>nc_set_chunk_cache</code> for
details.</p></div>
</dd>
<dt id="netCDF4.stringtoarr"><code class="name flex">
<span>def <span class="ident">stringtoarr</span></span>(<span>string, NUMCHARS, dtype='S')</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>stringtoarr(a, NUMCHARS,dtype='S')</code></strong></p>
<p>convert a string to a character array of length <code>NUMCHARS</code></p>
<p><strong><code>a</code></strong>:
Input python string.</p>
<p><strong><code>NUMCHARS</code></strong>:
number of characters used to represent string
(if len(a) < <code>NUMCHARS</code>, it will be padded on the right with blanks).</p>
<p><strong><code>dtype</code></strong>:
type of numpy array to return.
Default is <code>'S'</code>, which
means an array of dtype <code>'S1'</code> will be returned.
If dtype=<code>'U'</code>, a
unicode array (dtype = <code>'U1'</code>) will be returned.</p>
<p>returns a rank 1 numpy character array of length NUMCHARS with datatype <code>'S1'</code>
(default) or <code>'U1'</code> (if dtype=<code>'U'</code>)</p></div>
</dd>
<dt id="netCDF4.stringtochar"><code class="name flex">
<span>def <span class="ident">stringtochar</span></span>(<span>a, encoding='utf-8')</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>stringtochar(a,encoding='utf-8')</code></strong></p>
<p>convert a string array to a character array with one extra dimension</p>
<p><strong><code>a</code></strong>:
Input numpy string array with numpy datatype <code>'SN'</code> or <code>'UN'</code>, where N
is the number of characters in each string.
Will be converted to
an array of characters (datatype <code>'S1'</code> or <code>'U1'</code>) of shape <code>a.shape + (N,)</code>.</p>
<p>optional kwarg <code>encoding</code> can be used to specify character encoding (default
<code>utf-8</code>). If <code>encoding</code> is 'none' or 'bytes', a <code>numpy.string_</code> the input array
is treated a raw byte strings (<code>numpy.string_</code>).</p>
<p>returns a numpy character array with datatype <code>'S1'</code> or <code>'U1'</code>
and shape <code>a.shape + (N,)</code>, where N is the length of each string in a.</p></div>
</dd>
</dl>
</section>
<section>
<h2 class="section-title" id="header-classes">Classes</h2>
<dl>
<dt id="netCDF4.CompoundType"><code class="flex name class">
<span>class <span class="ident">CompoundType</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>A <code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> instance is used to describe a compound data
type, and can be passed to the the <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method of
a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance.
Compound data types map to numpy structured arrays.
See <code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> for more details.</p>
<p>The instance variables <code>dtype</code> and <code>name</code> should not be modified by
the user.</p>
<p><strong><em><code>__init__(group, datatype, datatype_name)</code></em></strong></p>
<p>CompoundType constructor.</p>
<p><strong><code>grp</code></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance to associate with the compound datatype.</p>
<p><strong><code>dt</code></strong>: A numpy dtype object describing a structured (a.k.a record)
array.
Can be composed of homogeneous numeric or character data types, or
other structured array data types.</p>
<p><strong><code>dtype_name</code></strong>: a Python string containing a description of the
compound data type.</p>
<p><strong><em>Note 1</em></strong>: When creating nested compound data types,
the inner compound data types must already be associated with CompoundType
instances (so create CompoundType instances for the innermost structures
first).</p>
<p><strong><em>Note 2</em></strong>: <code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> instances should be created using the
<code><a title="netCDF4.Dataset.createCompoundType" href="#netCDF4.Dataset.createCompoundType">Dataset.createCompoundType()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance, not using this class directly.</p></div>
<h3>Instance variables</h3>
<dl>
<dt id="netCDF4.CompoundType.dtype"><code class="name">var <span class="ident">dtype</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.CompoundType.dtype_view"><code class="name">var <span class="ident">dtype_view</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.CompoundType.name"><code class="name">var <span class="ident">name</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.Dataset"><code class="flex name class">
<span>class <span class="ident">Dataset</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>A netCDF <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> is a collection of dimensions, groups, variables and
attributes. Together they describe the meaning of data and relations among
data fields stored in a netCDF file. See <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> for more
details.</p>
<p>A list of attribute names corresponding to global netCDF attributes
defined for the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> can be obtained with the
<code><a title="netCDF4.Dataset.ncattrs" href="#netCDF4.Dataset.ncattrs">Dataset.ncattrs()</a></code> method.
These attributes can be created by assigning to an attribute of the
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance. A dictionary containing all the netCDF attribute
name/value pairs is provided by the <code>__dict__</code> attribute of a
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance.</p>
<p>The following class variables are read-only and should not be
modified by the user.</p>
<p><strong><code>dimensions</code></strong>: The <code>dimensions</code> dictionary maps the names of
dimensions defined for the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> to instances of the
<code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> class.</p>
<p><strong><code>variables</code></strong>: The <code>variables</code> dictionary maps the names of variables
defined for this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> to instances of the
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> class.</p>
<p><strong><code>groups</code></strong>: The groups dictionary maps the names of groups created for
this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> to instances of the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> class (the
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> class is simply a special case of the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> class which
describes the root group in the netCDF4 file).</p>
<p><strong><code>cmptypes</code></strong>: The <code>cmptypes</code> dictionary maps the names of
compound types defined for the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> to instances of the
<code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> class.</p>
<p><strong><code>vltypes</code></strong>: The <code>vltypes</code> dictionary maps the names of
variable-length types defined for the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> to instances
of the <code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> class.</p>
<p><strong><code>enumtypes</code></strong>: The <code>enumtypes</code> dictionary maps the names of
Enum types defined for the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> to instances
of the <code><a title="netCDF4.EnumType" href="#netCDF4.EnumType">EnumType</a></code> class.</p>
<p><strong><code>data_model</code></strong>: <code>data_model</code> describes the netCDF
data model version, one of <code>NETCDF3_CLASSIC</code>, <code>NETCDF4</code>,
<code>NETCDF4_CLASSIC</code>, <code>NETCDF3_64BIT_OFFSET</code> or <code>NETCDF3_64BIT_DATA</code>.</p>
<p><strong><code>file_format</code></strong>: same as <code>data_model</code>, retained for backwards compatibility.</p>
<p><strong><code>disk_format</code></strong>: <code>disk_format</code> describes the underlying
file format, one of <code>NETCDF3</code>, <code>HDF5</code>, <code>HDF4</code>,
<code>PNETCDF</code>, <code>DAP2</code>, <code>DAP4</code> or <code>UNDEFINED</code>. Only available if using
netcdf C library version >= 4.3.1, otherwise will always return
<code>UNDEFINED</code>.</p>
<p><strong><code>parent</code></strong>: <code>parent</code> is a reference to the parent
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance. <code>None</code> for the root group or <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>
instance.</p>
<p><strong><code>path</code></strong>: <code>path</code> shows the location of the <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> in
the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> in a unix directory format (the names of groups in the
hierarchy separated by backslashes). A <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance is the root
group, so the path is simply <code>'/'</code>.</p>
<p><strong><code>keepweakref</code></strong>: If <code>True</code>, child Dimension and Variables objects only keep weak
references to the parent Dataset or Group.</p>
<p><strong><code>_ncstring_attrs__</code></strong>: If <code>True</code>, all text attributes will be written as variable-length
strings.</p>
<p><strong><code>__init__(self, filename, mode="r", clobber=True, diskless=False,
persist=False, keepweakref=False, memory=None, encoding=None,
parallel=False, comm=None, info=None, format='NETCDF4')</code></strong></p>
<p><code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> constructor.</p>
<p><strong><code>filename</code></strong>: Name of netCDF file to hold dataset. Can also
be a python 3 pathlib instance or the URL of an OpenDAP dataset.
When memory is
set this is just used to set the <code>filepath()</code>.</p>
<p><strong><code>mode</code></strong>: access mode. <code>r</code> means read-only; no data can be
modified. <code>w</code> means write; a new file is created, an existing file with
the same name is deleted. <code>x</code> means write, but fail if an existing
file with the same name already exists. <code>a</code> and <code>r+</code> mean append;
an existing file is opened for reading and writing, if
file does not exist already, one is created.
Appending <code>s</code> to modes <code>r</code>, <code>w</code>, <code>r+</code> or <code>a</code> will enable unbuffered shared
access to <code>NETCDF3_CLASSIC</code>, <code>NETCDF3_64BIT_OFFSET</code> or
<code>NETCDF3_64BIT_DATA</code> formatted files.
Unbuffered access may be useful even if you don't need shared
access, since it may be faster for programs that don't access data
sequentially. This option is ignored for <code>NETCDF4</code> and <code>NETCDF4_CLASSIC</code>
formatted files.</p>
<p><strong><code>clobber</code></strong>: if <code>True</code> (default), opening a file with <code>mode='w'</code>
will clobber an existing file with the same name.
if <code>False</code>, an
exception will be raised if a file with the same name already exists.
mode=<code>x</code> is identical to mode=<code>w</code> with clobber=False.</p>
<p><strong><code>format</code></strong>: underlying file format (one of <code>'NETCDF4',
'NETCDF4_CLASSIC', 'NETCDF3_CLASSIC'<code>, </code>'NETCDF3_64BIT_OFFSET'</code> or
<code>'NETCDF3_64BIT_DATA'</code>.
Only relevant if <code>mode = 'w'</code> (if <code>mode = 'r','a'</code> or <code>'r+'</code> the file format
is automatically detected). Default <code>'NETCDF4'</code>, which means the data is
stored in an HDF5 file, using netCDF 4 API features.
Setting
<code>format='NETCDF4_CLASSIC'</code> will create an HDF5 file, using only netCDF 3
compatible API features. netCDF 3 clients must be recompiled and linked
against the netCDF 4 library to read files in <code>NETCDF4_CLASSIC</code> format.
<code>'NETCDF3_CLASSIC'</code> is the classic netCDF 3 file format that does not
handle 2+ Gb files. <code>'NETCDF3_64BIT_OFFSET'</code> is the 64-bit offset
version of the netCDF 3 file format, which fully supports 2+ GB files, but
is only compatible with clients linked against netCDF version 3.6.0 or
later. <code>'NETCDF3_64BIT_DATA'</code> is the 64-bit data version of the netCDF 3
file format, which supports 64-bit dimension sizes plus unsigned and
64 bit integer data types, but is only compatible with clients linked against
netCDF version 4.4.0 or later.</p>
<p><strong><code>diskless</code></strong>: If <code>True</code>, create diskless (in-core) file.
This is a feature added to the C library after the
netcdf-4.2 release. If you need to access the memory buffer directly,
use the in-memory feature instead (see <code>memory</code> kwarg).</p>
<p><strong><code>persist</code></strong>: if <code>diskless=True</code>, persist file to disk when closed
(default <code>False</code>).</p>
<p><strong><code>keepweakref</code></strong>: if <code>True</code>, child Dimension and Variable instances will keep weak
references to the parent Dataset or Group object.
Default is <code>False</code>, which
means strong references will be kept.
Having Dimension and Variable instances
keep a strong reference to the parent Dataset instance, which in turn keeps a
reference to child Dimension and Variable instances, creates circular references.
Circular references complicate garbage collection, which may mean increased
memory usage for programs that create may Dataset instances with lots of
Variables. It also will result in the Dataset object never being deleted, which
means it may keep open files alive as well. Setting <code>keepweakref=True</code> allows
Dataset instances to be garbage collected as soon as they go out of scope, potentially
reducing memory usage and open file handles.
However, in many cases this is not
desirable, since the associated Variable instances may still be needed, but are
rendered unusable when the parent Dataset instance is garbage collected.</p>
<p><strong><code>memory</code></strong>: if not <code>None</code>, create or open an in-memory Dataset.
If mode = <code>r</code>, the memory kwarg must contain a memory buffer object
(an object that supports the python buffer interface).
The Dataset will then be created with contents taken from this block of memory.
If mode = <code>w</code>, the memory kwarg should contain the anticipated size
of the Dataset in bytes (used only for NETCDF3 files).
A memory
buffer containing a copy of the Dataset is returned by the
<code><a title="netCDF4.Dataset.close" href="#netCDF4.Dataset.close">Dataset.close()</a></code> method. Requires netcdf-c version 4.4.1 for mode=<code>r</code>
netcdf-c 4.6.2 for mode=<code>w</code>. To persist the file to disk, the raw
bytes from the returned buffer can be written into a binary file.
The Dataset can also be re-opened using this memory buffer.</p>
<p><strong><code>encoding</code></strong>: encoding used to encode filename string into bytes.
Default is None (<code>sys.getdefaultfileencoding()</code> is used).</p>
<p><strong><code>parallel</code></strong>: open for parallel access using MPI (requires mpi4py and
parallel-enabled netcdf-c and hdf5 libraries).
Default is <code>False</code>. If
<code>True</code>, <code>comm</code> and <code>info</code> kwargs may also be specified.</p>
<p><strong><code>comm</code></strong>: MPI_Comm object for parallel access. Default <code>None</code>, which
means MPI_COMM_WORLD will be used.
Ignored if <code>parallel=False</code>.</p>
<p><strong><code>info</code></strong>: MPI_Info object for parallel access. Default <code>None</code>, which
means MPI_INFO_NULL will be used.
Ignored if <code>parallel=False</code>.</p>
<p><strong><code>auto_complex</code></strong>: if <code>True</code>, then automatically convert complex number types</p></div>
<h3>Subclasses</h3>
<ul class="hlist">
<li>netCDF4._netCDF4.Group</li>
<li>netCDF4._netCDF4.MFDataset</li>
</ul>
<h3>Static methods</h3>
<dl>
<dt id="netCDF4.Dataset.fromcdl"><code class="name flex">
<span>def <span class="ident">fromcdl</span></span>(<span>cdlfilename, ncfilename=None, mode='a', format='NETCDF4')</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>fromcdl(cdlfilename, ncfilename=None, mode='a',format='NETCDF4')</code></strong></p>
<p>call <a href="https://www.unidata.ucar.edu/software/netcdf/docs/netcdf_utilities_guide.html#ncgen_guide">ncgen</a> via subprocess to create Dataset from <a href="https://www.unidata.ucar.edu/software/netcdf/docs/netcdf_utilities_guide.html#cdl_guide">CDL</a>
text representation. Requires <a href="https://www.unidata.ucar.edu/software/netcdf/docs/netcdf_utilities_guide.html#ncgen_guide">ncgen</a> to be installed and in <code>$PATH</code>.</p>
<p><strong><code>cdlfilename</code></strong>:
CDL file.</p>
<p><strong><code>ncfilename</code></strong>: netCDF file to create. If not given, CDL filename with
suffix replaced by <code>.nc</code> is used..</p>
<p><strong><code>mode</code></strong>:
Access mode to open Dataset (Default <code>'a'</code>).</p>
<p><strong><code>format</code></strong>: underlying file format to use (one of <code>'NETCDF4',
'NETCDF4_CLASSIC', 'NETCDF3_CLASSIC'<code>, </code>'NETCDF3_64BIT_OFFSET'</code> or
<code>'NETCDF3_64BIT_DATA'</code>. Default <code>'NETCDF4'</code>.</p>
<p>Dataset instance for <code>ncfilename</code> is returned.</p></div>
</dd>
</dl>
<h3>Instance variables</h3>
<dl>
<dt id="netCDF4.Dataset.auto_complex"><code class="name">var <span class="ident">auto_complex</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.cmptypes"><code class="name">var <span class="ident">cmptypes</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.data_model"><code class="name">var <span class="ident">data_model</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.dimensions"><code class="name">var <span class="ident">dimensions</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.disk_format"><code class="name">var <span class="ident">disk_format</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.enumtypes"><code class="name">var <span class="ident">enumtypes</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.file_format"><code class="name">var <span class="ident">file_format</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.groups"><code class="name">var <span class="ident">groups</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.keepweakref"><code class="name">var <span class="ident">keepweakref</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.name"><code class="name">var <span class="ident">name</span></code></dt>
<dd>
<div class="desc"><p>string name of Group instance</p></div>
</dd>
<dt id="netCDF4.Dataset.parent"><code class="name">var <span class="ident">parent</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.path"><code class="name">var <span class="ident">path</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.variables"><code class="name">var <span class="ident">variables</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Dataset.vltypes"><code class="name">var <span class="ident">vltypes</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
</dl>
<h3>Methods</h3>
<dl>
<dt id="netCDF4.Dataset.close"><code class="name flex">
<span>def <span class="ident">close</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>close(self)</code></strong></p>
<p>Close the Dataset.</p></div>
</dd>
<dt id="netCDF4.Dataset.createCompoundType"><code class="name flex">
<span>def <span class="ident">createCompoundType</span></span>(<span>self, datatype, datatype_name)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>createCompoundType(self, datatype, datatype_name)</code></strong></p>
<p>Creates a new compound data type named <code>datatype_name</code> from the numpy
dtype object <code>datatype</code>.</p>
<p><strong><em>Note</em></strong>: If the new compound data type contains other compound data types
(i.e. it is a 'nested' compound type, where not all of the elements
are homogeneous numeric data types), then the 'inner' compound types <strong>must</strong> be
created first.</p>
<p>The return value is the <code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> class instance describing the new
datatype.</p></div>
</dd>
<dt id="netCDF4.Dataset.createDimension"><code class="name flex">
<span>def <span class="ident">createDimension</span></span>(<span>self, dimname, size=None)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>createDimension(self, dimname, size=None)</code></strong></p>
<p>Creates a new dimension with the given <code>dimname</code> and <code>size</code>.</p>
<p><code>size</code> must be a positive integer or <code>None</code>, which stands for
"unlimited" (default is <code>None</code>). Specifying a size of 0 also
results in an unlimited dimension. The return value is the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code>
class instance describing the new dimension.
To determine the current
maximum size of the dimension, use the <code>len</code> function on the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code>
instance. To determine if a dimension is 'unlimited', use the
<code><a title="netCDF4.Dimension.isunlimited" href="#netCDF4.Dimension.isunlimited">Dimension.isunlimited()</a></code> method of the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance.</p></div>
</dd>
<dt id="netCDF4.Dataset.createEnumType"><code class="name flex">
<span>def <span class="ident">createEnumType</span></span>(<span>self, datatype, datatype_name, enum_dict)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>createEnumType(self, datatype, datatype_name, enum_dict)</code></strong></p>
<p>Creates a new Enum data type named <code>datatype_name</code> from a numpy
integer dtype object <code>datatype</code>, and a python dictionary
defining the enum fields and values.</p>
<p>The return value is the <code><a title="netCDF4.EnumType" href="#netCDF4.EnumType">EnumType</a></code> class instance describing the new
datatype.</p></div>
</dd>
<dt id="netCDF4.Dataset.createGroup"><code class="name flex">
<span>def <span class="ident">createGroup</span></span>(<span>self, groupname)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>createGroup(self, groupname)</code></strong></p>
<p>Creates a new <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> with the given <code>groupname</code>.</p>
<p>If <code>groupname</code> is specified as a path, using forward slashes as in unix to
separate components, then intermediate groups will be created as necessary
(analogous to <code>mkdir -p</code> in unix).
For example,
<code>createGroup('/GroupA/GroupB/GroupC')</code> will create <code>GroupA</code>,
<code>GroupA/GroupB</code>, and <code>GroupA/GroupB/GroupC</code>, if they don't already exist.
If the specified path describes a group that already exists, no error is
raised.</p>
<p>The return value is a <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> class instance.</p></div>
</dd>
<dt id="netCDF4.Dataset.createVLType"><code class="name flex">
<span>def <span class="ident">createVLType</span></span>(<span>self, datatype, datatype_name)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>createVLType(self, datatype, datatype_name)</code></strong></p>
<p>Creates a new VLEN data type named <code>datatype_name</code> from a numpy
dtype object <code>datatype</code>.</p>
<p>The return value is the <code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> class instance describing the new
datatype.</p></div>
</dd>
<dt id="netCDF4.Dataset.createVariable"><code class="name flex">
<span>def <span class="ident">createVariable</span></span>(<span>self, varname, datatype, dimensions=(), compression=None, zlib=False, complevel=4, shuffle=True, szip_coding='nn', szip_pixels_per_block=8, blosc_shuffle=1, fletcher32=False, contiguous=False, chunksizes=None, endian='native', least_significant_digit=None, significant_digits=None, quantize_mode='BitGroom', fill_value=None, chunk_cache=None)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>createVariable(self, varname, datatype, dimensions=(), compression=None, zlib=False,
complevel=4, shuffle=True, fletcher32=False, contiguous=False, chunksizes=None,
szip_coding='nn', szip_pixels_per_block=8, blosc_shuffle=1,
endian='native', least_significant_digit=None, significant_digits=None, quantize_mode='BitGroom',
fill_value=None, chunk_cache=None)</code></strong></p>
<p>Creates a new variable with the given <code>varname</code>, <code>datatype</code>, and
<code>dimensions</code>. If dimensions are not given, the variable is assumed to be
a scalar.</p>
<p>If <code>varname</code> is specified as a path, using forward slashes as in unix to
separate components, then intermediate groups will be created as necessary
For example, <code>createVariable('/GroupA/GroupB/VarC', float, ('x','y'))</code> will create groups <code>GroupA</code>
and <code>GroupA/GroupB</code>, plus the variable <code>GroupA/GroupB/VarC</code>, if the preceding
groups don't already exist.</p>
<p>The <code>datatype</code> can be a numpy datatype object, or a string that describes
a numpy dtype object (like the <code>dtype.str</code> attribute of a numpy array).
Supported specifiers include: <code>'S1' or 'c' (NC_CHAR), 'i1' or 'b' or 'B'
(NC_BYTE), 'u1' (NC_UBYTE), 'i2' or 'h' or 's' (NC_SHORT), 'u2'
(NC_USHORT), 'i4' or 'i' or 'l' (NC_INT), 'u4' (NC_UINT), 'i8' (NC_INT64),
'u8' (NC_UINT64), 'f4' or 'f' (NC_FLOAT), 'f8' or 'd' (NC_DOUBLE)</code>.
<code>datatype</code> can also be a <code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> instance
(for a structured, or compound array), a <code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> instance
(for a variable-length array), or the python <code>str</code> builtin
(for a variable-length string array). Numpy string and unicode datatypes with
length greater than one are aliases for <code>str</code>.</p>
<p>Data from netCDF variables is presented to python as numpy arrays with
the corresponding data type.</p>
<p><code>dimensions</code> must be a tuple containing <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instances and/or
dimension names (strings) that have been defined
previously using <code><a title="netCDF4.Dataset.createDimension" href="#netCDF4.Dataset.createDimension">Dataset.createDimension()</a></code>. The default value
is an empty tuple, which means the variable is a scalar.</p>
<p>If the optional keyword argument <code>compression</code> is set, the data will be
compressed in the netCDF file using the specified compression algorithm.
Currently <code>zlib</code>,<code>szip</code>,<code>zstd</code>,<code>bzip2</code>,<code>blosc_lz</code>,<code>blosc_lz4</code>,<code>blosc_lz4hc</code>,
<code>blosc_zlib</code> and <code>blosc_zstd</code> are supported.
Default is <code>None</code> (no compression).
All of the compressors except
<code>zlib</code> and <code>szip</code> use the HDF5 plugin architecture.</p>
<p>If the optional keyword <code>zlib</code> is <code>True</code>, the data will be compressed in
the netCDF file using zlib compression (default <code>False</code>).
The use of this option is
deprecated in favor of <code>compression='zlib'</code>.</p>
<p>The optional keyword <code>complevel</code> is an integer between 0 and 9 describing
the level of compression desired (default 4). Ignored if <code>compression=None</code>.
A value of zero disables compression.</p>
<p>If the optional keyword <code>shuffle</code> is <code>True</code>, the HDF5 shuffle filter
will be applied before compressing the data with zlib (default <code>True</code>).
This
significantly improves compression. Default is <code>True</code>. Ignored if
<code>zlib=False</code>.</p>
<p>The optional kwarg <code>blosc_shuffle</code>is
ignored
unless the blosc compressor is used. <code>blosc_shuffle</code> can be 0 (no shuffle),
1 (byte-wise shuffle) or 2 (bit-wise shuffle). Default is 1.</p>
<p>The optional kwargs <code>szip_coding</code> and <code>szip_pixels_per_block</code> are ignored
unless the szip compressor is used. <code>szip_coding</code> can be <code>ec</code> (entropy coding)
or <code>nn</code> (nearest neighbor coding). Default is <code>nn</code>. <code>szip_pixels_per_block</code>
can be 4, 8, 16 or 32 (default 8).</p>
<p>If the optional keyword <code>fletcher32</code> is <code>True</code>, the Fletcher32 HDF5
checksum algorithm is activated to detect errors. Default <code>False</code>.</p>
<p>If the optional keyword <code>contiguous</code> is <code>True</code>, the variable data is
stored contiguously on disk.
Default <code>False</code>. Setting to <code>True</code> for
a variable with an unlimited dimension will trigger an error.
Fixed size variables (with no unlimited dimension) with no compression filters
are contiguous by default.</p>
<p>The optional keyword <code>chunksizes</code> can be used to manually specify the
HDF5 chunksizes for each dimension of the variable.
A detailed discussion of HDF chunking and I/O performance is available
<a href="https://support.hdfgroup.org/HDF5/doc/Advanced/Chunking">here</a>.
The default chunking scheme in the netcdf-c library is discussed
<a href="https://www.unidata.ucar.edu/software/netcdf/documentation/NUG/netcdf_perf_chunking.html">here</a>.
Basically, you want the chunk size for each dimension to match as
closely as possible the size of the data block that users will read
from the file. <code>chunksizes</code> cannot be set if <code>contiguous=True</code>.</p>
<p>The optional keyword <code>endian</code> can be used to control whether the
data is stored in little or big endian format on disk. Possible
values are <code>little, big</code> or <code>native</code> (default). The library
will automatically handle endian conversions when the data is read,
but if the data is always going to be read on a computer with the
opposite format as the one used to create the file, there may be
some performance advantage to be gained by setting the endian-ness.</p>
<p>The optional keyword <code>fill_value</code> can be used to override the default
netCDF <code>_FillValue</code> (the value that the variable gets filled with before
any data is written to it, defaults given in the dict <code>netCDF4.default_fillvals</code>).
If fill_value is set to <code>False</code>, then the variable is not pre-filled.</p>
<p>If the optional keyword parameters <code>least_significant_digit</code> or <code>significant_digits</code> are
specified, variable data will be truncated (quantized). In conjunction
with <code>compression='zlib'</code> this produces 'lossy', but significantly more
efficient compression. For example, if <code>least_significant_digit=1</code>,
data will be quantized using <code>numpy.around(scale*data)/scale</code>, where
scale = 2**bits, and bits is determined so that a precision of 0.1 is
retained (in this case bits=4). From the
<a href="http://www.esrl.noaa.gov/psl/data/gridded/conventions/cdc_netcdf_standard.shtml">PSL metadata conventions</a>:
"least_significant_digit – power of ten of the smallest decimal place
in unpacked data that is a reliable value." Default is <code>None</code>, or no
quantization, or 'lossless' compression.
If <code>significant_digits=3</code>
then the data will be quantized so that three significant digits are retained, independent
of the floating point exponent. The keyword argument <code>quantize_mode</code> controls
the quantization algorithm (default 'BitGroom', 'BitRound' and
'GranularBitRound' also available).
The 'GranularBitRound'
algorithm may result in better compression for typical geophysical datasets.
This <code>significant_digits</code> kwarg is only available
with netcdf-c >= 4.9.0, and
only works with <code>NETCDF4</code> or <code>NETCDF4_CLASSIC</code> formatted files.</p>
<p>When creating variables in a <code>NETCDF4</code> or <code>NETCDF4_CLASSIC</code> formatted file,
HDF5 creates something called a 'chunk cache' for each variable.
The
default size of the chunk cache may be large enough to completely fill
available memory when creating thousands of variables.
The optional
keyword <code>chunk_cache</code> allows you to reduce (or increase) the size of
the default chunk cache when creating a variable.
The setting only
persists as long as the Dataset is open - you can use the set_var_chunk_cache
method to change it the next time the Dataset is opened.
Warning - messing with this parameter can seriously degrade performance.</p>
<p>The return value is the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> class instance describing the new
variable.</p>
<p>A list of names corresponding to netCDF variable attributes can be
obtained with the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> method <code><a title="netCDF4.Variable.ncattrs" href="#netCDF4.Variable.ncattrs">Variable.ncattrs()</a></code>. A dictionary
containing all the netCDF attribute name/value pairs is provided by
the <code>__dict__</code> attribute of a <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance.</p>
<p><code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instances behave much like array objects. Data can be
assigned to or retrieved from a variable with indexing and slicing
operations on the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance. A <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance has six
Dataset standard attributes: <code>dimensions, dtype, shape, ndim, name</code> and
<code>least_significant_digit</code>. Application programs should never modify
these attributes. The <code>dimensions</code> attribute is a tuple containing the
names of the dimensions associated with this variable. The <code>dtype</code>
attribute is a string describing the variable's data type (<code>i4, f8,
S1,<code> etc). The </code>shape</code> attribute is a tuple describing the current
sizes of all the variable's dimensions. The <code>name</code> attribute is a
string containing the name of the Variable instance.
The <code>least_significant_digit</code>
attributes describes the power of ten of the smallest decimal place in
the data the contains a reliable value.
assigned to the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>
instance. The <code>ndim</code> attribute
is the number of variable dimensions.</p></div>
</dd>
<dt id="netCDF4.Dataset.delncattr"><code class="name flex">
<span>def <span class="ident">delncattr</span></span>(<span>self, name)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>delncattr(self,name,value)</code></strong></p>
<p>delete a netCDF dataset or group attribute.
Use if you need to delete a
netCDF attribute with the same name as one of the reserved python
attributes.</p></div>
</dd>
<dt id="netCDF4.Dataset.filepath"><code class="name flex">
<span>def <span class="ident">filepath</span></span>(<span>self, encoding=None)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>filepath(self,encoding=None)</code></strong></p>
<p>Get the file system path (or the opendap URL) which was used to
open/create the Dataset. Requires netcdf >= 4.1.2.
The path
is decoded into a string using <code>sys.getfilesystemencoding()</code> by default, this can be
changed using the <code>encoding</code> kwarg.</p></div>
</dd>
<dt id="netCDF4.Dataset.get_variables_by_attributes"><code class="name flex">
<span>def <span class="ident">get_variables_by_attributes</span></span>(<span>self, **kwargs)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>get_variables_by_attributes(self, **kwargs)</code></strong></p>
<p>Returns a list of variables that match specific conditions.</p>
<p>Can pass in key=value parameters and variables are returned that
contain all of the matches. For example,</p>
<pre><code class="language-python">>>> # Get variables with x-axis attribute.
>>> vs = nc.get_variables_by_attributes(axis='X')
>>> # Get variables with matching "standard_name" attribute
>>> vs = nc.get_variables_by_attributes(standard_name='northward_sea_water_velocity')
</code></pre>
<p>Can pass in key=callable parameter and variables are returned if the
callable returns True.
The callable should accept a single parameter,
the attribute value.
None is given as the attribute value when the
attribute does not exist on the variable. For example,</p>
<pre><code class="language-python">>>> # Get Axis variables
>>> vs = nc.get_variables_by_attributes(axis=lambda v: v in ['X', 'Y', 'Z', 'T'])
>>> # Get variables that don't have an "axis" attribute
>>> vs = nc.get_variables_by_attributes(axis=lambda v: v is None)
>>> # Get variables that have a "grid_mapping" attribute
>>> vs = nc.get_variables_by_attributes(grid_mapping=lambda v: v is not None)
</code></pre></div>
</dd>
<dt id="netCDF4.Dataset.getncattr"><code class="name flex">
<span>def <span class="ident">getncattr</span></span>(<span>self, name, encoding='utf-8')</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>getncattr(self,name)</code></strong></p>
<p>retrieve a netCDF dataset or group attribute.
Use if you need to get a netCDF attribute with the same
name as one of the reserved python attributes.</p>
<p>option kwarg <code>encoding</code> can be used to specify the
character encoding of a string attribute (default is <code>utf-8</code>).</p></div>
</dd>
<dt id="netCDF4.Dataset.has_blosc_filter"><code class="name flex">
<span>def <span class="ident">has_blosc_filter</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>has_blosc_filter(self)</code></strong>
returns True if blosc compression filter is available</p></div>
</dd>
<dt id="netCDF4.Dataset.has_bzip2_filter"><code class="name flex">
<span>def <span class="ident">has_bzip2_filter</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>has_bzip2_filter(self)</code></strong>
returns True if bzip2 compression filter is available</p></div>
</dd>
<dt id="netCDF4.Dataset.has_szip_filter"><code class="name flex">
<span>def <span class="ident">has_szip_filter</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>has_szip_filter(self)</code></strong>
returns True if szip compression filter is available</p></div>
</dd>
<dt id="netCDF4.Dataset.has_zstd_filter"><code class="name flex">
<span>def <span class="ident">has_zstd_filter</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>has_zstd_filter(self)</code></strong>
returns True if zstd compression filter is available</p></div>
</dd>
<dt id="netCDF4.Dataset.isopen"><code class="name flex">
<span>def <span class="ident">isopen</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>isopen(self)</code></strong></p>
<p>Is the Dataset open or closed?</p></div>
</dd>
<dt id="netCDF4.Dataset.ncattrs"><code class="name flex">
<span>def <span class="ident">ncattrs</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>ncattrs(self)</code></strong></p>
<p>return netCDF global attribute names for this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> in a list.</p></div>
</dd>
<dt id="netCDF4.Dataset.renameAttribute"><code class="name flex">
<span>def <span class="ident">renameAttribute</span></span>(<span>self, oldname, newname)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>renameAttribute(self, oldname, newname)</code></strong></p>
<p>rename a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> attribute named <code>oldname</code> to <code>newname</code>.</p></div>
</dd>
<dt id="netCDF4.Dataset.renameDimension"><code class="name flex">
<span>def <span class="ident">renameDimension</span></span>(<span>self, oldname, newname)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>renameDimension(self, oldname, newname)</code></strong></p>
<p>rename a <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> named <code>oldname</code> to <code>newname</code>.</p></div>
</dd>
<dt id="netCDF4.Dataset.renameGroup"><code class="name flex">
<span>def <span class="ident">renameGroup</span></span>(<span>self, oldname, newname)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>renameGroup(self, oldname, newname)</code></strong></p>
<p>rename a <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> named <code>oldname</code> to <code>newname</code> (requires netcdf >= 4.3.1).</p></div>
</dd>
<dt id="netCDF4.Dataset.renameVariable"><code class="name flex">
<span>def <span class="ident">renameVariable</span></span>(<span>self, oldname, newname)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>renameVariable(self, oldname, newname)</code></strong></p>
<p>rename a <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> named <code>oldname</code> to <code>newname</code></p></div>
</dd>
<dt id="netCDF4.Dataset.set_always_mask"><code class="name flex">
<span>def <span class="ident">set_always_mask</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_always_mask(self, True_or_False)</code></strong></p>
<p>Call <code><a title="netCDF4.Variable.set_always_mask" href="#netCDF4.Variable.set_always_mask">Variable.set_always_mask()</a></code> for all variables contained in
this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>, as well as for all
variables in all its subgroups.</p>
<p><strong><code>True_or_False</code></strong>: Boolean determining if automatic conversion of
masked arrays with no missing values to regular numpy arrays shall be
applied for all variables. Default True. Set to False to restore the default behaviour
in versions prior to 1.4.1 (numpy array returned unless missing values are present,
otherwise masked array returned).</p>
<p><strong><em>Note</em></strong>: Calling this function only affects existing
variables. Variables created after calling this function will follow
the default behaviour.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_auto_chartostring"><code class="name flex">
<span>def <span class="ident">set_auto_chartostring</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_chartostring(self, True_or_False)</code></strong></p>
<p>Call <code><a title="netCDF4.Variable.set_auto_chartostring" href="#netCDF4.Variable.set_auto_chartostring">Variable.set_auto_chartostring()</a></code> for all variables contained in this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>, as well as for all variables in all its subgroups.</p>
<p><strong><code>True_or_False</code></strong>: Boolean determining if automatic conversion of
all character arrays <–> string arrays should be performed for
character variables (variables of type <code>NC_CHAR</code> or <code>S1</code>) with the
<code>_Encoding</code> attribute set.</p>
<p><strong><em>Note</em></strong>: Calling this function only affects existing variables. Variables created
after calling this function will follow the default behaviour.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_auto_mask"><code class="name flex">
<span>def <span class="ident">set_auto_mask</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_mask(self, True_or_False)</code></strong></p>
<p>Call <code><a title="netCDF4.Variable.set_auto_mask" href="#netCDF4.Variable.set_auto_mask">Variable.set_auto_mask()</a></code> for all variables contained in this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>, as well as for all variables in all its subgroups. Only affects
Variables with primitive or enum types (not compound or vlen Variables).</p>
<p><strong><code>True_or_False</code></strong>: Boolean determining if automatic conversion to masked arrays
shall be applied for all variables.</p>
<p><strong><em>Note</em></strong>: Calling this function only affects existing variables. Variables created
after calling this function will follow the default behaviour.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_auto_maskandscale"><code class="name flex">
<span>def <span class="ident">set_auto_maskandscale</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_maskandscale(self, True_or_False)</code></strong></p>
<p>Call <code><a title="netCDF4.Variable.set_auto_maskandscale" href="#netCDF4.Variable.set_auto_maskandscale">Variable.set_auto_maskandscale()</a></code> for all variables contained in this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>, as well as for all variables in all its subgroups.</p>
<p><strong><code>True_or_False</code></strong>: Boolean determining if automatic conversion to masked arrays
and variable scaling shall be applied for all variables.</p>
<p><strong><em>Note</em></strong>: Calling this function only affects existing variables. Variables created
after calling this function will follow the default behaviour.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_auto_scale"><code class="name flex">
<span>def <span class="ident">set_auto_scale</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_scale(self, True_or_False)</code></strong></p>
<p>Call <code><a title="netCDF4.Variable.set_auto_scale" href="#netCDF4.Variable.set_auto_scale">Variable.set_auto_scale()</a></code> for all variables contained in this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>, as well as for all variables in all its subgroups.</p>
<p><strong><code>True_or_False</code></strong>: Boolean determining if automatic variable scaling
shall be applied for all variables.</p>
<p><strong><em>Note</em></strong>: Calling this function only affects existing variables. Variables created
after calling this function will follow the default behaviour.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_fill_off"><code class="name flex">
<span>def <span class="ident">set_fill_off</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_fill_off(self)</code></strong></p>
<p>Sets the fill mode for a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> open for writing to <code>off</code>.</p>
<p>This will prevent the data from being pre-filled with fill values, which
may result in some performance improvements. However, you must then make
sure the data is actually written before being read.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_fill_on"><code class="name flex">
<span>def <span class="ident">set_fill_on</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_fill_on(self)</code></strong></p>
<p>Sets the fill mode for a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> open for writing to <code>on</code>.</p>
<p>This causes data to be pre-filled with fill values. The fill values can be
controlled by the variable's <code>_Fill_Value</code> attribute, but is usually
sufficient to the use the netCDF default <code>_Fill_Value</code> (defined
separately for each variable type). The default behavior of the netCDF
library corresponds to <code>set_fill_on</code>.
Data which are equal to the
<code>_Fill_Value</code> indicate that the variable was created, but never written
to.</p></div>
</dd>
<dt id="netCDF4.Dataset.set_ncstring_attrs"><code class="name flex">
<span>def <span class="ident">set_ncstring_attrs</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_ncstring_attrs(self, True_or_False)</code></strong></p>
<p>Call <code><a title="netCDF4.Variable.set_ncstring_attrs" href="#netCDF4.Variable.set_ncstring_attrs">Variable.set_ncstring_attrs()</a></code> for all variables contained in
this <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>, as well as for all its
subgroups and their variables.</p>
<p><strong><code>True_or_False</code></strong>: Boolean determining if all string attributes are
created as variable-length NC_STRINGs, (if True), or if ascii text
attributes are stored as NC_CHARs (if False; default)</p>
<p><strong><em>Note</em></strong>: Calling this function only affects newly created attributes
of existing (sub-) groups and their variables.</p></div>
</dd>
<dt id="netCDF4.Dataset.setncattr"><code class="name flex">
<span>def <span class="ident">setncattr</span></span>(<span>self, name, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>setncattr(self,name,value)</code></strong></p>
<p>set a netCDF dataset or group attribute using name,value pair.
Use if you need to set a netCDF attribute with the
with the same name as one of the reserved python attributes.</p></div>
</dd>
<dt id="netCDF4.Dataset.setncattr_string"><code class="name flex">
<span>def <span class="ident">setncattr_string</span></span>(<span>self, name, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>setncattr_string(self,name,value)</code></strong></p>
<p>set a netCDF dataset or group string attribute using name,value pair.
Use if you need to ensure that a netCDF attribute is created with type
<code>NC_STRING</code> if the file format is <code>NETCDF4</code>.</p></div>
</dd>
<dt id="netCDF4.Dataset.setncatts"><code class="name flex">
<span>def <span class="ident">setncatts</span></span>(<span>self, attdict)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>setncatts(self,attdict)</code></strong></p>
<p>set a bunch of netCDF dataset or group attributes at once using a python dictionary.
This may be faster when setting a lot of attributes for a <code>NETCDF3</code>
formatted file, since nc_redef/nc_enddef is not called in between setting
each attribute</p></div>
</dd>
<dt id="netCDF4.Dataset.sync"><code class="name flex">
<span>def <span class="ident">sync</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>sync(self)</code></strong></p>
<p>Writes all buffered data in the <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> to the disk file.</p></div>
</dd>
<dt id="netCDF4.Dataset.tocdl"><code class="name flex">
<span>def <span class="ident">tocdl</span></span>(<span>self, coordvars=False, data=False, outfile=None)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>tocdl(self, coordvars=False, data=False, outfile=None)</code></strong></p>
<p>call <a href="https://www.unidata.ucar.edu/software/netcdf/docs/netcdf_utilities_guide.html#ncdump_guide">ncdump</a> via subprocess to create <a href="https://www.unidata.ucar.edu/software/netcdf/docs/netcdf_utilities_guide.html#cdl_guide">CDL</a>
text representation of Dataset. Requires <a href="https://www.unidata.ucar.edu/software/netcdf/docs/netcdf_utilities_guide.html#ncdump_guide">ncdump</a>
to be installed and in <code>$PATH</code>.</p>
<p><strong><code>coordvars</code></strong>: include coordinate variable data (via <code>ncdump -c</code>). Default False</p>
<p><strong><code>data</code></strong>: if True, write out variable data (Default False).</p>
<p><strong><code>outfile</code></strong>: If not None, file to output ncdump to. Default is to return a string.</p></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.Dimension"><code class="flex name class">
<span>class <span class="ident">Dimension</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>A netCDF <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> is used to describe the coordinates of a <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>.
See <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> for more details.</p>
<p>The current maximum size of a <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance can be obtained by
calling the python <code>len</code> function on the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance. The
<code><a title="netCDF4.Dimension.isunlimited" href="#netCDF4.Dimension.isunlimited">Dimension.isunlimited()</a></code> method of a <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance can be used to
determine if the dimension is unlimited.</p>
<p>Read-only class variables:</p>
<p><strong><code>name</code></strong>: String name, used when creating a <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> with
<code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code>.</p>
<p><strong><code>size</code></strong>: Current <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> size (same as <code>len(d)</code>, where <code>d</code> is a
<code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance).</p>
<p><strong><code>__init__(self, group, name, size=None)</code></strong></p>
<p><code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> constructor.</p>
<p><strong><code>group</code></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance to associate with dimension.</p>
<p><strong><code>name</code></strong>: Name of the dimension.</p>
<p><strong><code>size</code></strong>: Size of the dimension. <code>None</code> or 0 means unlimited. (Default <code>None</code>).</p>
<p><strong><em>Note</em></strong>: <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instances should be created using the
<code><a title="netCDF4.Dataset.createDimension" href="#netCDF4.Dataset.createDimension">Dataset.createDimension()</a></code> method of a <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance, not using <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> directly.</p></div>
<h3>Instance variables</h3>
<dl>
<dt id="netCDF4.Dimension.name"><code class="name">var <span class="ident">name</span></code></dt>
<dd>
<div class="desc"><p>string name of Dimension instance</p></div>
</dd>
<dt id="netCDF4.Dimension.size"><code class="name">var <span class="ident">size</span></code></dt>
<dd>
<div class="desc"><p>current size of Dimension (calls <code>len</code> on Dimension instance)</p></div>
</dd>
</dl>
<h3>Methods</h3>
<dl>
<dt id="netCDF4.Dimension.group"><code class="name flex">
<span>def <span class="ident">group</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>group(self)</code></strong></p>
<p>return the group that this <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> is a member of.</p></div>
</dd>
<dt id="netCDF4.Dimension.isunlimited"><code class="name flex">
<span>def <span class="ident">isunlimited</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>isunlimited(self)</code></strong></p>
<p>returns <code>True</code> if the <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instance is unlimited, <code>False</code> otherwise.</p></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.EnumType"><code class="flex name class">
<span>class <span class="ident">EnumType</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>A <code><a title="netCDF4.EnumType" href="#netCDF4.EnumType">EnumType</a></code> instance is used to describe an Enum data
type, and can be passed to the the <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method of
a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance. See
<code><a title="netCDF4.EnumType" href="#netCDF4.EnumType">EnumType</a></code> for more details.</p>
<p>The instance variables <code>dtype</code>, <code>name</code> and <code>enum_dict</code> should not be modified by
the user.</p>
<p><strong><code>__init__(group, datatype, datatype_name, enum_dict)</code></strong></p>
<p>EnumType constructor.</p>
<p><strong><code>group</code></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance to associate with the VLEN datatype.</p>
<p><strong><code>datatype</code></strong>: An numpy integer dtype object describing the base type
for the Enum.</p>
<p><strong><code>datatype_name</code></strong>: a Python string containing a description of the
Enum data type.</p>
<p><strong><code>enum_dict</code></strong>: a Python dictionary containing the Enum field/value
pairs.</p>
<p><strong><em><code>Note</code></em></strong>: <code><a title="netCDF4.EnumType" href="#netCDF4.EnumType">EnumType</a></code> instances should be created using the
<code><a title="netCDF4.Dataset.createEnumType" href="#netCDF4.Dataset.createEnumType">Dataset.createEnumType()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance, not using this class directly.</p></div>
<h3>Instance variables</h3>
<dl>
<dt id="netCDF4.EnumType.dtype"><code class="name">var <span class="ident">dtype</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.EnumType.enum_dict"><code class="name">var <span class="ident">enum_dict</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.EnumType.name"><code class="name">var <span class="ident">name</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.Group"><code class="flex name class">
<span>class <span class="ident">Group</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>Groups define a hierarchical namespace within a netCDF file. They are
analogous to directories in a unix filesystem. Each <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> behaves like
a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> within a Dataset, and can contain it's own variables,
dimensions and attributes (and other Groups). See <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>
for more details.</p>
<p><code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> inherits from <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code>, so all the
<code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> class methods and variables are available
to a <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance (except the <code>close</code> method).</p>
<p>Additional read-only class variables:</p>
<p><strong><code>name</code></strong>: String describing the group name.</p>
<p><strong><code>__init__(self, parent, name)</code></strong>
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> constructor.</p>
<p><strong><code>parent</code></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance for the parent group.
If being created
in the root group, use a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance.</p>
<p><strong><code>name</code></strong>: - Name of the group.</p>
<p><strong><em>Note</em></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instances should be created using the
<code><a title="netCDF4.Dataset.createGroup" href="#netCDF4.Dataset.createGroup">Dataset.createGroup()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance, or
another <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance, not using this class directly.</p></div>
<h3>Ancestors</h3>
<ul class="hlist">
<li>netCDF4._netCDF4.Dataset</li>
</ul>
<h3>Methods</h3>
<dl>
<dt id="netCDF4.Group.close"><code class="name flex">
<span>def <span class="ident">close</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>close(self)</code></strong></p>
<p>overrides <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> close method which does not apply to <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code>
instances, raises OSError.</p></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.MFDataset"><code class="flex name class">
<span>class <span class="ident">MFDataset</span></span>
<span>(</span><span>files, check=False, aggdim=None, exclude=[], master_file=None)</span>
</code></dt>
<dd>
<div class="desc"><p>Class for reading multi-file netCDF Datasets, making variables
spanning multiple files appear as if they were in one file.
Datasets must be in <code>NETCDF4_CLASSIC, NETCDF3_CLASSIC, NETCDF3_64BIT_OFFSET
or NETCDF3_64BIT_DATA<code> format (</code>NETCDF4</code> Datasets won't work).</p>
<p>Adapted from <a href="http://pysclint.sourceforge.net/pycdf">pycdf</a> by Andre Gosselin.</p>
<p>Example usage (See <code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code> for more details):</p>
<pre><code class="language-python">>>> import numpy as np
>>> # create a series of netCDF files with a variable sharing
>>> # the same unlimited dimension.
>>> for nf in range(10):
... with Dataset("mftest%s.nc" % nf, "w", format='NETCDF4_CLASSIC') as f:
... f.createDimension("x",None)
... x = f.createVariable("x","i",("x",))
... x[0:10] = np.arange(nf*10,10*(nf+1))
>>> # now read all those files in at once, in one Dataset.
>>> f = MFDataset("mftest*nc")
>>> print(f.variables["x"][:])
[ 0 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]
</code></pre>
<p><strong><code>__init__(self, files, check=False, aggdim=None, exclude=[],
master_file=None)</code></strong></p>
<p>Open a Dataset spanning multiple files, making it look as if it was a
single file. Variables in the list of files that share the same
dimension (specified with the keyword <code>aggdim</code>) are aggregated. If
<code>aggdim</code> is not specified, the unlimited is aggregated.
Currently,
<code>aggdim</code> must be the leftmost (slowest varying) dimension of each
of the variables to be aggregated.</p>
<p><strong><code>files</code></strong>: either a sequence of netCDF files or a string with a
wildcard (converted to a sorted list of files using glob)
If
the <code>master_file</code> kwarg is not specified, the first file
in the list will become the "master" file, defining all the
variables with an aggregation dimension which may span
subsequent files. Attribute access returns attributes only from "master"
file. The files are always opened in read-only mode.</p>
<p><strong><code>check</code></strong>: True if you want to do consistency checking to ensure the
correct variables structure for all of the netcdf files.
Checking makes
the initialization of the MFDataset instance much slower. Default is
False.</p>
<p><strong><code>aggdim</code></strong>: The name of the dimension to aggregate over (must
be the leftmost dimension of each of the variables to be aggregated).
If None (default), aggregate over the unlimited dimension.</p>
<p><strong><code>exclude</code></strong>: A list of variable names to exclude from aggregation.
Default is an empty list.</p>
<p><strong><code>master_file</code></strong>: file to use as "master file", defining all the
variables with an aggregation dimension and all global attributes.</p></div>
<h3>Ancestors</h3>
<ul class="hlist">
<li>netCDF4._netCDF4.Dataset</li>
</ul>
<h3>Methods</h3>
<dl>
<dt id="netCDF4.MFDataset.close"><code class="name flex">
<span>def <span class="ident">close</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>close(self)</code></strong></p>
<p>close all the open files.</p></div>
</dd>
<dt id="netCDF4.MFDataset.isopen"><code class="name flex">
<span>def <span class="ident">isopen</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>isopen(self)</code></strong></p>
<p>True if all files are open, False otherwise.</p></div>
</dd>
<dt id="netCDF4.MFDataset.ncattrs"><code class="name flex">
<span>def <span class="ident">ncattrs</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>ncattrs(self)</code></strong></p>
<p>return the netcdf attribute names from the master file.</p></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.MFTime"><code class="flex name class">
<span>class <span class="ident">MFTime</span></span>
<span>(</span><span>time, units=None, calendar=None)</span>
</code></dt>
<dd>
<div class="desc"><p>Class providing an interface to a MFDataset time Variable by imposing a unique common
time unit and/or calendar to all files.</p>
<p>Example usage (See <code><a title="netCDF4.MFTime" href="#netCDF4.MFTime">MFTime</a></code> for more details):</p>
<pre><code class="language-python">>>> import numpy as np
>>> f1 = Dataset("mftest_1.nc","w", format="NETCDF4_CLASSIC")
>>> f2 = Dataset("mftest_2.nc","w", format="NETCDF4_CLASSIC")
>>> f1.createDimension("time",None)
>>> f2.createDimension("time",None)
>>> t1 = f1.createVariable("time","i",("time",))
>>> t2 = f2.createVariable("time","i",("time",))
>>> t1.units = "days since 2000-01-01"
>>> t2.units = "days since 2000-02-01"
>>> t1.calendar = "standard"
>>> t2.calendar = "standard"
>>> t1[:] = np.arange(31)
>>> t2[:] = np.arange(30)
>>> f1.close()
>>> f2.close()
>>> # Read the two files in at once, in one Dataset.
>>> f = MFDataset("mftest_*nc")
>>> t = f.variables["time"]
>>> print(t.units)
days since 2000-01-01
>>> print(t[32]) # The value written in the file, inconsistent with the MF time units.
1
>>> T = MFTime(t)
>>> print(T[32])
32
</code></pre>
<p><strong><code>__init__(self, time, units=None, calendar=None)</code></strong></p>
<p>Create a time Variable with units consistent across a multifile
dataset.</p>
<p><strong><code>time</code></strong>: Time variable from a <code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code>.</p>
<p><strong><code>units</code></strong>: Time units, for example, <code>'days since 1979-01-01'</code>. If <code>None</code>,
use the units from the master variable.</p>
<p><strong><code>calendar</code></strong>: Calendar overload to use across all files, for example,
<code>'standard'</code> or <code>'gregorian'</code>. If <code>None</code>, check that the calendar attribute
is present on each variable and values are unique across files raising a
<code>ValueError</code> otherwise.</p></div>
<h3>Ancestors</h3>
<ul class="hlist">
<li>netCDF4._netCDF4._Variable</li>
</ul>
</dd>
<dt id="netCDF4.VLType"><code class="flex name class">
<span>class <span class="ident">VLType</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>A <code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> instance is used to describe a variable length (VLEN) data
type, and can be passed to the the <code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method of
a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance. See
<code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> for more details.</p>
<p>The instance variables <code>dtype</code> and <code>name</code> should not be modified by
the user.</p>
<p><strong><code>__init__(group, datatype, datatype_name)</code></strong></p>
<p>VLType constructor.</p>
<p><strong><code>group</code></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance to associate with the VLEN datatype.</p>
<p><strong><code>datatype</code></strong>: An numpy dtype object describing the component type for the
variable length array.</p>
<p><strong><code>datatype_name</code></strong>: a Python string containing a description of the
VLEN data type.</p>
<p><strong><em><code>Note</code></em></strong>: <code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> instances should be created using the
<code><a title="netCDF4.Dataset.createVLType" href="#netCDF4.Dataset.createVLType">Dataset.createVLType()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance, not using this class directly.</p></div>
<h3>Instance variables</h3>
<dl>
<dt id="netCDF4.VLType.dtype"><code class="name">var <span class="ident">dtype</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.VLType.name"><code class="name">var <span class="ident">name</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
</dl>
</dd>
<dt id="netCDF4.Variable"><code class="flex name class">
<span>class <span class="ident">Variable</span></span>
<span>(</span><span>...)</span>
</code></dt>
<dd>
<div class="desc"><p>A netCDF <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> is used to read and write netCDF data.
They are
analogous to numpy array objects. See <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> for more
details.</p>
<p>A list of attribute names corresponding to netCDF attributes defined for
the variable can be obtained with the <code><a title="netCDF4.Variable.ncattrs" href="#netCDF4.Variable.ncattrs">Variable.ncattrs()</a></code> method. These
attributes can be created by assigning to an attribute of the
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance. A dictionary containing all the netCDF attribute
name/value pairs is provided by the <code>__dict__</code> attribute of a
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance.</p>
<p>The following class variables are read-only:</p>
<p><strong><code>dimensions</code></strong>: A tuple containing the names of the
dimensions associated with this variable.</p>
<p><strong><code>dtype</code></strong>: A numpy dtype object describing the
variable's data type.</p>
<p><strong><code>ndim</code></strong>: The number of variable dimensions.</p>
<p><strong><code>shape</code></strong>: A tuple with the current shape (length of all dimensions).</p>
<p><strong><code>scale</code></strong>: If True, <code>scale_factor</code> and <code>add_offset</code> are
applied, and signed integer data is automatically converted to
unsigned integer data if the <code>_Unsigned</code> attribute is set to "true" or "True".
Default is <code>True</code>, can be reset using <code><a title="netCDF4.Variable.set_auto_scale" href="#netCDF4.Variable.set_auto_scale">Variable.set_auto_scale()</a></code> and
<code><a title="netCDF4.Variable.set_auto_maskandscale" href="#netCDF4.Variable.set_auto_maskandscale">Variable.set_auto_maskandscale()</a></code> methods.</p>
<p><strong><code>mask</code></strong>: If True, data is automatically converted to/from masked
arrays when missing values or fill values are present. Default is <code>True</code>, can be
reset using <code><a title="netCDF4.Variable.set_auto_mask" href="#netCDF4.Variable.set_auto_mask">Variable.set_auto_mask()</a></code> and <code><a title="netCDF4.Variable.set_auto_maskandscale" href="#netCDF4.Variable.set_auto_maskandscale">Variable.set_auto_maskandscale()</a></code>
methods. Only relevant for Variables with primitive or enum types (ignored
for compound and vlen Variables).</p>
<p><strong><code><a title="netCDF4.chartostring" href="#netCDF4.chartostring">chartostring()</a></code></strong>: If True, data is automatically converted to/from character
arrays to string arrays when the <code>_Encoding</code> variable attribute is set.
Default is <code>True</code>, can be reset using
<code><a title="netCDF4.Variable.set_auto_chartostring" href="#netCDF4.Variable.set_auto_chartostring">Variable.set_auto_chartostring()</a></code> method.</p>
<p><strong><code>least_significant_digit</code></strong>: Describes the power of ten of the
smallest decimal place in the data the contains a reliable value.
Data is
truncated to this decimal place when it is assigned to the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>
instance. If <code>None</code>, the data is not truncated.</p>
<p><strong><code>significant_digits</code></strong>: New in version 1.6.0. Describes the number of significant
digits in the data the contains a reliable value.
Data is
truncated to retain this number of significant digits when it is assigned to the
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instance. If <code>None</code>, the data is not truncated.
Only available with netcdf-c >= 4.9.0,
and only works with <code>NETCDF4</code> or <code>NETCDF4_CLASSIC</code> formatted files.
The number of significant digits used in the quantization of variable data can be
obtained using the <code>Variable.significant_digits</code> method. Default <code>None</code> -
no quantization done.</p>
<p><strong><code>quantize_mode</code></strong>: New in version 1.6.0. Controls
the quantization algorithm (default 'BitGroom', 'BitRound' and
'GranularBitRound' also available).
The 'GranularBitRound'
algorithm may result in better compression for typical geophysical datasets.
Ignored if <code>significant_digits</code> not specified. If 'BitRound' is used, then
<code>significant_digits</code> is interpreted as binary (not decimal) digits.</p>
<p><strong><code>__orthogonal_indexing__</code></strong>: Always <code>True</code>.
Indicates to client code
that the object supports 'orthogonal indexing', which means that slices
that are 1d arrays or lists slice along each dimension independently.
This
behavior is similar to Fortran or Matlab, but different than numpy.</p>
<p><strong><code>datatype</code></strong>: numpy data type (for primitive data types) or VLType/CompoundType
instance (for compound or vlen data types).</p>
<p><strong><code>name</code></strong>: String name.</p>
<p><strong><code>size</code></strong>: The number of stored elements.</p>
<p><strong><code>__init__(self, group, name, datatype, dimensions=(), compression=None, zlib=False,
complevel=4, shuffle=True, szip_coding='nn', szip_pixels_per_block=8,
blosc_shuffle=1, fletcher32=False, contiguous=False,
chunksizes=None, endian='native',
least_significant_digit=None,fill_value=None,chunk_cache=None)</code></strong></p>
<p><code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> constructor.</p>
<p><strong><code>group</code></strong>: <code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> or <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> instance to associate with variable.</p>
<p><strong><code>name</code></strong>: Name of the variable.</p>
<p><strong><code>datatype</code></strong>: <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> data type. Can be specified by providing a
numpy dtype object, or a string that describes a numpy dtype object.
Supported values, corresponding to <code>str</code> attribute of numpy dtype
objects, include <code>'f4'</code> (32-bit floating point), <code>'f8'</code> (64-bit floating
point), <code>'i4'</code> (32-bit signed integer), <code>'i2'</code> (16-bit signed integer),
<code>'i8'</code> (64-bit signed integer), <code>'i4'</code> (8-bit signed integer), <code>'i1'</code>
(8-bit signed integer), <code>'u1'</code> (8-bit unsigned integer), <code>'u2'</code> (16-bit
unsigned integer), <code>'u4'</code> (32-bit unsigned integer), <code>'u8'</code> (64-bit
unsigned integer), or <code>'S1'</code> (single-character string).
From
compatibility with Scientific.IO.NetCDF, the old Numeric single character
typecodes can also be used (<code>'f'</code> instead of <code>'f4'</code>, <code>'d'</code> instead of
<code>'f8'</code>, <code>'h'</code> or <code>'s'</code> instead of <code>'i2'</code>, <code>'b'</code> or <code>'B'</code> instead of
<code>'i1'</code>, <code>'c'</code> instead of <code>'S1'</code>, and <code>'i'</code> or <code>'l'</code> instead of
<code>'i4'</code>). <code>datatype</code> can also be a <code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code> instance
(for a structured, or compound array), a <code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code> instance
(for a variable-length array), or the python <code>str</code> builtin
(for a variable-length string array). Numpy string and unicode datatypes with
length greater than one are aliases for <code>str</code>.</p>
<p><strong><code>dimensions</code></strong>: a tuple containing the variable's Dimension instances
(defined previously with <code>createDimension</code>). Default is an empty tuple
which means the variable is a scalar (and therefore has no dimensions).</p>
<p><strong><code>compression</code></strong>: compression algorithm to use.
Currently <code>zlib</code>,<code>szip</code>,<code>zstd</code>,<code>bzip2</code>,<code>blosc_lz</code>,<code>blosc_lz4</code>,<code>blosc_lz4hc</code>,
<code>blosc_zlib</code> and <code>blosc_zstd</code> are supported.
Default is <code>None</code> (no compression).
All of the compressors except
<code>zlib</code> and <code>szip</code> use the HDF5 plugin architecture.</p>
<p><strong><code>zlib</code></strong>: if <code>True</code>, data assigned to the <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>
instance is compressed on disk. Default <code>False</code>. Deprecated - use
<code>compression='zlib'</code> instead.</p>
<p><strong><code>complevel</code></strong>: the level of compression to use (1 is the fastest,
but poorest compression, 9 is the slowest but best compression). Default 4.
Ignored if <code>compression=None</code> or <code>szip</code>. A value of 0 disables compression.</p>
<p><strong><code>shuffle</code></strong>: if <code>True</code>, the HDF5 shuffle filter is applied
to improve zlib compression. Default <code>True</code>. Ignored unless <code>compression = 'zlib'</code>.</p>
<p><strong><code>blosc_shuffle</code></strong>: shuffle filter inside blosc compressor (only
relevant if compression kwarg set to one of the blosc compressors).
Can be 0 (no blosc shuffle), 1 (bytewise shuffle) or 2 (bitwise
shuffle)). Default is 1. Ignored if blosc compressor not used.</p>
<p><strong><code>szip_coding</code></strong>: szip coding method. Can be <code>ec</code> (entropy coding)
or <code>nn</code> (nearest neighbor coding). Default is <code>nn</code>.
Ignored if szip compressor not used.</p>
<p><strong><code>szip_pixels_per_block</code></strong>: Can be 4,8,16 or 32 (Default 8).
Ignored if szip compressor not used.</p>
<p><strong><code>fletcher32</code></strong>: if <code>True</code> (default <code>False</code>), the Fletcher32 checksum
algorithm is used for error detection.</p>
<p><strong><code>contiguous</code></strong>: if <code>True</code> (default <code>False</code>), the variable data is
stored contiguously on disk.
Default <code>False</code>. Setting to <code>True</code> for
a variable with an unlimited dimension will trigger an error. Fixed
size variables (with no unlimited dimension) with no compression
filters are contiguous by default.</p>
<p><strong><code>chunksizes</code></strong>: Can be used to specify the HDF5 chunksizes for each
dimension of the variable. A detailed discussion of HDF chunking and I/O
performance is available
<a href="https://support.hdfgroup.org/HDF5/doc/Advanced/Chunking">here</a>.
The default chunking scheme in the netcdf-c library is discussed
<a href="https://www.unidata.ucar.edu/software/netcdf/documentation/NUG/netcdf_perf_chunking.html">here</a>.
Basically, you want the chunk size for each dimension to match as
closely as possible the size of the data block that users will read
from the file. <code>chunksizes</code> cannot be set if <code>contiguous=True</code>.</p>
<p><strong><code>endian</code></strong>: Can be used to control whether the
data is stored in little or big endian format on disk. Possible
values are <code>little, big</code> or <code>native</code> (default). The library
will automatically handle endian conversions when the data is read,
but if the data is always going to be read on a computer with the
opposite format as the one used to create the file, there may be
some performance advantage to be gained by setting the endian-ness.
For netCDF 3 files (that don't use HDF5), only <code>endian='native'</code> is allowed.</p>
<p>The <code>compression, zlib, complevel, shuffle, fletcher32, contiguous</code> and <code>chunksizes</code>
keywords are silently ignored for netCDF 3 files that do not use HDF5.</p>
<p><strong><code>least_significant_digit</code></strong>: If this or <code>significant_digits</code> are specified,
variable data will be truncated (quantized).
In conjunction with <code>compression='zlib'</code> this produces
'lossy', but significantly more efficient compression. For example, if
<code>least_significant_digit=1</code>, data will be quantized using
around(scale<em>data)/scale, where scale = 2</em>*bits, and bits is determined
so that a precision of 0.1 is retained (in this case bits=4). Default is
<code>None</code>, or no quantization.</p>
<p><strong><code>significant_digits</code></strong>: New in version 1.6.0.
As described for <code>least_significant_digit</code>
except the number of significant digits retained is prescribed independent
of the floating point exponent. Default <code>None</code> - no quantization done.</p>
<p><strong><code>quantize_mode</code></strong>: New in version 1.6.0. Controls
the quantization algorithm (default 'BitGroom', 'BitRound' and
'GranularBitRound' also available).
The 'GranularBitRound'
algorithm may result in better compression for typical geophysical datasets.
Ignored if <code>significant_digts</code> not specified. If 'BitRound' is used, then
<code>significant_digits</code> is interpreted as binary (not decimal) digits.</p>
<p><strong><code>fill_value</code></strong>:
If specified, the default netCDF fill value (the
value that the variable gets filled with before any data is written to it)
is replaced with this value, and the <code>_FillValue</code> attribute is set.
If fill_value is set to <code>False</code>, then the variable is not pre-filled.
The default netCDF fill values can be found in the dictionary <code>netCDF4.default_fillvals</code>.
If not set, the default fill value will be used but no <code>_FillValue</code> attribute will be created
(this is the default behavior of the netcdf-c library). If you want to use the
default fill value, but have the <code>_FillValue</code> attribute set, use
<code>fill_value='default'</code> (note - this only works for primitive data types). <code><a title="netCDF4.Variable.get_fill_value" href="#netCDF4.Variable.get_fill_value">Variable.get_fill_value()</a></code>
can be used to retrieve the fill value, even if the <code>_FillValue</code> attribute is not set.</p>
<p><strong><code>chunk_cache</code></strong>: If specified, sets the chunk cache size for this variable.
Persists as long as Dataset is open. Use <code>set_var_chunk_cache</code> to
change it when Dataset is re-opened.</p>
<p><strong><em>Note</em></strong>: <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> instances should be created using the
<code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">Dataset.createVariable()</a></code> method of a <code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code> or
<code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code> instance, not using this class directly.</p></div>
<h3>Instance variables</h3>
<dl>
<dt id="netCDF4.Variable.always_mask"><code class="name">var <span class="ident">always_mask</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.auto_complex"><code class="name">var <span class="ident">auto_complex</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.chartostring"><code class="name">var <span class="ident">chartostring</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.datatype"><code class="name">var <span class="ident">datatype</span></code></dt>
<dd>
<div class="desc"><p>numpy data type (for primitive data types) or
VLType/CompoundType/EnumType instance
(for compound, vlen
or enum data types)</p></div>
</dd>
<dt id="netCDF4.Variable.dimensions"><code class="name">var <span class="ident">dimensions</span></code></dt>
<dd>
<div class="desc"><p>get variables's dimension names</p></div>
</dd>
<dt id="netCDF4.Variable.dtype"><code class="name">var <span class="ident">dtype</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.mask"><code class="name">var <span class="ident">mask</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.name"><code class="name">var <span class="ident">name</span></code></dt>
<dd>
<div class="desc"><p>string name of Variable instance</p></div>
</dd>
<dt id="netCDF4.Variable.ndim"><code class="name">var <span class="ident">ndim</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.scale"><code class="name">var <span class="ident">scale</span></code></dt>
<dd>
<div class="desc"></div>
</dd>
<dt id="netCDF4.Variable.shape"><code class="name">var <span class="ident">shape</span></code></dt>
<dd>
<div class="desc"><p>find current sizes of all variable dimensions</p></div>
</dd>
<dt id="netCDF4.Variable.size"><code class="name">var <span class="ident">size</span></code></dt>
<dd>
<div class="desc"><p>Return the number of stored elements.</p></div>
</dd>
</dl>
<h3>Methods</h3>
<dl>
<dt id="netCDF4.Variable.assignValue"><code class="name flex">
<span>def <span class="ident">assignValue</span></span>(<span>self, val)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>assignValue(self, val)</code></strong></p>
<p>assign a value to a scalar variable.
Provided for compatibility with
Scientific.IO.NetCDF, can also be done by assigning to an Ellipsis slice ([…]).</p></div>
</dd>
<dt id="netCDF4.Variable.chunking"><code class="name flex">
<span>def <span class="ident">chunking</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>chunking(self)</code></strong></p>
<p>return variable chunking information.
If the dataset is
defined to be contiguous (and hence there is no chunking) the word 'contiguous'
is returned.
Otherwise, a sequence with the chunksize for
each dimension is returned.</p></div>
</dd>
<dt id="netCDF4.Variable.delncattr"><code class="name flex">
<span>def <span class="ident">delncattr</span></span>(<span>self, name)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>delncattr(self,name,value)</code></strong></p>
<p>delete a netCDF variable attribute.
Use if you need to delete a
netCDF attribute with the same name as one of the reserved python
attributes.</p></div>
</dd>
<dt id="netCDF4.Variable.endian"><code class="name flex">
<span>def <span class="ident">endian</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>endian(self)</code></strong></p>
<p>return endian-ness (<code>little,big,native</code>) of variable (as stored in HDF5 file).</p></div>
</dd>
<dt id="netCDF4.Variable.filters"><code class="name flex">
<span>def <span class="ident">filters</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>filters(self)</code></strong></p>
<p>return dictionary containing HDF5 filter parameters.</p></div>
</dd>
<dt id="netCDF4.Variable.getValue"><code class="name flex">
<span>def <span class="ident">getValue</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>getValue(self)</code></strong></p>
<p>get the value of a scalar variable.
Provided for compatibility with
Scientific.IO.NetCDF, can also be done by slicing with an Ellipsis ([…]).</p></div>
</dd>
<dt id="netCDF4.Variable.get_dims"><code class="name flex">
<span>def <span class="ident">get_dims</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>get_dims(self)</code></strong></p>
<p>return a tuple of <code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code> instances associated with this
<code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code>.</p></div>
</dd>
<dt id="netCDF4.Variable.get_fill_value"><code class="name flex">
<span>def <span class="ident">get_fill_value</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>get_fill_value(self)</code></strong></p>
<p>return the fill value associated with this <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> (returns <code>None</code> if data is not
pre-filled). Works even if default fill value was used, and <code>_FillValue</code> attribute
does not exist.</p></div>
</dd>
<dt id="netCDF4.Variable.get_var_chunk_cache"><code class="name flex">
<span>def <span class="ident">get_var_chunk_cache</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>get_var_chunk_cache(self)</code></strong></p>
<p>return variable chunk cache information in a tuple (size,nelems,preemption).
See netcdf C library documentation for <code>nc_get_var_chunk_cache</code> for
details.</p></div>
</dd>
<dt id="netCDF4.Variable.getncattr"><code class="name flex">
<span>def <span class="ident">getncattr</span></span>(<span>self, name, encoding='utf-8')</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>getncattr(self,name)</code></strong></p>
<p>retrieve a netCDF variable attribute.
Use if you need to set a
netCDF attribute with the same name as one of the reserved python
attributes.</p>
<p>option kwarg <code>encoding</code> can be used to specify the
character encoding of a string attribute (default is <code>utf-8</code>).</p></div>
</dd>
<dt id="netCDF4.Variable.group"><code class="name flex">
<span>def <span class="ident">group</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>group(self)</code></strong></p>
<p>return the group that this <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> is a member of.</p></div>
</dd>
<dt id="netCDF4.Variable.ncattrs"><code class="name flex">
<span>def <span class="ident">ncattrs</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>ncattrs(self)</code></strong></p>
<p>return netCDF attribute names for this <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> in a list.</p></div>
</dd>
<dt id="netCDF4.Variable.quantization"><code class="name flex">
<span>def <span class="ident">quantization</span></span>(<span>self)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>quantization(self)</code></strong></p>
<p>return number of significant digits and the algorithm used in quantization.
Returns None if quantization not active.</p></div>
</dd>
<dt id="netCDF4.Variable.renameAttribute"><code class="name flex">
<span>def <span class="ident">renameAttribute</span></span>(<span>self, oldname, newname)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>renameAttribute(self, oldname, newname)</code></strong></p>
<p>rename a <code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code> attribute named <code>oldname</code> to <code>newname</code>.</p></div>
</dd>
<dt id="netCDF4.Variable.set_always_mask"><code class="name flex">
<span>def <span class="ident">set_always_mask</span></span>(<span>self, always_mask)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_always_mask(self,always_mask)</code></strong></p>
<p>turn on or off conversion of data without missing values to regular
numpy arrays.</p>
<p><code>always_mask</code> is a Boolean determining if automatic conversion of
masked arrays with no missing values to regular numpy arrays shall be
applied. Default is True. Set to False to restore the default behaviour
in versions prior to 1.4.1 (numpy array returned unless missing values are present,
otherwise masked array returned).</p></div>
</dd>
<dt id="netCDF4.Variable.set_auto_chartostring"><code class="name flex">
<span>def <span class="ident">set_auto_chartostring</span></span>(<span>self, chartostring)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_chartostring(self,<a title="netCDF4.chartostring" href="#netCDF4.chartostring">chartostring()</a>)</code></strong></p>
<p>turn on or off automatic conversion of character variable data to and
from numpy fixed length string arrays when the <code>_Encoding</code> variable attribute
is set.</p>
<p>If <code><a title="netCDF4.chartostring" href="#netCDF4.chartostring">chartostring()</a></code> is set to <code>True</code>, when data is read from a character variable
(dtype = <code>S1</code>) that has an <code>_Encoding</code> attribute, it is converted to a numpy
fixed length unicode string array (dtype = <code>UN</code>, where <code>N</code> is the length
of the the rightmost dimension of the variable).
The value of <code>_Encoding</code>
is the unicode encoding that is used to decode the bytes into strings.</p>
<p>When numpy string data is written to a variable it is converted back to
individual bytes, with the number of bytes in each string equalling the
rightmost dimension of the variable.</p>
<p>The default value of <code><a title="netCDF4.chartostring" href="#netCDF4.chartostring">chartostring()</a></code> is <code>True</code>
(automatic conversions are performed).</p></div>
</dd>
<dt id="netCDF4.Variable.set_auto_mask"><code class="name flex">
<span>def <span class="ident">set_auto_mask</span></span>(<span>self, mask)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_mask(self,mask)</code></strong></p>
<p>turn on or off automatic conversion of variable data to and
from masked arrays .</p>
<p>If <code>mask</code> is set to <code>True</code>, when data is read from a variable
it is converted to a masked array if any of the values are exactly
equal to the either the netCDF _FillValue or the value specified by the
missing_value variable attribute. The fill_value of the masked array
is set to the missing_value attribute (if it exists), otherwise
the netCDF _FillValue attribute (which has a default value
for each data type). If the variable has no missing_value attribute, the
_FillValue is used instead. If the variable has valid_min/valid_max and
missing_value attributes, data outside the specified range will be masked.
When data is written to a variable, the masked
array is converted back to a regular numpy array by replacing all the
masked values by the missing_value attribute of the variable (if it
exists).
If the variable has no missing_value attribute, the _FillValue
is used instead.</p>
<p>The default value of <code>mask</code> is <code>True</code>
(automatic conversions are performed).</p></div>
</dd>
<dt id="netCDF4.Variable.set_auto_maskandscale"><code class="name flex">
<span>def <span class="ident">set_auto_maskandscale</span></span>(<span>self, maskandscale)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_maskandscale(self,maskandscale)</code></strong></p>
<p>turn on or off automatic conversion of variable data to and
from masked arrays, automatic packing/unpacking of variable
data using <code>scale_factor</code> and <code>add_offset</code> attributes and
automatic conversion of signed integer data to unsigned integer
data if the <code>_Unsigned</code> attribute exists and is set to "true" (or "True").</p>
<p>If <code>maskandscale</code> is set to <code>True</code>, when data is read from a variable
it is converted to a masked array if any of the values are exactly
equal to the either the netCDF _FillValue or the value specified by the
missing_value variable attribute. The fill_value of the masked array
is set to the missing_value attribute (if it exists), otherwise
the netCDF _FillValue attribute (which has a default value
for each data type). If the variable has no missing_value attribute, the
_FillValue is used instead. If the variable has valid_min/valid_max and
missing_value attributes, data outside the specified range will be masked.
When data is written to a variable, the masked
array is converted back to a regular numpy array by replacing all the
masked values by the missing_value attribute of the variable (if it
exists).
If the variable has no missing_value attribute, the _FillValue
is used instead.</p>
<p>If <code>maskandscale</code> is set to <code>True</code>, and the variable has a
<code>scale_factor</code> or an <code>add_offset</code> attribute, then data read
from that variable is unpacked using::</p>
<pre><code>data = self.scale_factor*data + self.add_offset
</code></pre>
<p>When data is written to a variable it is packed using::</p>
<pre><code>data = (data - self.add_offset)/self.scale_factor
</code></pre>
<p>If either scale_factor is present, but add_offset is missing, add_offset
is assumed zero.
If add_offset is present, but scale_factor is missing,
scale_factor is assumed to be one.
For more information on how <code>scale_factor</code> and <code>add_offset</code> can be
used to provide simple compression, see the
<a href="http://www.esrl.noaa.gov/psl/data/gridded/conventions/cdc_netcdf_standard.shtml">PSL metadata conventions</a>.</p>
<p>In addition, if <code>maskandscale</code> is set to <code>True</code>, and if the variable has an
attribute <code>_Unsigned</code> set to "true", and the variable has a signed integer data type,
a view to the data is returned with the corresponding unsigned integer data type.
This convention is used by the netcdf-java library to save unsigned integer
data in <code>NETCDF3</code> or <code>NETCDF4_CLASSIC</code> files (since the <code>NETCDF3</code>
data model does not have unsigned integer data types).</p>
<p>The default value of <code>maskandscale</code> is <code>True</code>
(automatic conversions are performed).</p></div>
</dd>
<dt id="netCDF4.Variable.set_auto_scale"><code class="name flex">
<span>def <span class="ident">set_auto_scale</span></span>(<span>self, scale)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_auto_scale(self,scale)</code></strong></p>
<p>turn on or off automatic packing/unpacking of variable
data using <code>scale_factor</code> and <code>add_offset</code> attributes.
Also turns on and off automatic conversion of signed integer data
to unsigned integer data if the variable has an <code>_Unsigned</code>
attribute set to "true" or "True".</p>
<p>If <code>scale</code> is set to <code>True</code>, and the variable has a
<code>scale_factor</code> or an <code>add_offset</code> attribute, then data read
from that variable is unpacked using::</p>
<pre><code>data = self.scale_factor*data + self.add_offset
</code></pre>
<p>When data is written to a variable it is packed using::</p>
<pre><code>data = (data - self.add_offset)/self.scale_factor
</code></pre>
<p>If either scale_factor is present, but add_offset is missing, add_offset
is assumed zero.
If add_offset is present, but scale_factor is missing,
scale_factor is assumed to be one.
For more information on how <code>scale_factor</code> and <code>add_offset</code> can be
used to provide simple compression, see the
<a href="http://www.esrl.noaa.gov/psl/data/gridded/conventions/cdc_netcdf_standard.shtml">PSL metadata conventions</a>.</p>
<p>In addition, if <code>scale</code> is set to <code>True</code>, and if the variable has an
attribute <code>_Unsigned</code> set to "true", and the variable has a signed integer data type,
a view to the data is returned with the corresponding unsigned integer datatype.
This convention is used by the netcdf-java library to save unsigned integer
data in <code>NETCDF3</code> or <code>NETCDF4_CLASSIC</code> files (since the <code>NETCDF3</code>
data model does not have unsigned integer data types).</p>
<p>The default value of <code>scale</code> is <code>True</code>
(automatic conversions are performed).</p></div>
</dd>
<dt id="netCDF4.Variable.set_collective"><code class="name flex">
<span>def <span class="ident">set_collective</span></span>(<span>self, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_collective(self,True_or_False)</code></strong></p>
<p>turn on or off collective parallel IO access. Ignored if file is not
open for parallel access.</p></div>
</dd>
<dt id="netCDF4.Variable.set_ncstring_attrs"><code class="name flex">
<span>def <span class="ident">set_ncstring_attrs</span></span>(<span>self, ncstring_attrs)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_always_mask(self,ncstring_attrs)</code></strong></p>
<p>turn on or off creating NC_STRING string attributes.</p>
<p>If <code>ncstring_attrs</code> is set to <code>True</code> then text attributes will be variable-length
NC_STRINGs.</p>
<p>The default value of <code>ncstring_attrs</code> is <code>False</code> (writing ascii text attributes as
NC_CHAR).</p></div>
</dd>
<dt id="netCDF4.Variable.set_var_chunk_cache"><code class="name flex">
<span>def <span class="ident">set_var_chunk_cache</span></span>(<span>self, size=None, nelems=None, preemption=None)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>set_var_chunk_cache(self,size=None,nelems=None,preemption=None)</code></strong></p>
<p>change variable chunk cache settings.
See netcdf C library documentation for <code>nc_set_var_chunk_cache</code> for
details.</p></div>
</dd>
<dt id="netCDF4.Variable.setncattr"><code class="name flex">
<span>def <span class="ident">setncattr</span></span>(<span>self, name, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>setncattr(self,name,value)</code></strong></p>
<p>set a netCDF variable attribute using name,value pair.
Use if you need to set a
netCDF attribute with the same name as one of the reserved python
attributes.</p></div>
</dd>
<dt id="netCDF4.Variable.setncattr_string"><code class="name flex">
<span>def <span class="ident">setncattr_string</span></span>(<span>self, name, value)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>setncattr_string(self,name,value)</code></strong></p>
<p>set a netCDF variable string attribute using name,value pair.
Use if you need to ensure that a netCDF attribute is created with type
<code>NC_STRING</code> if the file format is <code>NETCDF4</code>.
Use if you need to set an attribute to an array of variable-length strings.</p></div>
</dd>
<dt id="netCDF4.Variable.setncatts"><code class="name flex">
<span>def <span class="ident">setncatts</span></span>(<span>self, attdict)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>setncatts(self,attdict)</code></strong></p>
<p>set a bunch of netCDF variable attributes at once using a python dictionary.
This may be faster when setting a lot of attributes for a <code>NETCDF3</code>
formatted file, since nc_redef/nc_enddef is not called in between setting
each attribute</p></div>
</dd>
<dt id="netCDF4.Variable.use_nc_get_vars"><code class="name flex">
<span>def <span class="ident">use_nc_get_vars</span></span>(<span>self, use_nc_get_vars)</span>
</code></dt>
<dd>
<div class="desc"><p><strong><code>use_nc_get_vars(self,_use_get_vars)</code></strong></p>
<p>enable the use of netcdf library routine <code>nc_get_vars</code>
to retrieve strided variable slices.
By default,
<code>nc_get_vars</code> may not used by default (depending on the
version of the netcdf-c library being used) since it may be
slower than multiple calls to the unstrided read routine <code>nc_get_vara</code>.</p></div>
</dd>
</dl>
</dd>
</dl>
</section>
</article>
<nav id="sidebar">
<div class="toc">
<ul>
<li><a href="#version-172">Version 1.7.2</a></li>
<li><a href="#introduction">Introduction</a><ul>
<li><a href="#quick-install">Quick Install</a></li>
<li><a href="#developer-install">Developer Install</a></li>
</ul>
</li>
<li><a href="#tutorial">Tutorial</a><ul>
<li><a href="#creatingopeningclosing-a-netcdf-file">Creating/Opening/Closing a netCDF file</a></li>
<li><a href="#groups-in-a-netcdf-file">Groups in a netCDF file</a></li>
<li><a href="#dimensions-in-a-netcdf-file">Dimensions in a netCDF file</a></li>
<li><a href="#variables-in-a-netcdf-file">Variables in a netCDF file</a></li>
<li><a href="#attributes-in-a-netcdf-file">Attributes in a netCDF file</a></li>
<li><a href="#writing-data-to-and-retrieving-data-from-a-netcdf-variable">Writing data to and retrieving data from a netCDF variable</a></li>
<li><a href="#dealing-with-time-coordinates">Dealing with time coordinates</a></li>
<li><a href="#reading-data-from-a-multi-file-netcdf-dataset">Reading data from a multi-file netCDF dataset</a></li>
<li><a href="#efficient-compression-of-netcdf-variables">Efficient compression of netCDF variables</a></li>
<li><a href="#beyond-homogeneous-arrays-of-a-fixed-type-compound-data-types">Beyond homogeneous arrays of a fixed type - compound data types</a></li>
<li><a href="#variable-length-vlen-data-types">Variable-length (vlen) data types</a></li>
<li><a href="#enum-data-type">Enum data type</a></li>
<li><a href="#parallel-io">Parallel IO</a></li>
<li><a href="#dealing-with-strings">Dealing with strings</a></li>
<li><a href="#in-memory-diskless-datasets">In-memory (diskless) Datasets</a></li>
<li><a href="#support-for-complex-numbers">Support for complex numbers</a></li>
</ul>
</li>
</ul>
</div>
<ul id="index">
<li><h3><a href="#header-functions">Functions</a></h3>
<ul class="two-column">
<li><code><a title="netCDF4.chartostring" href="#netCDF4.chartostring">chartostring</a></code></li>
<li><code><a title="netCDF4.date2index" href="#netCDF4.date2index">date2index</a></code></li>
<li><code><a title="netCDF4.date2num" href="#netCDF4.date2num">date2num</a></code></li>
<li><code><a title="netCDF4.get_alignment" href="#netCDF4.get_alignment">get_alignment</a></code></li>
<li><code><a title="netCDF4.get_chunk_cache" href="#netCDF4.get_chunk_cache">get_chunk_cache</a></code></li>
<li><code><a title="netCDF4.getlibversion" href="#netCDF4.getlibversion">getlibversion</a></code></li>
<li><code><a title="netCDF4.num2date" href="#netCDF4.num2date">num2date</a></code></li>
<li><code><a title="netCDF4.rc_get" href="#netCDF4.rc_get">rc_get</a></code></li>
<li><code><a title="netCDF4.rc_set" href="#netCDF4.rc_set">rc_set</a></code></li>
<li><code><a title="netCDF4.set_alignment" href="#netCDF4.set_alignment">set_alignment</a></code></li>
<li><code><a title="netCDF4.set_chunk_cache" href="#netCDF4.set_chunk_cache">set_chunk_cache</a></code></li>
<li><code><a title="netCDF4.stringtoarr" href="#netCDF4.stringtoarr">stringtoarr</a></code></li>
<li><code><a title="netCDF4.stringtochar" href="#netCDF4.stringtochar">stringtochar</a></code></li>
</ul>
</li>
<li><h3><a href="#header-classes">Classes</a></h3>
<ul>
<li>
<h4><code><a title="netCDF4.CompoundType" href="#netCDF4.CompoundType">CompoundType</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.CompoundType.dtype" href="#netCDF4.CompoundType.dtype">dtype</a></code></li>
<li><code><a title="netCDF4.CompoundType.dtype_view" href="#netCDF4.CompoundType.dtype_view">dtype_view</a></code></li>
<li><code><a title="netCDF4.CompoundType.name" href="#netCDF4.CompoundType.name">name</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.Dataset" href="#netCDF4.Dataset">Dataset</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.Dataset.auto_complex" href="#netCDF4.Dataset.auto_complex">auto_complex</a></code></li>
<li><code><a title="netCDF4.Dataset.close" href="#netCDF4.Dataset.close">close</a></code></li>
<li><code><a title="netCDF4.Dataset.cmptypes" href="#netCDF4.Dataset.cmptypes">cmptypes</a></code></li>
<li><code><a title="netCDF4.Dataset.createCompoundType" href="#netCDF4.Dataset.createCompoundType">createCompoundType</a></code></li>
<li><code><a title="netCDF4.Dataset.createDimension" href="#netCDF4.Dataset.createDimension">createDimension</a></code></li>
<li><code><a title="netCDF4.Dataset.createEnumType" href="#netCDF4.Dataset.createEnumType">createEnumType</a></code></li>
<li><code><a title="netCDF4.Dataset.createGroup" href="#netCDF4.Dataset.createGroup">createGroup</a></code></li>
<li><code><a title="netCDF4.Dataset.createVLType" href="#netCDF4.Dataset.createVLType">createVLType</a></code></li>
<li><code><a title="netCDF4.Dataset.createVariable" href="#netCDF4.Dataset.createVariable">createVariable</a></code></li>
<li><code><a title="netCDF4.Dataset.data_model" href="#netCDF4.Dataset.data_model">data_model</a></code></li>
<li><code><a title="netCDF4.Dataset.delncattr" href="#netCDF4.Dataset.delncattr">delncattr</a></code></li>
<li><code><a title="netCDF4.Dataset.dimensions" href="#netCDF4.Dataset.dimensions">dimensions</a></code></li>
<li><code><a title="netCDF4.Dataset.disk_format" href="#netCDF4.Dataset.disk_format">disk_format</a></code></li>
<li><code><a title="netCDF4.Dataset.enumtypes" href="#netCDF4.Dataset.enumtypes">enumtypes</a></code></li>
<li><code><a title="netCDF4.Dataset.file_format" href="#netCDF4.Dataset.file_format">file_format</a></code></li>
<li><code><a title="netCDF4.Dataset.filepath" href="#netCDF4.Dataset.filepath">filepath</a></code></li>
<li><code><a title="netCDF4.Dataset.fromcdl" href="#netCDF4.Dataset.fromcdl">fromcdl</a></code></li>
<li><code><a title="netCDF4.Dataset.get_variables_by_attributes" href="#netCDF4.Dataset.get_variables_by_attributes">get_variables_by_attributes</a></code></li>
<li><code><a title="netCDF4.Dataset.getncattr" href="#netCDF4.Dataset.getncattr">getncattr</a></code></li>
<li><code><a title="netCDF4.Dataset.groups" href="#netCDF4.Dataset.groups">groups</a></code></li>
<li><code><a title="netCDF4.Dataset.has_blosc_filter" href="#netCDF4.Dataset.has_blosc_filter">has_blosc_filter</a></code></li>
<li><code><a title="netCDF4.Dataset.has_bzip2_filter" href="#netCDF4.Dataset.has_bzip2_filter">has_bzip2_filter</a></code></li>
<li><code><a title="netCDF4.Dataset.has_szip_filter" href="#netCDF4.Dataset.has_szip_filter">has_szip_filter</a></code></li>
<li><code><a title="netCDF4.Dataset.has_zstd_filter" href="#netCDF4.Dataset.has_zstd_filter">has_zstd_filter</a></code></li>
<li><code><a title="netCDF4.Dataset.isopen" href="#netCDF4.Dataset.isopen">isopen</a></code></li>
<li><code><a title="netCDF4.Dataset.keepweakref" href="#netCDF4.Dataset.keepweakref">keepweakref</a></code></li>
<li><code><a title="netCDF4.Dataset.name" href="#netCDF4.Dataset.name">name</a></code></li>
<li><code><a title="netCDF4.Dataset.ncattrs" href="#netCDF4.Dataset.ncattrs">ncattrs</a></code></li>
<li><code><a title="netCDF4.Dataset.parent" href="#netCDF4.Dataset.parent">parent</a></code></li>
<li><code><a title="netCDF4.Dataset.path" href="#netCDF4.Dataset.path">path</a></code></li>
<li><code><a title="netCDF4.Dataset.renameAttribute" href="#netCDF4.Dataset.renameAttribute">renameAttribute</a></code></li>
<li><code><a title="netCDF4.Dataset.renameDimension" href="#netCDF4.Dataset.renameDimension">renameDimension</a></code></li>
<li><code><a title="netCDF4.Dataset.renameGroup" href="#netCDF4.Dataset.renameGroup">renameGroup</a></code></li>
<li><code><a title="netCDF4.Dataset.renameVariable" href="#netCDF4.Dataset.renameVariable">renameVariable</a></code></li>
<li><code><a title="netCDF4.Dataset.set_always_mask" href="#netCDF4.Dataset.set_always_mask">set_always_mask</a></code></li>
<li><code><a title="netCDF4.Dataset.set_auto_chartostring" href="#netCDF4.Dataset.set_auto_chartostring">set_auto_chartostring</a></code></li>
<li><code><a title="netCDF4.Dataset.set_auto_mask" href="#netCDF4.Dataset.set_auto_mask">set_auto_mask</a></code></li>
<li><code><a title="netCDF4.Dataset.set_auto_maskandscale" href="#netCDF4.Dataset.set_auto_maskandscale">set_auto_maskandscale</a></code></li>
<li><code><a title="netCDF4.Dataset.set_auto_scale" href="#netCDF4.Dataset.set_auto_scale">set_auto_scale</a></code></li>
<li><code><a title="netCDF4.Dataset.set_fill_off" href="#netCDF4.Dataset.set_fill_off">set_fill_off</a></code></li>
<li><code><a title="netCDF4.Dataset.set_fill_on" href="#netCDF4.Dataset.set_fill_on">set_fill_on</a></code></li>
<li><code><a title="netCDF4.Dataset.set_ncstring_attrs" href="#netCDF4.Dataset.set_ncstring_attrs">set_ncstring_attrs</a></code></li>
<li><code><a title="netCDF4.Dataset.setncattr" href="#netCDF4.Dataset.setncattr">setncattr</a></code></li>
<li><code><a title="netCDF4.Dataset.setncattr_string" href="#netCDF4.Dataset.setncattr_string">setncattr_string</a></code></li>
<li><code><a title="netCDF4.Dataset.setncatts" href="#netCDF4.Dataset.setncatts">setncatts</a></code></li>
<li><code><a title="netCDF4.Dataset.sync" href="#netCDF4.Dataset.sync">sync</a></code></li>
<li><code><a title="netCDF4.Dataset.tocdl" href="#netCDF4.Dataset.tocdl">tocdl</a></code></li>
<li><code><a title="netCDF4.Dataset.variables" href="#netCDF4.Dataset.variables">variables</a></code></li>
<li><code><a title="netCDF4.Dataset.vltypes" href="#netCDF4.Dataset.vltypes">vltypes</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.Dimension" href="#netCDF4.Dimension">Dimension</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.Dimension.group" href="#netCDF4.Dimension.group">group</a></code></li>
<li><code><a title="netCDF4.Dimension.isunlimited" href="#netCDF4.Dimension.isunlimited">isunlimited</a></code></li>
<li><code><a title="netCDF4.Dimension.name" href="#netCDF4.Dimension.name">name</a></code></li>
<li><code><a title="netCDF4.Dimension.size" href="#netCDF4.Dimension.size">size</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.EnumType" href="#netCDF4.EnumType">EnumType</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.EnumType.dtype" href="#netCDF4.EnumType.dtype">dtype</a></code></li>
<li><code><a title="netCDF4.EnumType.enum_dict" href="#netCDF4.EnumType.enum_dict">enum_dict</a></code></li>
<li><code><a title="netCDF4.EnumType.name" href="#netCDF4.EnumType.name">name</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.Group" href="#netCDF4.Group">Group</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.Group.close" href="#netCDF4.Group.close">close</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.MFDataset" href="#netCDF4.MFDataset">MFDataset</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.MFDataset.close" href="#netCDF4.MFDataset.close">close</a></code></li>
<li><code><a title="netCDF4.MFDataset.isopen" href="#netCDF4.MFDataset.isopen">isopen</a></code></li>
<li><code><a title="netCDF4.MFDataset.ncattrs" href="#netCDF4.MFDataset.ncattrs">ncattrs</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.MFTime" href="#netCDF4.MFTime">MFTime</a></code></h4>
</li>
<li>
<h4><code><a title="netCDF4.VLType" href="#netCDF4.VLType">VLType</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.VLType.dtype" href="#netCDF4.VLType.dtype">dtype</a></code></li>
<li><code><a title="netCDF4.VLType.name" href="#netCDF4.VLType.name">name</a></code></li>
</ul>
</li>
<li>
<h4><code><a title="netCDF4.Variable" href="#netCDF4.Variable">Variable</a></code></h4>
<ul class="">
<li><code><a title="netCDF4.Variable.always_mask" href="#netCDF4.Variable.always_mask">always_mask</a></code></li>
<li><code><a title="netCDF4.Variable.assignValue" href="#netCDF4.Variable.assignValue">assignValue</a></code></li>
<li><code><a title="netCDF4.Variable.auto_complex" href="#netCDF4.Variable.auto_complex">auto_complex</a></code></li>
<li><code><a title="netCDF4.Variable.chartostring" href="#netCDF4.Variable.chartostring">chartostring</a></code></li>
<li><code><a title="netCDF4.Variable.chunking" href="#netCDF4.Variable.chunking">chunking</a></code></li>
<li><code><a title="netCDF4.Variable.datatype" href="#netCDF4.Variable.datatype">datatype</a></code></li>
<li><code><a title="netCDF4.Variable.delncattr" href="#netCDF4.Variable.delncattr">delncattr</a></code></li>
<li><code><a title="netCDF4.Variable.dimensions" href="#netCDF4.Variable.dimensions">dimensions</a></code></li>
<li><code><a title="netCDF4.Variable.dtype" href="#netCDF4.Variable.dtype">dtype</a></code></li>
<li><code><a title="netCDF4.Variable.endian" href="#netCDF4.Variable.endian">endian</a></code></li>
<li><code><a title="netCDF4.Variable.filters" href="#netCDF4.Variable.filters">filters</a></code></li>
<li><code><a title="netCDF4.Variable.getValue" href="#netCDF4.Variable.getValue">getValue</a></code></li>
<li><code><a title="netCDF4.Variable.get_dims" href="#netCDF4.Variable.get_dims">get_dims</a></code></li>
<li><code><a title="netCDF4.Variable.get_fill_value" href="#netCDF4.Variable.get_fill_value">get_fill_value</a></code></li>
<li><code><a title="netCDF4.Variable.get_var_chunk_cache" href="#netCDF4.Variable.get_var_chunk_cache">get_var_chunk_cache</a></code></li>
<li><code><a title="netCDF4.Variable.getncattr" href="#netCDF4.Variable.getncattr">getncattr</a></code></li>
<li><code><a title="netCDF4.Variable.group" href="#netCDF4.Variable.group">group</a></code></li>
<li><code><a title="netCDF4.Variable.mask" href="#netCDF4.Variable.mask">mask</a></code></li>
<li><code><a title="netCDF4.Variable.name" href="#netCDF4.Variable.name">name</a></code></li>
<li><code><a title="netCDF4.Variable.ncattrs" href="#netCDF4.Variable.ncattrs">ncattrs</a></code></li>
<li><code><a title="netCDF4.Variable.ndim" href="#netCDF4.Variable.ndim">ndim</a></code></li>
<li><code><a title="netCDF4.Variable.quantization" href="#netCDF4.Variable.quantization">quantization</a></code></li>
<li><code><a title="netCDF4.Variable.renameAttribute" href="#netCDF4.Variable.renameAttribute">renameAttribute</a></code></li>
<li><code><a title="netCDF4.Variable.scale" href="#netCDF4.Variable.scale">scale</a></code></li>
<li><code><a title="netCDF4.Variable.set_always_mask" href="#netCDF4.Variable.set_always_mask">set_always_mask</a></code></li>
<li><code><a title="netCDF4.Variable.set_auto_chartostring" href="#netCDF4.Variable.set_auto_chartostring">set_auto_chartostring</a></code></li>
<li><code><a title="netCDF4.Variable.set_auto_mask" href="#netCDF4.Variable.set_auto_mask">set_auto_mask</a></code></li>
<li><code><a title="netCDF4.Variable.set_auto_maskandscale" href="#netCDF4.Variable.set_auto_maskandscale">set_auto_maskandscale</a></code></li>
<li><code><a title="netCDF4.Variable.set_auto_scale" href="#netCDF4.Variable.set_auto_scale">set_auto_scale</a></code></li>
<li><code><a title="netCDF4.Variable.set_collective" href="#netCDF4.Variable.set_collective">set_collective</a></code></li>
<li><code><a title="netCDF4.Variable.set_ncstring_attrs" href="#netCDF4.Variable.set_ncstring_attrs">set_ncstring_attrs</a></code></li>
<li><code><a title="netCDF4.Variable.set_var_chunk_cache" href="#netCDF4.Variable.set_var_chunk_cache">set_var_chunk_cache</a></code></li>
<li><code><a title="netCDF4.Variable.setncattr" href="#netCDF4.Variable.setncattr">setncattr</a></code></li>
<li><code><a title="netCDF4.Variable.setncattr_string" href="#netCDF4.Variable.setncattr_string">setncattr_string</a></code></li>
<li><code><a title="netCDF4.Variable.setncatts" href="#netCDF4.Variable.setncatts">setncatts</a></code></li>
<li><code><a title="netCDF4.Variable.shape" href="#netCDF4.Variable.shape">shape</a></code></li>
<li><code><a title="netCDF4.Variable.size" href="#netCDF4.Variable.size">size</a></code></li>
<li><code><a title="netCDF4.Variable.use_nc_get_vars" href="#netCDF4.Variable.use_nc_get_vars">use_nc_get_vars</a></code></li>
</ul>
</li>
</ul>
</li>
</ul>
</nav>
</main>
<footer id="footer">
<p>Generated by <a href="https://pdoc3.github.io/pdoc" title="pdoc: Python API documentation generator"><cite>pdoc</cite> 0.11.1</a>.</p>
</footer>
</body>
</html>
|