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
|
# Copyright Cartopy Contributors
#
# This file is part of Cartopy and is released under the LGPL license.
# See COPYING and COPYING.LESSER in the root of the repository for full
# licensing details.
"""
The crs module defines Coordinate Reference Systems and the transformations
between them.
"""
from abc import ABCMeta
from collections import OrderedDict
from functools import lru_cache
import io
import json
import math
import warnings
import numpy as np
import shapely.geometry as sgeom
from pyproj import Transformer
from pyproj.exceptions import ProjError
from shapely.prepared import prep
import cartopy.trace
try:
# https://github.com/pyproj4/pyproj/pull/912
from pyproj.crs import CustomConstructorCRS as _CRS
except ImportError:
from pyproj import CRS as _CRS
__document_these__ = ['CRS', 'Geocentric', 'Geodetic', 'Globe']
WGS84_SEMIMAJOR_AXIS = 6378137.0
WGS84_SEMIMINOR_AXIS = 6356752.3142
# Cache the transformer creation method
@lru_cache
def _get_transformer_from_crs(src_crs, tgt_crs):
return Transformer.from_crs(src_crs, tgt_crs, always_xy=True)
def _safe_pj_transform(src_crs, tgt_crs, x, y, z=None, trap=True):
transformer = _get_transformer_from_crs(src_crs, tgt_crs)
# if a projection is essentially 2d there
# should be no harm in setting its z to 0
if z is None:
z = np.zeros_like(x)
return transformer.transform(x, y, z, errcheck=trap)
class Globe:
"""
Define an ellipsoid and, optionally, how to relate it to the real world.
"""
def __init__(self, datum=None, ellipse='WGS84',
semimajor_axis=None, semiminor_axis=None,
flattening=None, inverse_flattening=None,
towgs84=None, nadgrids=None):
"""
Parameters
----------
datum
Proj "datum" definition. Defaults to None.
ellipse
Proj "ellps" definition. Defaults to 'WGS84'.
semimajor_axis
Semimajor axis of the spheroid / ellipsoid. Defaults to None.
semiminor_axis
Semiminor axis of the ellipsoid. Defaults to None.
flattening
Flattening of the ellipsoid. Defaults to None.
inverse_flattening
Inverse flattening of the ellipsoid. Defaults to None.
towgs84
Passed through to the Proj definition. Defaults to None.
nadgrids
Passed through to the Proj definition. Defaults to None.
"""
self.datum = datum
self.ellipse = ellipse
self.semimajor_axis = semimajor_axis
self.semiminor_axis = semiminor_axis
self.flattening = flattening
self.inverse_flattening = inverse_flattening
self.towgs84 = towgs84
self.nadgrids = nadgrids
def to_proj4_params(self):
"""
Create an OrderedDict of key value pairs which represents this globe
in terms of proj params.
"""
proj4_params = (
['datum', self.datum],
['ellps', self.ellipse],
['a', self.semimajor_axis],
['b', self.semiminor_axis],
['f', self.flattening],
['rf', self.inverse_flattening],
['towgs84', self.towgs84],
['nadgrids', self.nadgrids]
)
return OrderedDict((k, v) for k, v in proj4_params if v is not None)
class CRS(_CRS):
"""
Define a Coordinate Reference System using proj.
"""
#: Whether this projection can handle ellipses.
_handles_ellipses = True
def __init__(self, proj4_params, globe=None):
"""
Parameters
----------
proj4_params: iterable of key-value pairs
The proj4 parameters required to define the
desired CRS. The parameters should not describe
the desired elliptic model, instead create an
appropriate Globe instance. The ``proj4_params``
parameters will override any parameters that the
Globe defines.
globe: :class:`~cartopy.crs.Globe` instance, optional
If omitted, the default Globe instance will be created.
See :class:`~cartopy.crs.Globe` for details.
"""
# for compatibility with pyproj.CRS and rasterio.crs.CRS
try:
proj4_params = proj4_params.to_wkt()
except AttributeError:
pass
# handle PROJ JSON
if (
isinstance(proj4_params, dict) and
"proj" not in proj4_params and
"init" not in proj4_params
):
proj4_params = json.dumps(proj4_params)
if globe is not None and isinstance(proj4_params, str):
raise ValueError("Cannot have 'globe' with string params.")
if globe is None and not isinstance(proj4_params, str):
if self._handles_ellipses:
globe = Globe()
else:
globe = Globe(semimajor_axis=WGS84_SEMIMAJOR_AXIS,
ellipse=None)
if globe is not None and not self._handles_ellipses:
a = globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS
b = globe.semiminor_axis or a
if a != b or globe.ellipse is not None:
warnings.warn(f'The {self.__class__.__name__!r} projection '
'does not handle elliptical globes.')
self.globe = globe
if isinstance(proj4_params, str):
self._proj4_params = {}
self.proj4_init = proj4_params
else:
self._proj4_params = self.globe.to_proj4_params()
self._proj4_params.update(proj4_params)
init_items = []
for k, v in self._proj4_params.items():
if v is not None:
if isinstance(v, float):
init_items.append(f'+{k}={v:.16}')
elif isinstance(v, np.float32):
init_items.append(f'+{k}={v:.8}')
else:
init_items.append(f'+{k}={v}')
else:
init_items.append(f'+{k}')
self.proj4_init = ' '.join(init_items) + ' +no_defs'
super().__init__(self.proj4_init)
def __eq__(self, other):
if isinstance(other, CRS) and self.proj4_init == other.proj4_init:
# Fast path Cartopy's CRS
return True
# For everything else, we let pyproj handle the comparison
return super().__eq__(other)
def __hash__(self):
"""Hash the CRS based on its proj4_init string."""
return hash(self.proj4_init)
def __reduce__(self):
"""
Implement the __reduce__ API so that unpickling produces a stateless
instance of this class (e.g. an empty tuple). The state will then be
added via __getstate__ and __setstate__.
We are forced to this approach because a CRS does not store
the constructor keyword arguments in its state.
"""
return self.__class__, (), self.__getstate__()
def __getstate__(self):
"""Return the full state of this instance for reconstruction
in ``__setstate__``.
"""
state = self.__dict__.copy()
# remove pyproj specific attrs
state.pop('srs', None)
state.pop('_local', None)
# Remove the proj4 instance and the proj4_init string, which can
# be re-created (in __setstate__) from the other arguments.
state.pop('proj4', None)
state.pop('proj4_init', None)
state['proj4_params'] = self.proj4_params
return state
def __setstate__(self, state):
"""
Take the dictionary created by ``__getstate__`` and passes it
through to this implementation's __init__ method.
"""
# Strip out the key state items for a CRS instance.
CRS_state = {key: state.pop(key) for key in ['proj4_params', 'globe']}
# Put everything else directly into the dict of the instance.
self.__dict__.update(state)
# Call the init of this class to ensure that the projection is
# properly initialised with proj4.
CRS.__init__(self, **CRS_state)
def _as_mpl_transform(self, axes=None):
"""
Cast this CRS instance into a :class:`matplotlib.axes.Axes` using
the Matplotlib ``_as_mpl_transform`` interface.
"""
# lazy import mpl.geoaxes (and therefore matplotlib) as mpl
# is only an optional dependency
import cartopy.mpl.geoaxes as geoaxes
if not isinstance(axes, geoaxes.GeoAxes):
raise ValueError(
'Axes should be an instance of GeoAxes, got %s' % type(axes)
)
return (
geoaxes.InterProjectionTransform(self, axes.projection) +
axes.transData
)
@property
def proj4_params(self):
return dict(self._proj4_params)
def as_geocentric(self):
"""
Return a new Geocentric CRS with the same ellipse/datum as this
CRS.
"""
return CRS(
{
"$schema": (
"https://proj.org/schemas/v0.2/projjson.schema.json"
),
"type": "GeodeticCRS",
"name": "unknown",
"datum": self.datum.to_json_dict(),
"coordinate_system": {
"subtype": "Cartesian",
"axis": [
{
"name": "Geocentric X",
"abbreviation": "X",
"direction": "geocentricX",
"unit": "metre"
},
{
"name": "Geocentric Y",
"abbreviation": "Y",
"direction": "geocentricY",
"unit": "metre"
},
{
"name": "Geocentric Z",
"abbreviation": "Z",
"direction": "geocentricZ",
"unit": "metre"
}
]
}
}
)
def as_geodetic(self):
"""
Return a new Geodetic CRS with the same ellipse/datum as this
CRS.
"""
return CRS(self.geodetic_crs.srs)
def is_geodetic(self):
return self.is_geographic
def transform_point(self, x, y, src_crs, trap=True):
"""
transform_point(x, y, src_crs)
Transform the given float64 coordinate pair, in the given source
coordinate system (``src_crs``), to this coordinate system.
Parameters
----------
x
the x coordinate, in ``src_crs`` coordinates, to transform
y
the y coordinate, in ``src_crs`` coordinates, to transform
src_crs
instance of :class:`CRS` that represents the coordinate
system of ``x`` and ``y``.
trap
Whether proj errors for "latitude or longitude exceeded limits" and
"tolerance condition error" should be trapped.
Returns
-------
(x, y) in this coordinate system
"""
result = self.transform_points(
src_crs, np.array([x]), np.array([y]), trap=trap,
).reshape((1, 3))
return result[0, 0], result[0, 1]
def transform_points(self, src_crs, x, y, z=None, trap=False):
"""
transform_points(src_crs, x, y[, z])
Transform the given coordinates, in the given source
coordinate system (``src_crs``), to this coordinate system.
Parameters
----------
src_crs
instance of :class:`CRS` that represents the
coordinate system of ``x``, ``y`` and ``z``.
x
the x coordinates (array), in ``src_crs`` coordinates,
to transform. May be 1 or 2 dimensional.
y
the y coordinates (array), in ``src_crs`` coordinates,
to transform. Its shape must match that of x.
z: optional
the z coordinates (array), in ``src_crs`` coordinates, to
transform. Defaults to None.
If supplied, its shape must match that of x.
trap
Whether proj errors for "latitude or longitude exceeded limits" and
"tolerance condition error" should be trapped.
Returns
-------
Array of shape ``x.shape + (3, )`` in this coordinate system.
"""
result_shape = tuple(x.shape[i] for i in range(x.ndim)) + (3, )
if z is None:
if x.ndim > 2 or y.ndim > 2:
raise ValueError('x and y arrays must be 1 or 2 dimensional')
elif x.ndim != 1 or y.ndim != 1:
x, y = x.flatten(), y.flatten()
if x.shape[0] != y.shape[0]:
raise ValueError('x and y arrays must have the same length')
else:
if x.ndim > 2 or y.ndim > 2 or z.ndim > 2:
raise ValueError('x, y and z arrays must be 1 or 2 '
'dimensional')
elif x.ndim != 1 or y.ndim != 1 or z.ndim != 1:
x, y, z = x.flatten(), y.flatten(), z.flatten()
if not x.shape[0] == y.shape[0] == z.shape[0]:
raise ValueError('x, y, and z arrays must have the same '
'length')
npts = x.shape[0]
result = np.empty([npts, 3], dtype=np.double)
if npts:
if self == src_crs and (
isinstance(src_crs, _CylindricalProjection) or
self.is_geodetic()):
# convert from [0,360] to [-180,180]
x = np.array(x, copy=True)
to_180 = (x > 180) | (x < -180)
x[to_180] = (((x[to_180] + 180) % 360) - 180)
try:
result[:, 0], result[:, 1], result[:, 2] = \
_safe_pj_transform(src_crs, self, x, y, z, trap=trap)
except ProjError as err:
msg = str(err).lower()
if (
"latitude" in msg or
"longitude" in msg or
"outside of projection domain" in msg or
"tolerance condition error" in msg
):
result[:] = np.nan
else:
raise
if not trap:
result[np.isinf(result)] = np.nan
if len(result_shape) > 2:
return result.reshape(result_shape)
return result
def transform_vectors(self, src_proj, x, y, u, v):
"""
transform_vectors(src_proj, x, y, u, v)
Transform the given vector components, with coordinates in the
given source coordinate system (``src_proj``), to this coordinate
system. The vector components must be given relative to the
source projection's coordinate reference system (grid eastward and
grid northward).
Parameters
----------
src_proj
The :class:`CRS.Projection` that represents the coordinate system
the vectors are defined in.
x
The x coordinates of the vectors in the source projection.
y
The y coordinates of the vectors in the source projection.
u
The grid-eastward components of the vectors.
v
The grid-northward components of the vectors.
Note
----
x, y, u and v may be 1 or 2 dimensional, but must all have matching
shapes.
Returns
-------
ut, vt: The transformed vector components.
Note
----
The algorithm used to transform vectors is an approximation
rather than an exact transform, but the accuracy should be
good enough for visualization purposes.
"""
if not (x.shape == y.shape == u.shape == v.shape):
raise ValueError('x, y, u and v arrays must be the same shape')
if x.ndim not in (1, 2):
raise ValueError('x, y, u and v must be 1 or 2 dimensional')
# Transform the coordinates to the target projection.
proj_xyz = self.transform_points(src_proj, x, y)
target_x, target_y = proj_xyz[..., 0], proj_xyz[..., 1]
# Rotate the input vectors to the projection.
#
# 1: Find the magnitude and direction of the input vectors.
vector_magnitudes = (u**2 + v**2)**0.5
vector_angles = np.arctan2(v, u)
# 2: Find a point in the direction of the original vector that is
# a small distance away from the base point of the vector (near
# the poles the point may have to be in the opposite direction
# to be valid).
factor = 360000.
delta = (src_proj.x_limits[1] - src_proj.x_limits[0]) / factor
x_perturbations = delta * np.cos(vector_angles)
y_perturbations = delta * np.sin(vector_angles)
# 3: Handle points that are invalid. These come from picking a new
# point that is outside the domain of the CRS. The first step is
# to apply the native transform to the input coordinates to make
# sure they are in the appropriate range. Then detect all the
# coordinates where the perturbation takes the point out of the
# valid x-domain and fix them. After that do the same for points
# that are outside the valid y-domain, which may reintroduce some
# points outside of the valid x-domain
proj_xyz = src_proj.transform_points(src_proj, x, y)
source_x, source_y = proj_xyz[..., 0], proj_xyz[..., 1]
# Detect all the coordinates where the perturbation takes the point
# outside of the valid x-domain, and reverse the direction of the
# perturbation to fix this.
eps = 1e-9
invalid_x = np.logical_or(
source_x + x_perturbations < src_proj.x_limits[0] - eps,
source_x + x_perturbations > src_proj.x_limits[1] + eps)
if invalid_x.any():
x_perturbations[invalid_x] *= -1
y_perturbations[invalid_x] *= -1
# Do the same for coordinates where the perturbation takes the point
# outside of the valid y-domain. This may reintroduce some points
# that will be outside the x-domain when the perturbation is
# applied.
invalid_y = np.logical_or(
source_y + y_perturbations < src_proj.y_limits[0] - eps,
source_y + y_perturbations > src_proj.y_limits[1] + eps)
if invalid_y.any():
x_perturbations[invalid_y] *= -1
y_perturbations[invalid_y] *= -1
# Keep track of the points where the perturbation direction was
# reversed.
reversed_vectors = np.logical_xor(invalid_x, invalid_y)
# See if there were any points where we cannot reverse the direction
# of the perturbation to get the perturbed point within the valid
# domain of the projection, and issue a warning if there are.
problem_points = np.logical_or(
source_x + x_perturbations < src_proj.x_limits[0] - eps,
source_x + x_perturbations > src_proj.x_limits[1] + eps)
if problem_points.any():
warnings.warn('Some vectors at source domain corners '
'may not have been transformed correctly')
# 4: Transform this set of points to the projection coordinates and
# find the angle between the base point and the perturbed point
# in the projection coordinates (reversing the direction at any
# points where the original was reversed in step 3).
proj_xyz = self.transform_points(src_proj,
source_x + x_perturbations,
source_y + y_perturbations)
target_x_perturbed = proj_xyz[..., 0]
target_y_perturbed = proj_xyz[..., 1]
projected_angles = np.arctan2(target_y_perturbed - target_y,
target_x_perturbed - target_x)
if reversed_vectors.any():
projected_angles[reversed_vectors] += np.pi
# 5: Form the projected vector components, preserving the magnitude
# of the original vectors.
projected_u = vector_magnitudes * np.cos(projected_angles)
projected_v = vector_magnitudes * np.sin(projected_angles)
return projected_u, projected_v
class Geodetic(CRS):
"""
Define a latitude/longitude coordinate system with spherical topology,
geographical distance and coordinates are measured in degrees.
"""
def __init__(self, globe=None):
"""
Parameters
----------
globe: A :class:`cartopy.crs.Globe`, optional
Defaults to a "WGS84" datum.
"""
proj4_params = [('proj', 'lonlat')]
globe = globe or Globe(datum='WGS84')
super().__init__(proj4_params, globe)
# XXX Implement fwd such as Basemap's Geod.
# Would be used in the tissot example.
# Could come from https://geographiclib.sourceforge.io
class Geocentric(CRS):
"""
Define a Geocentric coordinate system, where x, y, z are Cartesian
coordinates from the center of the Earth.
"""
def __init__(self, globe=None):
"""
Parameters
----------
globe: A :class:`cartopy.crs.Globe`, optional
Defaults to a "WGS84" datum.
"""
proj4_params = [('proj', 'geocent')]
globe = globe or Globe(datum='WGS84')
super().__init__(proj4_params, globe)
class RotatedGeodetic(CRS):
"""
Define a rotated latitude/longitude coordinate system with spherical
topology and geographical distance.
Coordinates are measured in degrees.
The class uses proj to perform an ob_tran operation, using the
pole_longitude to set a lon_0 then performing two rotations based on
pole_latitude and central_rotated_longitude.
This is equivalent to setting the new pole to a location defined by
the pole_latitude and pole_longitude values in the GeogCRS defined by
globe, then rotating this new CRS about it's pole using the
central_rotated_longitude value.
"""
def __init__(self, pole_longitude, pole_latitude,
central_rotated_longitude=0.0, globe=None):
"""
Parameters
----------
pole_longitude
Pole longitude position, in unrotated degrees.
pole_latitude
Pole latitude position, in unrotated degrees.
central_rotated_longitude: optional
Longitude rotation about the new pole, in degrees. Defaults to 0.
globe: optional
A :class:`cartopy.crs.Globe`. Defaults to a "WGS84" datum.
"""
globe = globe or Globe(datum='WGS84')
proj4_params = [('proj', 'ob_tran'), ('o_proj', 'latlon'),
('o_lon_p', central_rotated_longitude),
('o_lat_p', pole_latitude),
('lon_0', 180 + pole_longitude),
('to_meter', math.radians(1) * (
globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS))]
super().__init__(proj4_params, globe=globe)
class Projection(CRS, metaclass=ABCMeta):
"""
Define a projected coordinate system with flat topology and Euclidean
distance.
"""
_method_map = {
'Point': '_project_point',
'LineString': '_project_line_string',
'LinearRing': '_project_linear_ring',
'Polygon': '_project_polygon',
'MultiPoint': '_project_multipoint',
'MultiLineString': '_project_multiline',
'MultiPolygon': '_project_multipolygon',
}
# Whether or not this projection can handle wrapped coordinates
_wrappable = False
def __init__(self, *args, **kwargs):
super().__init__(*args, **kwargs)
self.bounds = None
if self.area_of_use:
# Convert lat/lon bounds to projected bounds.
# Geographic area of the entire dataset referenced to WGS 84
# NB. We can't use a polygon transform at this stage because
# that relies on the existence of the map boundary... the very
# thing we're trying to work out! ;-)
x0 = self.area_of_use.west
x1 = self.area_of_use.east
y0 = self.area_of_use.south
y1 = self.area_of_use.north
lons = np.array([x0, x0, x1, x1])
lats = np.array([y0, y1, y1, y0])
points = self.transform_points(self.as_geodetic(), lons, lats)
x = points[:, 0]
y = points[:, 1]
self.bounds = (x.min(), x.max(), y.min(), y.max())
x0, x1, y0, y1 = self.bounds
self.threshold = min(x1 - x0, y1 - y0) / 100.
@property
def boundary(self):
if self.bounds is None:
raise NotImplementedError
x0, x1, y0, y1 = self.bounds
return sgeom.LineString([(x0, y0), (x0, y1), (x1, y1), (x1, y0),
(x0, y0)])
@property
def x_limits(self):
if self.bounds is None:
raise NotImplementedError
x0, x1, y0, y1 = self.bounds
return (x0, x1)
@property
def y_limits(self):
if self.bounds is None:
raise NotImplementedError
x0, x1, y0, y1 = self.bounds
return (y0, y1)
@property
def threshold(self):
return getattr(self, '_threshold', 0.5)
@threshold.setter
def threshold(self, t):
self._threshold = t
@property
def cw_boundary(self):
try:
boundary = self._cw_boundary
except AttributeError:
boundary = sgeom.LinearRing(self.boundary)
self._cw_boundary = boundary
return boundary
@property
def ccw_boundary(self):
try:
boundary = self._ccw_boundary
except AttributeError:
boundary = sgeom.LinearRing(self.boundary.coords[::-1])
self._ccw_boundary = boundary
return boundary
@property
def domain(self):
try:
domain = self._domain
except AttributeError:
domain = self._domain = sgeom.Polygon(self.boundary)
return domain
def is_geodetic(self):
return False
def _determine_longitude_bounds(self, central_longitude):
# In new proj, using exact limits will wrap-around, so subtract a
# small epsilon:
epsilon = 1e-10
minlon = -180 + central_longitude
maxlon = 180 + central_longitude
if central_longitude > 0:
maxlon -= epsilon
elif central_longitude < 0:
minlon += epsilon
return minlon, maxlon
def _repr_html_(self):
from html import escape
try:
# As matplotlib is not a core cartopy dependency, don't error
# if it's not available.
import matplotlib.pyplot as plt
except ImportError:
# We can't return an SVG of the CRS, so let Jupyter fall back to
# a default repr by returning None.
return None
# Produce a visual repr of the Projection instance.
fig, ax = plt.subplots(figsize=(5, 3),
subplot_kw={'projection': self})
ax.set_global()
ax.coastlines('auto')
ax.gridlines()
buf = io.StringIO()
fig.savefig(buf, format='svg', bbox_inches='tight')
plt.close(fig)
# "Rewind" the buffer to the start and return it as an svg string.
buf.seek(0)
svg = buf.read()
return f'{svg}<pre>{escape(object.__repr__(self))}</pre>'
def _as_mpl_axes(self):
import cartopy.mpl.geoaxes as geoaxes
return geoaxes.GeoAxes, {'projection': self}
def project_geometry(self, geometry, src_crs=None):
"""
Project the given geometry into this projection.
Parameters
----------
geometry
The geometry to (re-)project.
src_crs: optional
The source CRS. Defaults to None.
If src_crs is None, the source CRS is assumed to be a geodetic
version of the target CRS.
Returns
-------
geometry
The projected result (a shapely geometry).
"""
if src_crs is None:
src_crs = self.as_geodetic()
elif not isinstance(src_crs, CRS):
raise TypeError('Source CRS must be an instance of CRS'
' or one of its subclasses, or None.')
geom_type = geometry.geom_type
method_name = self._method_map.get(geom_type)
if not method_name:
raise ValueError(f'Unsupported geometry type {geom_type!r}')
return getattr(self, method_name)(geometry, src_crs)
def _project_point(self, point, src_crs):
return sgeom.Point(*self.transform_point(point.x, point.y, src_crs))
def _project_line_string(self, geometry, src_crs):
return cartopy.trace.project_linear(geometry, src_crs, self)
def _project_linear_ring(self, linear_ring, src_crs):
"""
Project the given LinearRing from the src_crs into this CRS and
returns a list of LinearRings and a single MultiLineString.
"""
debug = False
# 1) Resolve the initial lines into projected segments
# 1abc
# def23ghi
# jkl41
multi_line_string = cartopy.trace.project_linear(linear_ring,
src_crs, self)
# Threshold for whether a point is close enough to be the same
# point as another.
threshold = max(np.abs(self.x_limits + self.y_limits)) * 1e-5
# 2) Simplify the segments where appropriate.
if len(multi_line_string.geoms) > 1:
# Stitch together segments which are close to continuous.
# This is important when:
# 1) The first source point projects into the map and the
# ring has been cut by the boundary.
# Continuing the example from above this gives:
# def23ghi
# jkl41abc
# 2) The cut ends of segments are too close to reliably
# place into an order along the boundary.
line_strings = list(multi_line_string.geoms)
any_modified = False
i = 0
if debug:
first_coord = np.array([ls.coords[0] for ls in line_strings])
last_coord = np.array([ls.coords[-1] for ls in line_strings])
print('Distance matrix:')
np.set_printoptions(precision=2)
x = first_coord[:, np.newaxis, :]
y = last_coord[np.newaxis, :, :]
print(np.abs(x - y).max(axis=-1))
while i < len(line_strings):
modified = False
j = 0
while j < len(line_strings):
if i != j and np.allclose(line_strings[i].coords[0],
line_strings[j].coords[-1],
atol=threshold):
if debug:
print(f'Joining together {i} and {j}.')
last_coords = list(line_strings[j].coords)
first_coords = list(line_strings[i].coords)[1:]
combo = sgeom.LineString(last_coords + first_coords)
if j < i:
i, j = j, i
del line_strings[j], line_strings[i]
line_strings.append(combo)
modified = True
any_modified = True
break
else:
j += 1
if not modified:
i += 1
if any_modified:
multi_line_string = sgeom.MultiLineString(line_strings)
# 3) Check for rings that have been created by the projection stage.
rings = []
line_strings = []
for line in multi_line_string.geoms:
if len(line.coords) > 3 and np.allclose(line.coords[0],
line.coords[-1],
atol=threshold):
result_geometry = sgeom.LinearRing(line.coords[:-1])
rings.append(result_geometry)
else:
line_strings.append(line)
# If we found any rings, then we should re-create the multi-line str.
if rings:
multi_line_string = sgeom.MultiLineString(line_strings)
return rings, multi_line_string
def _project_multipoint(self, geometry, src_crs):
geoms = []
for geom in geometry.geoms:
geoms.append(self._project_point(geom, src_crs))
if geoms:
return sgeom.MultiPoint(geoms)
else:
return sgeom.MultiPoint()
def _project_multiline(self, geometry, src_crs):
geoms = []
for geom in geometry.geoms:
r = self._project_line_string(geom, src_crs)
if r:
geoms.extend(r.geoms)
if geoms:
return sgeom.MultiLineString(geoms)
else:
return []
def _project_multipolygon(self, geometry, src_crs):
geoms = []
for geom in geometry.geoms:
r = self._project_polygon(geom, src_crs)
if r:
geoms.extend(r.geoms)
if geoms:
result = sgeom.MultiPolygon(geoms)
else:
result = sgeom.MultiPolygon()
return result
def _project_polygon(self, polygon, src_crs):
"""
Return the projected polygon(s) derived from the given polygon.
"""
# Determine orientation of polygon.
# TODO: Consider checking the internal rings have the opposite
# orientation to the external rings?
if src_crs.is_geodetic():
is_ccw = True
else:
is_ccw = polygon.exterior.is_ccw
# Project the polygon exterior/interior rings.
# Each source ring will result in either a ring, or one or more
# lines.
rings = []
multi_lines = []
for src_ring in [polygon.exterior] + list(polygon.interiors):
p_rings, p_mline = self._project_linear_ring(src_ring, src_crs)
if p_rings:
rings.extend(p_rings)
if len(p_mline.geoms) > 0:
multi_lines.append(p_mline)
# Convert any lines to rings by attaching them to the boundary.
if multi_lines:
rings.extend(self._attach_lines_to_boundary(multi_lines, is_ccw))
# Resolve all the inside vs. outside rings, and convert to the
# final MultiPolygon.
return self._rings_to_multi_polygon(rings, is_ccw)
def _attach_lines_to_boundary(self, multi_line_strings, is_ccw):
"""
Return a list of LinearRings by attaching the ends of the given lines
to the boundary, paying attention to the traversal directions of the
lines and boundary.
"""
debug = False
debug_plot_edges = False
# Accumulate all the boundary and segment end points, along with
# their distance along the boundary.
edge_things = []
# Get the boundary as a LineString of the correct orientation
# so we can compute distances along it.
if is_ccw:
boundary = self.ccw_boundary
else:
boundary = self.cw_boundary
def boundary_distance(xy):
return boundary.project(sgeom.Point(*xy))
# Squash all the LineStrings into a single list.
line_strings = []
for multi_line_string in multi_line_strings:
line_strings.extend(multi_line_string.geoms)
# Record the positions of all the segment ends
for i, line_string in enumerate(line_strings):
first_dist = boundary_distance(line_string.coords[0])
thing = _BoundaryPoint(first_dist, False,
(i, 'first', line_string.coords[0]))
edge_things.append(thing)
last_dist = boundary_distance(line_string.coords[-1])
thing = _BoundaryPoint(last_dist, False,
(i, 'last', line_string.coords[-1]))
edge_things.append(thing)
# Record the positions of all the boundary vertices
for xy in boundary.coords[:-1]:
point = sgeom.Point(*xy)
dist = boundary.project(point)
thing = _BoundaryPoint(dist, True, point)
edge_things.append(thing)
if debug_plot_edges:
import matplotlib.pyplot as plt
current_fig = plt.gcf()
fig = plt.figure()
# Reset the current figure so we don't upset anything.
plt.figure(current_fig.number)
ax = fig.add_subplot(1, 1, 1)
# Order everything as if walking around the boundary.
# NB. We make line end-points take precedence over boundary points
# to ensure that end-points are still found and followed when they
# coincide.
edge_things.sort(key=lambda thing: (thing.distance, thing.kind))
remaining_ls = dict(enumerate(line_strings))
prev_thing = None
for edge_thing in edge_things[:]:
if (prev_thing is not None and
not edge_thing.kind and
not prev_thing.kind and
edge_thing.data[0] == prev_thing.data[0]):
j = edge_thing.data[0]
# Insert a edge boundary point in between this geometry.
mid_dist = (edge_thing.distance + prev_thing.distance) * 0.5
mid_point = boundary.interpolate(mid_dist)
new_thing = _BoundaryPoint(mid_dist, True, mid_point)
if debug:
print(f'Artificially insert boundary: {new_thing}')
ind = edge_things.index(edge_thing)
edge_things.insert(ind, new_thing)
prev_thing = None
else:
prev_thing = edge_thing
if debug:
print()
print('Edge things')
for thing in edge_things:
print(' ', thing)
if debug_plot_edges:
for thing in edge_things:
if isinstance(thing.data, sgeom.Point):
ax.plot(*thing.data.xy, marker='o')
else:
ax.plot(*thing.data[2], marker='o')
ls = line_strings[thing.data[0]]
coords = np.array(ls.coords)
ax.plot(coords[:, 0], coords[:, 1])
ax.text(coords[0, 0], coords[0, 1], thing.data[0])
ax.text(coords[-1, 0], coords[-1, 1],
f'{thing.data[0]}.')
def filter_last(t):
return t.kind or t.data[1] == 'first'
edge_things = list(filter(filter_last, edge_things))
processed_ls = []
while remaining_ls:
# Rename line_string to current_ls
i, current_ls = remaining_ls.popitem()
if debug:
import sys
sys.stdout.write('+')
sys.stdout.flush()
print()
print(f'Processing: {i}, {current_ls}')
added_linestring = set()
while True:
# Find out how far around this linestring's last
# point is on the boundary. We will use this to find
# the next point on the boundary.
d_last = boundary_distance(current_ls.coords[-1])
if debug:
print(f' d_last: {d_last!r}')
next_thing = _find_first_ge(edge_things, d_last)
# Remove this boundary point from the edge.
edge_things.remove(next_thing)
if debug:
print(' next_thing:', next_thing)
if next_thing.kind:
# We've just got a boundary point, add it, and keep going.
if debug:
print(' adding boundary point')
boundary_point = next_thing.data
combined_coords = (list(current_ls.coords) +
[(boundary_point.x, boundary_point.y)])
current_ls = sgeom.LineString(combined_coords)
elif next_thing.data[0] == i:
# We've gone all the way around and are now back at the
# first boundary thing.
if debug:
print(' close loop')
processed_ls.append(current_ls)
if debug_plot_edges:
coords = np.array(current_ls.coords)
ax.plot(coords[:, 0], coords[:, 1], color='black',
linestyle='--')
break
else:
if debug:
print(' adding line')
j = next_thing.data[0]
line_to_append = line_strings[j]
if j in remaining_ls:
remaining_ls.pop(j)
coords_to_append = list(line_to_append.coords)
# Build up the linestring.
current_ls = sgeom.LineString(list(current_ls.coords) +
coords_to_append)
# Catch getting stuck in an infinite loop by checking that
# linestring only added once.
if j not in added_linestring:
added_linestring.add(j)
else:
if debug_plot_edges:
plt.show()
raise RuntimeError('Unidentified problem with '
'geometry, linestring being '
're-added. Please raise an issue.')
# filter out any non-valid linear rings
def makes_valid_ring(line_string):
if len(line_string.coords) == 3:
# When sgeom.LinearRing is passed a LineString of length 3,
# if the first and last coordinate are equal, a LinearRing
# with 3 coordinates will be created. This object will cause
# a segfault when evaluated.
coords = list(line_string.coords)
return coords[0] != coords[-1] and line_string.is_valid
else:
return len(line_string.coords) > 3 and line_string.is_valid
linear_rings = [
sgeom.LinearRing(line_string)
for line_string in processed_ls
if makes_valid_ring(line_string)]
if debug:
print(' DONE')
return linear_rings
def _rings_to_multi_polygon(self, rings, is_ccw):
exterior_rings = []
interior_rings = []
for ring in rings:
if ring.is_ccw != is_ccw:
interior_rings.append(ring)
else:
exterior_rings.append(ring)
polygon_bits = []
# Turn all the exterior rings into polygon definitions,
# "slurping up" any interior rings they contain.
for exterior_ring in exterior_rings:
polygon = sgeom.Polygon(exterior_ring)
prep_polygon = prep(polygon)
holes = []
for interior_ring in interior_rings[:]:
if prep_polygon.contains(interior_ring):
holes.append(interior_ring)
interior_rings.remove(interior_ring)
elif polygon.crosses(interior_ring):
# Likely that we have an invalid geometry such as
# that from #509 or #537.
holes.append(interior_ring)
interior_rings.remove(interior_ring)
polygon_bits.append((exterior_ring.coords,
[ring.coords for ring in holes]))
# Any left over "interior" rings need "inverting" with respect
# to the boundary.
if interior_rings:
boundary_poly = self.domain
x3, y3, x4, y4 = boundary_poly.bounds
bx = (x4 - x3) * 0.1
by = (y4 - y3) * 0.1
x3 -= bx
y3 -= by
x4 += bx
y4 += by
for ring in interior_rings:
# Use shapely buffer in an attempt to fix invalid geometries
polygon = sgeom.Polygon(ring).buffer(0)
if not polygon.is_empty and polygon.is_valid:
x1, y1, x2, y2 = polygon.bounds
bx = (x2 - x1) * 0.1
by = (y2 - y1) * 0.1
x1 -= bx
y1 -= by
x2 += bx
y2 += by
box = sgeom.box(min(x1, x3), min(y1, y3),
max(x2, x4), max(y2, y4))
# Invert the polygon
polygon = box.difference(polygon)
# Intersect the inverted polygon with the boundary
polygon = boundary_poly.intersection(polygon)
if not polygon.is_empty:
polygon_bits.append(polygon)
if polygon_bits:
multi_poly = sgeom.MultiPolygon(polygon_bits)
else:
multi_poly = sgeom.MultiPolygon()
return multi_poly
def quick_vertices_transform(self, vertices, src_crs):
"""
Where possible, return a vertices array transformed to this CRS from
the given vertices array of shape ``(n, 2)`` and the source CRS.
Note
----
This method may return None to indicate that the vertices cannot
be transformed quickly, and a more complex geometry transformation
is required (see :meth:`cartopy.crs.Projection.project_geometry`).
"""
return_value = None
if self == src_crs:
x = vertices[:, 0]
y = vertices[:, 1]
# Extend the limits a tiny amount to allow for precision mistakes
epsilon = 1.e-10
x_limits = (self.x_limits[0] - epsilon, self.x_limits[1] + epsilon)
y_limits = (self.y_limits[0] - epsilon, self.y_limits[1] + epsilon)
if (x.min() >= x_limits[0] and x.max() <= x_limits[1] and
y.min() >= y_limits[0] and y.max() <= y_limits[1]):
return_value = vertices
return return_value
class _RectangularProjection(Projection, metaclass=ABCMeta):
"""
The abstract superclass of projections with a rectangular domain which
is symmetric about the origin.
"""
_wrappable = True
def __init__(self, proj4_params, half_width, half_height, globe=None):
self._half_width = half_width
self._half_height = half_height
super().__init__(proj4_params, globe=globe)
@property
def boundary(self):
w, h = self._half_width, self._half_height
return sgeom.LinearRing([(-w, -h), (-w, h), (w, h), (w, -h), (-w, -h)])
@property
def x_limits(self):
return (-self._half_width, self._half_width)
@property
def y_limits(self):
return (-self._half_height, self._half_height)
class _CylindricalProjection(_RectangularProjection, metaclass=ABCMeta):
"""
The abstract class which denotes cylindrical projections where we
want to allow x values to wrap around.
"""
_wrappable = True
def _ellipse_boundary(semimajor=2, semiminor=1, easting=0, northing=0, n=201):
"""
Define a projection boundary using an ellipse.
This type of boundary is used by several projections.
"""
t = np.linspace(0, -2 * np.pi, n) # Clockwise boundary.
coords = np.vstack([semimajor * np.cos(t), semiminor * np.sin(t)])
coords += ([easting], [northing])
return coords
class PlateCarree(_CylindricalProjection):
def __init__(self, central_longitude=0.0, globe=None):
globe = globe or Globe(semimajor_axis=WGS84_SEMIMAJOR_AXIS)
proj4_params = [('proj', 'eqc'), ('lon_0', central_longitude),
('to_meter', math.radians(1) * (
globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS)),
('vto_meter', 1)]
x_max = 180
y_max = 90
# Set the threshold around 0.5 if the x max is 180.
self.threshold = x_max / 360
super().__init__(proj4_params, x_max, y_max, globe=globe)
def _bbox_and_offset(self, other_plate_carree):
"""
Return a pair of (xmin, xmax) pairs and an offset which can be used
for identification of whether data in ``other_plate_carree`` needs
to be transformed to wrap appropriately.
>>> import cartopy.crs as ccrs
>>> src = ccrs.PlateCarree(central_longitude=10)
>>> bboxes, offset = ccrs.PlateCarree()._bbox_and_offset(src)
>>> print(bboxes)
[[-180, -170.0], [-170.0, 180]]
>>> print(offset)
10.0
The returned values are longitudes in ``other_plate_carree``'s
coordinate system.
Warning
-------
The two CRSs must be identical in every way, other than their
central longitudes. No checking of this is done.
"""
self_lon_0 = self.proj4_params['lon_0']
other_lon_0 = other_plate_carree.proj4_params['lon_0']
lon_0_offset = other_lon_0 - self_lon_0
lon_lower_bound_0 = self.x_limits[0]
lon_lower_bound_1 = (other_plate_carree.x_limits[0] + lon_0_offset)
if lon_lower_bound_1 < self.x_limits[0]:
lon_lower_bound_1 += np.diff(self.x_limits)[0]
lon_lower_bound_0, lon_lower_bound_1 = sorted(
[lon_lower_bound_0, lon_lower_bound_1])
bbox = [[lon_lower_bound_0, lon_lower_bound_1],
[lon_lower_bound_1, lon_lower_bound_0]]
bbox[1][1] += np.diff(self.x_limits)[0]
return bbox, lon_0_offset
def quick_vertices_transform(self, vertices, src_crs):
return_value = super().quick_vertices_transform(vertices, src_crs)
# Optimise the PlateCarree -> PlateCarree case where no
# wrapping or interpolation needs to take place.
if return_value is None and isinstance(src_crs, PlateCarree):
self_params = self.proj4_params.copy()
src_params = src_crs.proj4_params.copy()
self_params.pop('lon_0'), src_params.pop('lon_0')
xs, ys = vertices[:, 0], vertices[:, 1]
potential = (self_params == src_params and
self.y_limits[0] <= ys.min() and
self.y_limits[1] >= ys.max())
if potential:
mod = np.diff(src_crs.x_limits)[0]
bboxes, proj_offset = self._bbox_and_offset(src_crs)
x_lim = xs.min(), xs.max()
for poly in bboxes:
# Arbitrarily choose the number of moduli to look
# above and below the -180->180 range. If data is beyond
# this range, we're not going to transform it quickly.
for i in [-1, 0, 1, 2]:
offset = mod * i - proj_offset
if ((poly[0] + offset) <= x_lim[0] and
(poly[1] + offset) >= x_lim[1]):
return_value = vertices + [[-offset, 0]]
break
if return_value is not None:
break
return return_value
class TransverseMercator(Projection):
"""
A Transverse Mercator projection.
"""
_wrappable = True
def __init__(self, central_longitude=0.0, central_latitude=0.0,
false_easting=0.0, false_northing=0.0,
scale_factor=1.0, globe=None, approx=False):
"""
Parameters
----------
central_longitude: optional
The true longitude of the central meridian in degrees.
Defaults to 0.
central_latitude: optional
The true latitude of the planar origin in degrees. Defaults to 0.
false_easting: optional
X offset from the planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from the planar origin in metres. Defaults to 0.
scale_factor: optional
Scale factor at the central meridian. Defaults to 1.
globe: optional
An instance of :class:`cartopy.crs.Globe`. If omitted, a default
globe is created.
approx: optional
Whether to use Proj's approximate projection (True), or the new
Extended Transverse Mercator code (False). Defaults to True, but
will change to False in the next release.
"""
proj4_params = [('proj', 'tmerc'), ('lon_0', central_longitude),
('lat_0', central_latitude), ('k', scale_factor),
('x_0', false_easting), ('y_0', false_northing),
('units', 'm')]
if approx:
proj4_params += [('approx', None)]
super().__init__(proj4_params, globe=globe)
self.threshold = 1e4
@property
def boundary(self):
x0, x1 = self.x_limits
y0, y1 = self.y_limits
return sgeom.LinearRing([(x0, y0), (x0, y1),
(x1, y1), (x1, y0),
(x0, y0)])
@property
def x_limits(self):
return (-2e7, 2e7)
@property
def y_limits(self):
return (-1e7, 1e7)
class OSGB(TransverseMercator):
def __init__(self, approx=False):
super().__init__(central_longitude=-2, central_latitude=49,
scale_factor=0.9996012717,
false_easting=400000, false_northing=-100000,
globe=Globe(datum='OSGB36', ellipse='airy'),
approx=approx)
@property
def boundary(self):
w = self.x_limits[1] - self.x_limits[0]
h = self.y_limits[1] - self.y_limits[0]
return sgeom.LinearRing([(0, 0), (0, h), (w, h), (w, 0), (0, 0)])
@property
def x_limits(self):
return (0, 7e5)
@property
def y_limits(self):
return (0, 13e5)
class OSNI(TransverseMercator):
def __init__(self, approx=False):
globe = Globe(semimajor_axis=6377340.189,
semiminor_axis=6356034.447938534)
super().__init__(central_longitude=-8, central_latitude=53.5,
scale_factor=1.000035,
false_easting=200000, false_northing=250000,
globe=globe, approx=approx)
@property
def boundary(self):
w = self.x_limits[1] - self.x_limits[0]
h = self.y_limits[1] - self.y_limits[0]
return sgeom.LinearRing([(0, 0), (0, h), (w, h), (w, 0), (0, 0)])
@property
def x_limits(self):
return (18814.9667, 386062.3293)
@property
def y_limits(self):
return (11764.8481, 464720.9559)
class UTM(Projection):
"""
Universal Transverse Mercator projection.
"""
def __init__(self, zone, southern_hemisphere=False, globe=None):
"""
Parameters
----------
zone
The numeric zone of the UTM required.
southern_hemisphere: optional
Set to True if the zone is in the southern hemisphere. Defaults to
False.
globe: optional
An instance of :class:`cartopy.crs.Globe`. If omitted, a default
globe is created.
"""
proj4_params = [('proj', 'utm'),
('units', 'm'),
('zone', zone)]
if southern_hemisphere:
proj4_params.append(('south', None))
super().__init__(proj4_params, globe=globe)
self.threshold = 1e2
@property
def boundary(self):
x0, x1 = self.x_limits
y0, y1 = self.y_limits
return sgeom.LinearRing([(x0, y0), (x0, y1),
(x1, y1), (x1, y0),
(x0, y0)])
@property
def x_limits(self):
easting = 5e5
# allow 50% overflow
return (0 - easting / 2, 2 * easting + easting / 2)
@property
def y_limits(self):
northing = 1e7
# allow 50% overflow
return (0 - northing, 2 * northing + northing / 2)
class EuroPP(UTM):
"""
UTM Zone 32 projection for EuroPP domain.
Ellipsoid is International 1924, Datum is ED50.
"""
def __init__(self):
globe = Globe(ellipse='intl')
super().__init__(32, globe=globe)
@property
def x_limits(self):
return (-1.4e6, 2e6)
@property
def y_limits(self):
return (4e6, 7.9e6)
class Mercator(Projection):
"""
A Mercator projection.
"""
_wrappable = True
def __init__(self, central_longitude=0.0,
min_latitude=-80.0, max_latitude=84.0,
globe=None, latitude_true_scale=None,
false_easting=0.0, false_northing=0.0, scale_factor=None):
"""
Parameters
----------
central_longitude: optional
The central longitude. Defaults to 0.
min_latitude: optional
The maximum southerly extent of the projection. Defaults
to -80 degrees.
max_latitude: optional
The maximum northerly extent of the projection. Defaults
to 84 degrees.
globe: A :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
latitude_true_scale: optional
The latitude where the scale is 1. Defaults to 0 degrees.
false_easting: optional
X offset from the planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from the planar origin in metres. Defaults to 0.
scale_factor: optional
Scale factor at natural origin. Defaults to unused.
Notes
-----
Only one of ``latitude_true_scale`` and ``scale_factor`` should
be included.
"""
proj4_params = [('proj', 'merc'),
('lon_0', central_longitude),
('x_0', false_easting),
('y_0', false_northing),
('units', 'm')]
# If it's None, we don't pass it to Proj4, in which case its default
# of 0.0 will be used.
if latitude_true_scale is not None:
proj4_params.append(('lat_ts', latitude_true_scale))
if scale_factor is not None:
if latitude_true_scale is not None:
raise ValueError('It does not make sense to provide both '
'"scale_factor" and "latitude_true_scale". ')
else:
proj4_params.append(('k_0', scale_factor))
super().__init__(proj4_params, globe=globe)
# Need to have x/y limits defined for the initial hash which
# gets used within transform_points for caching
self._x_limits = self._y_limits = None
# Calculate limits.
minlon, maxlon = self._determine_longitude_bounds(central_longitude)
limits = self.transform_points(self.as_geodetic(),
np.array([minlon, maxlon]),
np.array([min_latitude, max_latitude]))
self._x_limits = tuple(limits[..., 0])
self._y_limits = tuple(limits[..., 1])
self.threshold = min(np.diff(self.x_limits)[0] / 720,
np.diff(self.y_limits)[0] / 360)
def __eq__(self, other):
res = super().__eq__(other)
if hasattr(other, "_y_limits") and hasattr(other, "_x_limits"):
res = res and self._y_limits == other._y_limits and \
self._x_limits == other._x_limits
return res
def __ne__(self, other):
return not self == other
def __hash__(self):
return hash((self.proj4_init, self._x_limits, self._y_limits))
@property
def boundary(self):
x0, x1 = self.x_limits
y0, y1 = self.y_limits
return sgeom.LinearRing([(x0, y0), (x0, y1),
(x1, y1), (x1, y0),
(x0, y0)])
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
# Define a specific instance of a Mercator projection, the Google mercator.
Mercator.GOOGLE = Mercator(min_latitude=-85.0511287798066,
max_latitude=85.0511287798066,
globe=Globe(ellipse=None,
semimajor_axis=WGS84_SEMIMAJOR_AXIS,
semiminor_axis=WGS84_SEMIMAJOR_AXIS,
nadgrids='@null'))
# Deprecated form
GOOGLE_MERCATOR = Mercator.GOOGLE
class LambertCylindrical(_RectangularProjection):
def __init__(self, central_longitude=0.0, globe=None):
globe = globe or Globe(semimajor_axis=WGS84_SEMIMAJOR_AXIS)
proj4_params = [('proj', 'cea'), ('lon_0', central_longitude),
('to_meter', math.radians(1) * (
globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS))]
super().__init__(proj4_params, 180, math.degrees(1), globe=globe)
class LambertConformal(Projection):
"""
A Lambert Conformal conic projection.
"""
def __init__(self, central_longitude=-96.0, central_latitude=39.0,
false_easting=0.0, false_northing=0.0,
standard_parallels=(33, 45),
globe=None, cutoff=-30):
"""
Parameters
----------
central_longitude: optional
The central longitude. Defaults to -96.
central_latitude: optional
The central latitude. Defaults to 39.
false_easting: optional
X offset from planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from planar origin in metres. Defaults to 0.
standard_parallels: optional
Standard parallel latitude(s). Defaults to (33, 45).
globe: optional
A :class:`cartopy.crs.Globe`. If omitted, a default globe is
created.
cutoff: optional
Latitude of map cutoff.
The map extends to infinity opposite the central pole
so we must cut off the map drawing before then.
A value of 0 will draw half the globe. Defaults to -30.
"""
proj4_params = [('proj', 'lcc'),
('lon_0', central_longitude),
('lat_0', central_latitude),
('x_0', false_easting),
('y_0', false_northing)]
n_parallels = len(standard_parallels)
if not 1 <= n_parallels <= 2:
raise ValueError('1 or 2 standard parallels must be specified. '
f'Got {n_parallels} ({standard_parallels})')
proj4_params.append(('lat_1', standard_parallels[0]))
if n_parallels == 2:
proj4_params.append(('lat_2', standard_parallels[1]))
super().__init__(proj4_params, globe=globe)
# Compute whether this projection is at the "north pole" or the
# "south pole" (after the central lon/lat have been taken into
# account).
if n_parallels == 1:
plat = 90 if standard_parallels[0] > 0 else -90
else:
# Which pole are the parallels closest to? That is the direction
# that the cone converges.
if abs(standard_parallels[0]) > abs(standard_parallels[1]):
poliest_sec = standard_parallels[0]
else:
poliest_sec = standard_parallels[1]
plat = 90 if poliest_sec > 0 else -90
self.cutoff = cutoff
n = 91
lons = np.empty(n + 2)
lats = np.full(n + 2, float(cutoff))
lons[0] = lons[-1] = 0
lats[0] = lats[-1] = plat
if plat == 90:
# Ensure clockwise
lons[1:-1] = np.linspace(central_longitude + 180 - 0.001,
central_longitude - 180 + 0.001, n)
else:
lons[1:-1] = np.linspace(central_longitude - 180 + 0.001,
central_longitude + 180 - 0.001, n)
points = self.transform_points(PlateCarree(), lons, lats)
self._boundary = sgeom.LinearRing(points)
mins = np.min(points, axis=0)
maxs = np.max(points, axis=0)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = 1e5
def __eq__(self, other):
res = super().__eq__(other)
if hasattr(other, "cutoff"):
res = res and self.cutoff == other.cutoff
return res
def __ne__(self, other):
return not self == other
def __hash__(self):
return hash((self.proj4_init, self.cutoff))
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class LambertAzimuthalEqualArea(Projection):
"""
A Lambert Azimuthal Equal-Area projection.
"""
_wrappable = True
def __init__(self, central_longitude=0.0, central_latitude=0.0,
false_easting=0.0, false_northing=0.0,
globe=None):
"""
Parameters
----------
central_longitude: optional
The central longitude. Defaults to 0.
central_latitude: optional
The central latitude. Defaults to 0.
false_easting: optional
X offset from planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from planar origin in metres. Defaults to 0.
globe: optional
A :class:`cartopy.crs.Globe`. If omitted, a default globe is
created.
"""
proj4_params = [('proj', 'laea'),
('lon_0', central_longitude),
('lat_0', central_latitude),
('x_0', false_easting),
('y_0', false_northing)]
super().__init__(proj4_params, globe=globe)
a = float(self.ellipsoid.semi_major_metre or WGS84_SEMIMAJOR_AXIS)
# Find the antipode, and shift it a small amount in latitude to
# approximate the extent of the projection:
lon = central_longitude + 180
sign = np.sign(central_latitude) or 1
lat = -central_latitude + sign * 0.01
x, max_y = self.transform_point(lon, lat, PlateCarree(globe=globe))
coords = _ellipse_boundary(a * 1.9999, max_y - false_northing,
false_easting, false_northing, 61)
self._boundary = sgeom.polygon.LinearRing(coords.T)
mins = np.min(coords, axis=1)
maxs = np.max(coords, axis=1)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = np.diff(self._x_limits)[0] * 1e-3
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class Miller(_RectangularProjection):
_handles_ellipses = False
def __init__(self, central_longitude=0.0, globe=None):
if globe is None:
globe = Globe(semimajor_axis=WGS84_SEMIMAJOR_AXIS, ellipse=None)
a = globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS
proj4_params = [('proj', 'mill'), ('lon_0', central_longitude)]
# See Snyder, 1987. Eqs (11-1) and (11-2) substituting maximums of
# (lambda-lambda0)=180 and phi=90 to get limits.
super().__init__(proj4_params, a * np.pi, a * 2.303412543376391,
globe=globe)
class RotatedPole(_CylindricalProjection):
"""
A rotated latitude/longitude projected coordinate system
with cylindrical topology and projected distance.
Coordinates are measured in projection metres.
The class uses proj to perform an ob_tran operation, using the
pole_longitude to set a lon_0 then performing two rotations based on
pole_latitude and central_rotated_longitude.
This is equivalent to setting the new pole to a location defined by
the pole_latitude and pole_longitude values in the GeogCRS defined by
globe, then rotating this new CRS about it's pole using the
central_rotated_longitude value.
"""
def __init__(self, pole_longitude=0.0, pole_latitude=90.0,
central_rotated_longitude=0.0, globe=None):
"""
Parameters
----------
pole_longitude: optional
Pole longitude position, in unrotated degrees. Defaults to 0.
pole_latitude: optional
Pole latitude position, in unrotated degrees. Defaults to 0.
central_rotated_longitude: optional
Longitude rotation about the new pole, in degrees. Defaults to 0.
globe: optional
An optional :class:`cartopy.crs.Globe`. Defaults to a "WGS84"
datum.
"""
globe = globe or Globe(semimajor_axis=WGS84_SEMIMAJOR_AXIS)
proj4_params = [('proj', 'ob_tran'), ('o_proj', 'latlon'),
('o_lon_p', central_rotated_longitude),
('o_lat_p', pole_latitude),
('lon_0', 180 + pole_longitude),
('to_meter', math.radians(1) * (
globe.semimajor_axis or WGS84_SEMIMAJOR_AXIS))]
super().__init__(proj4_params, 180, 90, globe=globe)
class Gnomonic(Projection):
_handles_ellipses = False
def __init__(self, central_latitude=0.0,
central_longitude=0.0, globe=None):
proj4_params = [('proj', 'gnom'), ('lat_0', central_latitude),
('lon_0', central_longitude)]
super().__init__(proj4_params, globe=globe)
self._max = 5e7
self.threshold = 1e5
@property
def boundary(self):
return sgeom.Point(0, 0).buffer(self._max).exterior
@property
def x_limits(self):
return (-self._max, self._max)
@property
def y_limits(self):
return (-self._max, self._max)
class Stereographic(Projection):
_wrappable = True
def __init__(self, central_latitude=0.0, central_longitude=0.0,
false_easting=0.0, false_northing=0.0,
true_scale_latitude=None,
scale_factor=None, globe=None):
proj4_params = [('proj', 'stere'), ('lat_0', central_latitude),
('lon_0', central_longitude),
('x_0', false_easting), ('y_0', false_northing)]
if true_scale_latitude is not None:
if central_latitude not in (-90., 90.):
warnings.warn('"true_scale_latitude" parameter is only used '
'for polar stereographic projections. Consider '
'the use of "scale_factor" instead.',
stacklevel=2)
proj4_params.append(('lat_ts', true_scale_latitude))
if scale_factor is not None:
if true_scale_latitude is not None:
raise ValueError('It does not make sense to provide both '
'"scale_factor" and "true_scale_latitude". '
'Ignoring "scale_factor".')
else:
proj4_params.append(('k_0', scale_factor))
super().__init__(proj4_params, globe=globe)
# TODO: Let the globe return the semimajor axis always.
a = float(self.ellipsoid.semi_major_metre or WGS84_SEMIMAJOR_AXIS)
b = float(self.ellipsoid.semi_minor_metre or WGS84_SEMIMINOR_AXIS)
# Note: The magic number has been picked to maintain consistent
# behaviour with a wgs84 globe. There is no guarantee that the scaling
# should even be linear.
x_axis_offset = 5e7 / WGS84_SEMIMAJOR_AXIS
y_axis_offset = 5e7 / WGS84_SEMIMINOR_AXIS
self._x_limits = (-a * x_axis_offset + false_easting,
a * x_axis_offset + false_easting)
self._y_limits = (-b * y_axis_offset + false_northing,
b * y_axis_offset + false_northing)
coords = _ellipse_boundary(self._x_limits[1], self._y_limits[1],
false_easting, false_northing, 91)
self._boundary = sgeom.LinearRing(coords.T)
self.threshold = np.diff(self._x_limits)[0] * 1e-3
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class NorthPolarStereo(Stereographic):
def __init__(self, central_longitude=0.0, true_scale_latitude=None,
globe=None):
super().__init__(
central_latitude=90,
central_longitude=central_longitude,
true_scale_latitude=true_scale_latitude, # None is +90
globe=globe)
class SouthPolarStereo(Stereographic):
def __init__(self, central_longitude=0.0, true_scale_latitude=None,
globe=None):
super().__init__(
central_latitude=-90,
central_longitude=central_longitude,
true_scale_latitude=true_scale_latitude, # None is -90
globe=globe)
class Orthographic(Projection):
_handles_ellipses = False
def __init__(self, central_longitude=0.0, central_latitude=0.0,
globe=None):
proj4_params = [('proj', 'ortho'), ('lon_0', central_longitude),
('lat_0', central_latitude)]
super().__init__(proj4_params, globe=globe)
# TODO: Let the globe return the semimajor axis always.
a = float(self.ellipsoid.semi_major_metre or WGS84_SEMIMAJOR_AXIS)
# To stabilise the projection of geometries, we reduce the boundary by
# a tiny fraction at the cost of the extreme edges.
coords = _ellipse_boundary(a * 0.99999, a * 0.99999, n=61)
self._boundary = sgeom.polygon.LinearRing(coords.T)
mins = np.min(coords, axis=1)
maxs = np.max(coords, axis=1)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = np.diff(self._x_limits)[0] * 0.02
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class _WarpedRectangularProjection(Projection, metaclass=ABCMeta):
_wrappable = True
def __init__(self, proj4_params, central_longitude,
false_easting=None, false_northing=None, globe=None):
if false_easting is not None:
proj4_params += [('x_0', false_easting)]
if false_northing is not None:
proj4_params += [('y_0', false_northing)]
super().__init__(proj4_params, globe=globe)
# Obtain boundary points
minlon, maxlon = self._determine_longitude_bounds(central_longitude)
n = 91
lon = np.empty(2 * n + 1)
lat = np.empty(2 * n + 1)
lon[:n] = minlon
lat[:n] = np.linspace(-90, 90, n)
lon[n:2 * n] = maxlon
lat[n:2 * n] = np.linspace(90, -90, n)
lon[-1] = minlon
lat[-1] = -90
points = self.transform_points(self.as_geodetic(), lon, lat)
self._boundary = sgeom.LinearRing(points)
mins = np.min(points, axis=0)
maxs = np.max(points, axis=0)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class _Eckert(_WarpedRectangularProjection, metaclass=ABCMeta):
"""
An Eckert projection.
This class implements all the methods common to the Eckert family of
projections.
"""
_handles_ellipses = False
def __init__(self, central_longitude=0, false_easting=None,
false_northing=None, globe=None):
"""
Parameters
----------
central_longitude: float, optional
The central longitude. Defaults to 0.
false_easting: float, optional
X offset from planar origin in metres. Defaults to 0.
false_northing: float, optional
Y offset from planar origin in metres. Defaults to 0.
globe: :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
.. note::
This projection does not handle elliptical globes.
"""
proj4_params = [('proj', self._proj_name),
('lon_0', central_longitude)]
super().__init__(proj4_params, central_longitude,
false_easting=false_easting,
false_northing=false_northing,
globe=globe)
self.threshold = 1e5
class EckertI(_Eckert):
"""
An Eckert I projection.
This projection is pseudocylindrical, but not equal-area. Both meridians
and parallels are straight lines. Its equal-area pair is :class:`EckertII`.
"""
_proj_name = 'eck1'
class EckertII(_Eckert):
"""
An Eckert II projection.
This projection is pseudocylindrical, and equal-area. Both meridians and
parallels are straight lines. Its non-equal-area pair with equally-spaced
parallels is :class:`EckertI`.
"""
_proj_name = 'eck2'
class EckertIII(_Eckert):
"""
An Eckert III projection.
This projection is pseudocylindrical, but not equal-area. Parallels are
equally-spaced straight lines, while meridians are elliptical arcs up to
semicircles on the edges. Its equal-area pair is :class:`EckertIV`.
"""
_proj_name = 'eck3'
class EckertIV(_Eckert):
"""
An Eckert IV projection.
This projection is pseudocylindrical, and equal-area. Parallels are
unequally-spaced straight lines, while meridians are elliptical arcs up to
semicircles on the edges. Its non-equal-area pair with equally-spaced
parallels is :class:`EckertIII`.
It is commonly used for world maps.
"""
_proj_name = 'eck4'
class EckertV(_Eckert):
"""
An Eckert V projection.
This projection is pseudocylindrical, but not equal-area. Parallels are
equally-spaced straight lines, while meridians are sinusoidal arcs. Its
equal-area pair is :class:`EckertVI`.
"""
_proj_name = 'eck5'
class EckertVI(_Eckert):
"""
An Eckert VI projection.
This projection is pseudocylindrical, and equal-area. Parallels are
unequally-spaced straight lines, while meridians are sinusoidal arcs. Its
non-equal-area pair with equally-spaced parallels is :class:`EckertV`.
It is commonly used for world maps.
"""
_proj_name = 'eck6'
class EqualEarth(_WarpedRectangularProjection):
"""
An Equal Earth projection.
This projection is pseudocylindrical, and equal area. Parallels are
unequally-spaced straight lines, while meridians are equally-spaced arcs.
It is intended for world maps.
Note
----
To use this projection, you must be using Proj 5.2.0 or newer.
References
----------
Bojan Šavrič, Tom Patterson & Bernhard Jenny (2018)
The Equal Earth map projection,
International Journal of Geographical Information Science,
DOI: 10.1080/13658816.2018.1504949
"""
def __init__(self, central_longitude=0, false_easting=None,
false_northing=None, globe=None):
"""
Parameters
----------
central_longitude: float, optional
The central longitude. Defaults to 0.
false_easting: float, optional
X offset from planar origin in metres. Defaults to 0.
false_northing: float, optional
Y offset from planar origin in metres. Defaults to 0.
globe: :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
"""
proj_params = [('proj', 'eqearth'), ('lon_0', central_longitude)]
super().__init__(proj_params, central_longitude,
false_easting=false_easting,
false_northing=false_northing,
globe=globe)
self.threshold = 1e5
class Mollweide(_WarpedRectangularProjection):
"""
A Mollweide projection.
This projection is pseudocylindrical, and equal area. Parallels are
unequally-spaced straight lines, while meridians are elliptical arcs up to
semicircles on the edges. Poles are points.
It is commonly used for world maps, or interrupted with several central
meridians.
"""
_handles_ellipses = False
def __init__(self, central_longitude=0, globe=None,
false_easting=None, false_northing=None):
"""
Parameters
----------
central_longitude: float, optional
The central longitude. Defaults to 0.
false_easting: float, optional
X offset from planar origin in metres. Defaults to 0.
false_northing: float, optional
Y offset from planar origin in metres. Defaults to 0.
globe: :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
.. note::
This projection does not handle elliptical globes.
"""
proj4_params = [('proj', 'moll'), ('lon_0', central_longitude)]
super().__init__(proj4_params, central_longitude,
false_easting=false_easting,
false_northing=false_northing,
globe=globe)
self.threshold = 1e5
class Robinson(_WarpedRectangularProjection):
"""
A Robinson projection.
This projection is pseudocylindrical, and a compromise that is neither
equal-area nor conformal. Parallels are unequally-spaced straight lines,
and meridians are curved lines of no particular form.
It is commonly used for "visually-appealing" world maps.
"""
_handles_ellipses = False
def __init__(self, central_longitude=0, globe=None,
false_easting=None, false_northing=None):
"""
Parameters
----------
central_longitude: float, optional
The central longitude. Defaults to 0.
false_easting: float, optional
X offset from planar origin in metres. Defaults to 0.
false_northing: float, optional
Y offset from planar origin in metres. Defaults to 0.
globe: :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
.. note::
This projection does not handle elliptical globes.
"""
proj4_params = [('proj', 'robin'), ('lon_0', central_longitude)]
super().__init__(proj4_params, central_longitude,
false_easting=false_easting,
false_northing=false_northing,
globe=globe)
self.threshold = 1e4
def transform_point(self, x, y, src_crs, trap=True):
"""
Capture and handle any input NaNs, else invoke parent function,
:meth:`_WarpedRectangularProjection.transform_point`.
Needed because input NaNs can trigger a fatal error in the underlying
implementation of the Robinson projection.
Note
----
Although the original can in fact translate (nan, lat) into
(nan, y-value), this patched version doesn't support that.
"""
if np.isnan(x) or np.isnan(y):
result = (np.nan, np.nan)
else:
result = super().transform_point(x, y, src_crs, trap=trap)
return result
def transform_points(self, src_crs, x, y, z=None, trap=False):
"""
Capture and handle NaNs in input points -- else as parent function,
:meth:`_WarpedRectangularProjection.transform_points`.
Needed because input NaNs can trigger a fatal error in the underlying
implementation of the Robinson projection.
Note
----
Although the original can in fact translate (nan, lat) into
(nan, y-value), this patched version doesn't support that.
Instead, we invalidate any of the points that contain a NaN.
"""
input_point_nans = np.isnan(x) | np.isnan(y)
if z is not None:
input_point_nans |= np.isnan(z)
handle_nans = np.any(input_point_nans)
if handle_nans:
# Remove NaN points from input data to avoid the error.
x[input_point_nans] = 0.0
y[input_point_nans] = 0.0
if z is not None:
z[input_point_nans] = 0.0
result = super().transform_points(src_crs, x, y, z, trap=trap)
if handle_nans:
# Result always has shape (N, 3).
# Blank out each (whole) point where we had a NaN in the input.
result[input_point_nans] = np.nan
return result
class InterruptedGoodeHomolosine(Projection):
"""
Composite equal-area projection empahsizing either land or
ocean features.
Original Reference:
Goode, J. P., 1925: The Homolosine Projection: A new device for
portraying the Earth's surface entire. Annals of the
Association of American Geographers, 15:3, 119-125,
DOI: 10.1080/00045602509356949
A central_longitude value of -160 is recommended for the oceanic view.
"""
_wrappable = True
def __init__(self, central_longitude=0, globe=None, emphasis='land'):
"""
Parameters
----------
central_longitude : float, optional
The central longitude, by default 0
globe : :class:`cartopy.crs.Globe`, optional
If omitted, a default Globe object is created, by default None
emphasis : str, optional
Options 'land' and 'ocean' are available, by default 'land'
"""
if emphasis == 'land':
proj4_params = [('proj', 'igh'), ('lon_0', central_longitude)]
super().__init__(proj4_params, globe=globe)
elif emphasis == 'ocean':
proj4_params = [('proj', 'igh_o'), ('lon_0', central_longitude)]
super().__init__(proj4_params, globe=globe)
else:
msg = '`emphasis` needs to be either \'land\' or \'ocean\''
raise ValueError(msg)
minlon, maxlon = self._determine_longitude_bounds(central_longitude)
epsilon = 1e-10
# Obtain boundary points
n = 31
if emphasis == 'land':
top_interrupted_lons = (-40.0,)
bottom_interrupted_lons = (80.0, -20.0, -100.0)
elif emphasis == 'ocean':
top_interrupted_lons = (-90.0, 60.0)
bottom_interrupted_lons = (90.0, -60.0)
lons = np.empty(
(2 + 2 * len(top_interrupted_lons + bottom_interrupted_lons)) * n +
1)
lats = np.empty(
(2 + 2 * len(top_interrupted_lons + bottom_interrupted_lons)) * n +
1)
end = 0
# Left boundary
lons[end:end + n] = minlon
lats[end:end + n] = np.linspace(-90, 90, n)
end += n
# Top boundary
for lon in top_interrupted_lons:
lons[end:end + n] = lon - epsilon + central_longitude
lats[end:end + n] = np.linspace(90, 0, n)
end += n
lons[end:end + n] = lon + epsilon + central_longitude
lats[end:end + n] = np.linspace(0, 90, n)
end += n
# Right boundary
lons[end:end + n] = maxlon
lats[end:end + n] = np.linspace(90, -90, n)
end += n
# Bottom boundary
for lon in bottom_interrupted_lons:
lons[end:end + n] = lon + epsilon + central_longitude
lats[end:end + n] = np.linspace(-90, 0, n)
end += n
lons[end:end + n] = lon - epsilon + central_longitude
lats[end:end + n] = np.linspace(0, -90, n)
end += n
# Close loop
lons[-1] = minlon
lats[-1] = -90
points = self.transform_points(self.as_geodetic(), lons, lats)
self._boundary = sgeom.LinearRing(points)
mins = np.min(points, axis=0)
maxs = np.max(points, axis=0)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = 2e4
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class _Satellite(Projection):
def __init__(self, projection, satellite_height=35785831,
central_longitude=0.0, central_latitude=0.0,
false_easting=0, false_northing=0, globe=None,
sweep_axis=None):
proj4_params = [('proj', projection), ('lon_0', central_longitude),
('lat_0', central_latitude), ('h', satellite_height),
('x_0', false_easting), ('y_0', false_northing),
('units', 'm')]
if sweep_axis:
proj4_params.append(('sweep', sweep_axis))
super().__init__(proj4_params, globe=globe)
def _set_boundary(self, coords):
self._boundary = sgeom.LinearRing(coords.T)
mins = np.min(coords, axis=1)
maxs = np.max(coords, axis=1)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = np.diff(self._x_limits)[0] * 0.02
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class Geostationary(_Satellite):
"""
A view appropriate for satellites in Geostationary Earth orbit.
Perspective view looking directly down from above a point on the equator.
In this projection, the projected coordinates are scanning angles measured
from the satellite looking directly downward, multiplied by the height of
the satellite.
"""
def __init__(self, central_longitude=0.0, satellite_height=35785831,
false_easting=0, false_northing=0, globe=None,
sweep_axis='y'):
"""
Parameters
----------
central_longitude: float, optional
The central longitude. Defaults to 0.
satellite_height: float, optional
The height of the satellite. Defaults to 35785831 metres
(true geostationary orbit).
false_easting:
X offset from planar origin in metres. Defaults to 0.
false_northing:
Y offset from planar origin in metres. Defaults to 0.
globe: :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
sweep_axis: 'x' or 'y', optional. Defaults to 'y'.
Controls which axis is scanned first, and thus which angle is
applied first. The default is appropriate for Meteosat, while
'x' should be used for GOES.
"""
super().__init__(
projection='geos',
satellite_height=satellite_height,
central_longitude=central_longitude,
central_latitude=0.0,
false_easting=false_easting,
false_northing=false_northing,
globe=globe,
sweep_axis=sweep_axis)
# TODO: Let the globe return the semimajor axis always.
a = float(self.ellipsoid.semi_major_metre or WGS84_SEMIMAJOR_AXIS)
b = float(self.ellipsoid.semi_minor_metre or WGS84_SEMIMINOR_AXIS)
h = float(satellite_height)
# To find the bound we trace around where the line from the satellite
# is tangent to the surface. This involves trigonometry on a sphere
# centered at the satellite. The two scanning angles form two legs of
# triangle on this sphere--the hypotenuse "c" (angle arc) is controlled
# by distance from center to the edge of the ellipse being seen.
# This is one of the angles in the spherical triangle and used to
# rotate around and "scan" the boundary
angleA = np.linspace(0, -2 * np.pi, 91) # Clockwise boundary.
# Convert the angle around center to the proper value to use in the
# parametric form of an ellipse
th = np.arctan(a / b * np.tan(angleA))
# Given the position on the ellipse, what is the distance from center
# to the ellipse--and thus the tangent point
r = np.hypot(a * np.cos(th), b * np.sin(th))
sat_dist = a + h
# Using this distance, solve for sin and tan of c in the triangle that
# includes the satellite, Earth center, and tangent point--we need to
# figure out the location of this tangent point on the elliptical
# cross-section through the Earth towards the satellite, where the
# major axis is a and the minor is r. With the ellipse centered on the
# Earth and the satellite on the y-axis (at y = a + h = sat_dist), the
# equation for an ellipse and some calculus gives us the tangent point
# (x0, y0) as:
# y0 = a**2 / sat_dist
# x0 = r * np.sqrt(1 - a**2 / sat_dist**2)
# which gives:
# sin_c = x0 / np.hypot(x0, sat_dist - y0)
# tan_c = x0 / (sat_dist - y0)
# A bit of algebra combines these to give directly:
sin_c = r / np.sqrt(sat_dist ** 2 - a ** 2 + r ** 2)
tan_c = r / np.sqrt(sat_dist ** 2 - a ** 2)
# Using Napier's rules for right spherical triangles R2 and R6,
# (See https://en.wikipedia.org/wiki/Spherical_trigonometry), we can
# solve for arc angles b and a, which are our x and y scanning angles,
# respectively.
coords = np.vstack([np.arctan(np.cos(angleA) * tan_c), # R6
np.arcsin(np.sin(angleA) * sin_c)]) # R2
# Need to multiply scanning angles by satellite height to get to the
# actual native coordinates for the projection.
coords *= h
coords += np.array([[false_easting], [false_northing]])
self._set_boundary(coords)
class NearsidePerspective(_Satellite):
"""
Perspective view looking directly down from above a point on the globe.
In this projection, the projected coordinates are x and y measured from
the origin of a plane tangent to the Earth directly below the perspective
point (e.g. a satellite).
"""
_handles_ellipses = False
def __init__(self, central_longitude=0.0, central_latitude=0.0,
satellite_height=35785831,
false_easting=0, false_northing=0, globe=None):
"""
Parameters
----------
central_longitude: float, optional
The central longitude. Defaults to 0.
central_latitude: float, optional
The central latitude. Defaults to 0.
satellite_height: float, optional
The height of the satellite. Defaults to 35785831 meters
(true geostationary orbit).
false_easting:
X offset from planar origin in metres. Defaults to 0.
false_northing:
Y offset from planar origin in metres. Defaults to 0.
globe: :class:`cartopy.crs.Globe`, optional
If omitted, a default globe is created.
.. note::
This projection does not handle elliptical globes.
"""
super().__init__(
projection='nsper',
satellite_height=satellite_height,
central_longitude=central_longitude,
central_latitude=central_latitude,
false_easting=false_easting,
false_northing=false_northing,
globe=globe)
# TODO: Let the globe return the semimajor axis always.
a = self.ellipsoid.semi_major_metre or WGS84_SEMIMAJOR_AXIS
h = float(satellite_height)
max_x = a * np.sqrt(h / (2 * a + h))
coords = _ellipse_boundary(max_x, max_x,
false_easting, false_northing, 61)
self._set_boundary(coords)
class AlbersEqualArea(Projection):
"""
An Albers Equal Area projection
This projection is conic and equal-area, and is commonly used for maps of
the conterminous United States.
"""
def __init__(self, central_longitude=0.0, central_latitude=0.0,
false_easting=0.0, false_northing=0.0,
standard_parallels=(20.0, 50.0), globe=None):
"""
Parameters
----------
central_longitude: optional
The central longitude. Defaults to 0.
central_latitude: optional
The central latitude. Defaults to 0.
false_easting: optional
X offset from planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from planar origin in metres. Defaults to 0.
standard_parallels: optional
The one or two latitudes of correct scale. Defaults to (20, 50).
globe: optional
A :class:`cartopy.crs.Globe`. If omitted, a default globe is
created.
"""
proj4_params = [('proj', 'aea'),
('lon_0', central_longitude),
('lat_0', central_latitude),
('x_0', false_easting),
('y_0', false_northing)]
if standard_parallels is not None:
try:
proj4_params.append(('lat_1', standard_parallels[0]))
try:
proj4_params.append(('lat_2', standard_parallels[1]))
except IndexError:
pass
except TypeError:
proj4_params.append(('lat_1', standard_parallels))
super().__init__(proj4_params, globe=globe)
# bounds
minlon, maxlon = self._determine_longitude_bounds(central_longitude)
n = 103
lons = np.empty(2 * n + 1)
lats = np.empty(2 * n + 1)
tmp = np.linspace(minlon, maxlon, n)
lons[:n] = tmp
lats[:n] = 90
lons[n:-1] = tmp[::-1]
lats[n:-1] = -90
lons[-1] = lons[0]
lats[-1] = lats[0]
points = self.transform_points(self.as_geodetic(), lons, lats)
self._boundary = sgeom.LinearRing(points)
mins = np.min(points, axis=0)
maxs = np.max(points, axis=0)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = 1e5
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class AzimuthalEquidistant(Projection):
"""
An Azimuthal Equidistant projection
This projection provides accurate angles about and distances through the
central position. Other angles, distances, or areas may be distorted.
"""
_wrappable = True
def __init__(self, central_longitude=0.0, central_latitude=0.0,
false_easting=0.0, false_northing=0.0,
globe=None):
"""
Parameters
----------
central_longitude: optional
The true longitude of the central meridian in degrees.
Defaults to 0.
central_latitude: optional
The true latitude of the planar origin in degrees.
Defaults to 0.
false_easting: optional
X offset from the planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from the planar origin in metres. Defaults to 0.
globe: optional
An instance of :class:`cartopy.crs.Globe`. If omitted, a default
globe is created.
"""
proj4_params = [('proj', 'aeqd'), ('lon_0', central_longitude),
('lat_0', central_latitude),
('x_0', false_easting), ('y_0', false_northing)]
super().__init__(proj4_params, globe=globe)
# TODO: Let the globe return the semimajor axis always.
a = float(self.ellipsoid.semi_major_metre or WGS84_SEMIMAJOR_AXIS)
b = float(self.ellipsoid.semi_minor_metre or a)
coords = _ellipse_boundary(a * np.pi, b * np.pi,
false_easting, false_northing, 61)
self._boundary = sgeom.LinearRing(coords.T)
mins = np.min(coords, axis=1)
maxs = np.max(coords, axis=1)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = 1e5
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class Sinusoidal(Projection):
"""
A Sinusoidal projection.
This projection is equal-area.
"""
def __init__(self, central_longitude=0.0, false_easting=0.0,
false_northing=0.0, globe=None):
"""
Parameters
----------
central_longitude: optional
The central longitude. Defaults to 0.
false_easting: optional
X offset from planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from planar origin in metres. Defaults to 0.
globe: optional
A :class:`cartopy.crs.Globe`. If omitted, a default globe is
created.
"""
proj4_params = [('proj', 'sinu'),
('lon_0', central_longitude),
('x_0', false_easting),
('y_0', false_northing)]
super().__init__(proj4_params, globe=globe)
# Obtain boundary points
minlon, maxlon = self._determine_longitude_bounds(central_longitude)
points = []
n = 91
lon = np.empty(2 * n + 1)
lat = np.empty(2 * n + 1)
lon[:n] = minlon
lat[:n] = np.linspace(-90, 90, n)
lon[n:2 * n] = maxlon
lat[n:2 * n] = np.linspace(90, -90, n)
lon[-1] = minlon
lat[-1] = -90
points = self.transform_points(self.as_geodetic(), lon, lat)
self._boundary = sgeom.LinearRing(points)
mins = np.min(points, axis=0)
maxs = np.max(points, axis=0)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = max(np.abs(self.x_limits + self.y_limits)) * 1e-5
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
# MODIS data products use a Sinusoidal projection of a spherical Earth
# https://modis-land.gsfc.nasa.gov/GCTP.html
Sinusoidal.MODIS = Sinusoidal(globe=Globe(ellipse=None,
semimajor_axis=6371007.181,
semiminor_axis=6371007.181))
class EquidistantConic(Projection):
"""
An Equidistant Conic projection.
This projection is conic and equidistant, and the scale is true along all
meridians and along one or two specified standard parallels.
"""
def __init__(self, central_longitude=0.0, central_latitude=0.0,
false_easting=0.0, false_northing=0.0,
standard_parallels=(20.0, 50.0), globe=None):
"""
Parameters
----------
central_longitude: optional
The central longitude. Defaults to 0.
central_latitude: optional
The true latitude of the planar origin in degrees. Defaults to 0.
false_easting: optional
X offset from planar origin in metres. Defaults to 0.
false_northing: optional
Y offset from planar origin in metres. Defaults to 0.
standard_parallels: optional
The one or two latitudes of correct scale. Defaults to (20, 50).
globe: optional
A :class:`cartopy.crs.Globe`. If omitted, a default globe is
created.
"""
proj4_params = [('proj', 'eqdc'),
('lon_0', central_longitude),
('lat_0', central_latitude),
('x_0', false_easting),
('y_0', false_northing)]
if standard_parallels is not None:
try:
proj4_params.append(('lat_1', standard_parallels[0]))
try:
proj4_params.append(('lat_2', standard_parallels[1]))
except IndexError:
pass
except TypeError:
proj4_params.append(('lat_1', standard_parallels))
super().__init__(proj4_params, globe=globe)
# bounds
n = 103
lons = np.empty(2 * n + 1)
lats = np.empty(2 * n + 1)
minlon, maxlon = self._determine_longitude_bounds(central_longitude)
tmp = np.linspace(minlon, maxlon, n)
lons[:n] = tmp
lats[:n] = 90
lons[n:-1] = tmp[::-1]
lats[n:-1] = -90
lons[-1] = lons[0]
lats[-1] = lats[0]
points = self.transform_points(self.as_geodetic(), lons, lats)
self._boundary = sgeom.LinearRing(points)
mins = np.min(points, axis=0)
maxs = np.max(points, axis=0)
self._x_limits = mins[0], maxs[0]
self._y_limits = mins[1], maxs[1]
self.threshold = 1e5
@property
def boundary(self):
return self._boundary
@property
def x_limits(self):
return self._x_limits
@property
def y_limits(self):
return self._y_limits
class _BoundaryPoint:
def __init__(self, distance, kind, data):
"""
A representation for a geometric object which is
connected to the boundary.
Parameters
----------
distance: float
The distance along the boundary that this object
can be found.
kind: bool
Whether this object represents a point from the
pre-computed boundary.
data: point or namedtuple
The actual data that this boundary object represents.
"""
self.distance = distance
self.kind = kind
self.data = data
def __repr__(self):
return f'_BoundaryPoint({self.distance!r}, {self.kind!r}, {self.data})'
def _find_first_ge(a, x):
for v in a:
if v.distance >= x:
return v
# We've gone all the way around, so pick the first point again.
return a[0]
def epsg(code):
"""
Return the projection which corresponds to the given EPSG code.
The EPSG code must correspond to a "projected coordinate system",
so EPSG codes such as 4326 (WGS-84) which define a "geodetic coordinate
system" will not work.
Note
----
The conversion is performed by pyproj.CRS.
"""
import cartopy._epsg
return cartopy._epsg._EPSGProjection(code)
|