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
|
/******************************************************************************
* Project: PROJ.4
* Purpose: Implement a (currently minimalistic) proj API based primarily
* on the PJ_COORD 4D geodetic spatiotemporal data type.
*
* Author: Thomas Knudsen, thokn@sdfe.dk, 2016-06-09/2016-11-06
*
******************************************************************************
* Copyright (c) 2016, 2017 Thomas Knudsen/SDFE
*
* Permission is hereby granted, free of charge, to any person obtaining a
* copy of this software and associated documentation files (the "Software"),
* to deal in the Software without restriction, including without limitation
* the rights to use, copy, modify, merge, publish, distribute, sublicense,
* and/or sell copies of the Software, and to permit persons to whom the
* Software is furnished to do so, subject to the following conditions:
*
* The above copyright notice and this permission notice shall be included
* in all copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
* OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
* FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
* THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
* FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
* DEALINGS IN THE SOFTWARE.
*****************************************************************************/
#define FROM_PROJ_CPP
#include <assert.h>
#include <errno.h>
#include <stddef.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#ifndef _MSC_VER
#include <strings.h>
#endif
#include <algorithm>
#include <limits>
#include "proj.h"
#include "proj_experimental.h"
#include "proj_internal.h"
#include <math.h>
#include "geodesic.h"
#include "grids.hpp"
#include "filemanager.hpp"
#include "proj/common.hpp"
#include "proj/coordinateoperation.hpp"
#include "proj/internal/internal.hpp"
#include "proj/internal/io_internal.hpp"
using namespace NS_PROJ::internal;
/* Initialize PJ_COORD struct */
PJ_COORD proj_coord (double x, double y, double z, double t) {
PJ_COORD res;
res.v[0] = x;
res.v[1] = y;
res.v[2] = z;
res.v[3] = t;
return res;
}
static PJ_DIRECTION opposite_direction(PJ_DIRECTION dir) {
return static_cast<PJ_DIRECTION>(-dir);
}
/*****************************************************************************/
int proj_angular_input (PJ *P, enum PJ_DIRECTION dir) {
/******************************************************************************
Returns 1 if the operator P expects angular input coordinates when
operating in direction dir, 0 otherwise.
dir: {PJ_FWD, PJ_INV}
******************************************************************************/
if (PJ_FWD==dir)
return pj_left (P)==PJ_IO_UNITS_RADIANS;
return pj_right (P)==PJ_IO_UNITS_RADIANS;
}
/*****************************************************************************/
int proj_angular_output (PJ *P, enum PJ_DIRECTION dir) {
/******************************************************************************
Returns 1 if the operator P provides angular output coordinates when
operating in direction dir, 0 otherwise.
dir: {PJ_FWD, PJ_INV}
******************************************************************************/
return proj_angular_input (P, opposite_direction(dir));
}
/*****************************************************************************/
int proj_degree_input (PJ *P, enum PJ_DIRECTION dir) {
/******************************************************************************
Returns 1 if the operator P expects degree input coordinates when
operating in direction dir, 0 otherwise.
dir: {PJ_FWD, PJ_INV}
******************************************************************************/
if (PJ_FWD==dir)
return pj_left (P)==PJ_IO_UNITS_DEGREES;
return pj_right (P)==PJ_IO_UNITS_DEGREES;
}
/*****************************************************************************/
int proj_degree_output (PJ *P, enum PJ_DIRECTION dir) {
/******************************************************************************
Returns 1 if the operator P provides degree output coordinates when
operating in direction dir, 0 otherwise.
dir: {PJ_FWD, PJ_INV}
******************************************************************************/
return proj_degree_input (P, opposite_direction(dir));
}
/* Geodesic distance (in meter) + fwd and rev azimuth between two points on the ellipsoid */
PJ_COORD proj_geod (const PJ *P, PJ_COORD a, PJ_COORD b) {
PJ_COORD c;
if( !P->geod ) {
return proj_coord_error();
}
/* Note: the geodesic code takes arguments in degrees */
geod_inverse (P->geod,
PJ_TODEG(a.lpz.phi), PJ_TODEG(a.lpz.lam),
PJ_TODEG(b.lpz.phi), PJ_TODEG(b.lpz.lam),
c.v, c.v+1, c.v+2
);
// cppcheck-suppress uninitvar
return c;
}
/* Geodesic distance (in meter) between two points with angular 2D coordinates */
double proj_lp_dist (const PJ *P, PJ_COORD a, PJ_COORD b) {
double s12, azi1, azi2;
/* Note: the geodesic code takes arguments in degrees */
if( !P->geod ) {
return HUGE_VAL;
}
geod_inverse (P->geod,
PJ_TODEG(a.lpz.phi), PJ_TODEG(a.lpz.lam),
PJ_TODEG(b.lpz.phi), PJ_TODEG(b.lpz.lam),
&s12, &azi1, &azi2
);
return s12;
}
/* The geodesic distance AND the vertical offset */
double proj_lpz_dist (const PJ *P, PJ_COORD a, PJ_COORD b) {
if (HUGE_VAL==a.lpz.lam || HUGE_VAL==b.lpz.lam)
return HUGE_VAL;
return hypot (proj_lp_dist (P, a, b), a.lpz.z - b.lpz.z);
}
/* Euclidean distance between two points with linear 2D coordinates */
double proj_xy_dist (PJ_COORD a, PJ_COORD b) {
return hypot (a.xy.x - b.xy.x, a.xy.y - b.xy.y);
}
/* Euclidean distance between two points with linear 3D coordinates */
double proj_xyz_dist (PJ_COORD a, PJ_COORD b) {
return hypot (proj_xy_dist (a, b), a.xyz.z - b.xyz.z);
}
/* Measure numerical deviation after n roundtrips fwd-inv (or inv-fwd) */
double proj_roundtrip (PJ *P, PJ_DIRECTION direction, int n, PJ_COORD *coord) {
int i;
PJ_COORD t, org;
if (nullptr==P)
return HUGE_VAL;
if (n < 1) {
proj_log_error(P, _("n should be >= 1"));
proj_errno_set (P, PROJ_ERR_OTHER_API_MISUSE);
return HUGE_VAL;
}
/* in the first half-step, we generate the output value */
org = *coord;
*coord = proj_trans (P, direction, org);
t = *coord;
/* now we take n-1 full steps in inverse direction: We are */
/* out of phase due to the half step already taken */
for (i = 0; i < n - 1; i++)
t = proj_trans (P, direction, proj_trans (P, opposite_direction(direction), t) );
/* finally, we take the last half-step */
t = proj_trans (P, opposite_direction(direction), t);
/* checking for angular *input* since we do a roundtrip, and end where we begin */
if (proj_angular_input (P, direction))
return proj_lpz_dist (P, org, t);
return proj_xyz_dist (org, t);
}
/**************************************************************************************/
int pj_get_suggested_operation(PJ_CONTEXT*,
const std::vector<PJCoordOperation>& opList,
const int iExcluded[2],
PJ_DIRECTION direction,
PJ_COORD coord)
/**************************************************************************************/
{
// Select the operations that match the area of use
// and has the best accuracy.
int iBest = -1;
double bestAccuracy = std::numeric_limits<double>::max();
const int nOperations = static_cast<int>(opList.size());
for( int i = 0; i < nOperations; i++ ) {
if( i == iExcluded[0] || i == iExcluded[1] ) {
continue;
}
const auto &alt = opList[i];
bool spatialCriterionOK = false;
if( direction == PJ_FWD ) {
if( coord.xyzt.x >= alt.minxSrc &&
coord.xyzt.y >= alt.minySrc &&
coord.xyzt.x <= alt.maxxSrc &&
coord.xyzt.y <= alt.maxySrc) {
spatialCriterionOK = true;
}
} else {
if( coord.xyzt.x >= alt.minxDst &&
coord.xyzt.y >= alt.minyDst &&
coord.xyzt.x <= alt.maxxDst &&
coord.xyzt.y <= alt.maxyDst ) {
spatialCriterionOK = true;
}
}
if( spatialCriterionOK ) {
// The offshore test is for the "Test bug 245 (use +datum=carthage)"
// of testvarious. The long=10 lat=34 point belongs both to the
// onshore and offshore Tunisia area of uses, but is slightly
// onshore. So in a general way, prefer a onshore area to a
// offshore one.
if( iBest < 0 ||
(alt.accuracy >= 0 &&
(alt.accuracy < bestAccuracy ||
// If two operations have the same accuracy, use the one that
// is contained within a larger one
(alt.accuracy == bestAccuracy &&
alt.minxSrc > opList[iBest].minxSrc &&
alt.minySrc > opList[iBest].minySrc &&
alt.maxxSrc < opList[iBest].maxxSrc &&
alt.maxySrc < opList[iBest].maxySrc)) &&
!alt.isOffshore) ) {
iBest = i;
bestAccuracy = alt.accuracy;
}
}
}
return iBest;
}
/**************************************************************************************/
PJ_COORD proj_trans (PJ *P, PJ_DIRECTION direction, PJ_COORD coord) {
/***************************************************************************************
Apply the transformation P to the coordinate coord, preferring the 4D interfaces if
available.
See also pj_approx_2D_trans and pj_approx_3D_trans in pj_internal.c, which work
similarly, but prefers the 2D resp. 3D interfaces if available.
***************************************************************************************/
if (nullptr==P || direction == PJ_IDENT)
return coord;
if (P->inverted)
direction = opposite_direction(direction);
if( !P->alternativeCoordinateOperations.empty() ) {
constexpr int N_MAX_RETRY = 2;
int iExcluded[N_MAX_RETRY] = {-1, -1};
const int nOperations = static_cast<int>(
P->alternativeCoordinateOperations.size());
// We may need several attempts. For example the point at
// lon=-111.5 lat=45.26 falls into the bounding box of the Canadian
// ntv2_0.gsb grid, except that it is not in any of the subgrids, being
// in the US. We thus need another retry that will select the conus
// grid.
for( int iRetry = 0; iRetry <= N_MAX_RETRY; iRetry++ )
{
// Do a first pass and select the operations that match the area of use
// and has the best accuracy.
int iBest = pj_get_suggested_operation(P->ctx,
P->alternativeCoordinateOperations,
iExcluded,
direction,
coord);
if( iBest < 0 ) {
break;
}
if( iRetry > 0 ) {
const int oldErrno = proj_errno_reset(P);
if (proj_log_level(P->ctx, PJ_LOG_TELL) >= PJ_LOG_DEBUG) {
pj_log(P->ctx, PJ_LOG_DEBUG, proj_context_errno_string(P->ctx, oldErrno));
}
pj_log(P->ctx, PJ_LOG_DEBUG,
"Did not result in valid result. "
"Attempting a retry with another operation.");
}
const auto& alt = P->alternativeCoordinateOperations[iBest];
if( P->iCurCoordOp != iBest ) {
if (proj_log_level(P->ctx, PJ_LOG_TELL) >= PJ_LOG_DEBUG) {
std::string msg("Using coordinate operation ");
msg += alt.name;
pj_log(P->ctx, PJ_LOG_DEBUG, msg.c_str());
}
P->iCurCoordOp = iBest;
}
PJ_COORD res = direction == PJ_FWD ?
pj_fwd4d( coord, alt.pj ) : pj_inv4d( coord, alt.pj );
if( proj_errno(alt.pj) == PROJ_ERR_OTHER_NETWORK_ERROR ) {
return proj_coord_error ();
}
if( res.xyzt.x != HUGE_VAL ) {
return res;
}
if( iRetry == N_MAX_RETRY ) {
break;
}
iExcluded[iRetry] = iBest;
}
// In case we did not find an operation whose area of use is compatible
// with the input coordinate, then goes through again the list, and
// use the first operation that does not require grids.
NS_PROJ::io::DatabaseContextPtr dbContext;
try
{
if( P->ctx->cpp_context ) {
dbContext = P->ctx->cpp_context->getDatabaseContext().as_nullable();
}
}
catch( const std::exception& ) {}
for( int i = 0; i < nOperations; i++ ) {
const auto &alt = P->alternativeCoordinateOperations[i];
auto coordOperation = dynamic_cast<
NS_PROJ::operation::CoordinateOperation*>(alt.pj->iso_obj.get());
if( coordOperation ) {
if( coordOperation->gridsNeeded(dbContext, true).empty() ) {
if( P->iCurCoordOp != i ) {
if (proj_log_level(P->ctx, PJ_LOG_TELL) >= PJ_LOG_DEBUG) {
std::string msg("Using coordinate operation ");
msg += alt.name;
msg += " as a fallback due to lack of more "
"appropriate operations";
pj_log(P->ctx, PJ_LOG_DEBUG, msg.c_str());
}
P->iCurCoordOp = i;
}
if( direction == PJ_FWD ) {
return pj_fwd4d( coord, alt.pj );
}
else {
return pj_inv4d( coord, alt.pj );
}
}
}
}
proj_errno_set (P, PROJ_ERR_COORD_TRANSFM_NO_OPERATION);
return proj_coord_error ();
}
if (direction == PJ_FWD)
return pj_fwd4d (coord, P);
else
return pj_inv4d (coord, P);
}
/*****************************************************************************/
int proj_trans_array (PJ *P, PJ_DIRECTION direction, size_t n, PJ_COORD *coord) {
/******************************************************************************
Batch transform an array of PJ_COORD.
Performs transformation on all points, even if errors occur on some points.
Individual points that fail to transform will have their components set to
HUGE_VAL
Returns 0 if all coordinates are transformed without error, otherwise
returns a precise error number if all coordinates that fail to transform
for the same reason, or a generic error code if they fail for different
reasons.
******************************************************************************/
size_t i;
int retErrno = 0;
bool hasSetRetErrno = false;
bool sameRetErrno = true;
for (i = 0; i < n; i++) {
proj_context_errno_set(P->ctx, 0);
coord[i] = proj_trans (P, direction, coord[i]);
int thisErrno = proj_errno(P);
if( thisErrno != 0 )
{
if( !hasSetRetErrno )
{
retErrno = thisErrno;
hasSetRetErrno = true;
}
else if( sameRetErrno && retErrno != thisErrno )
{
sameRetErrno = false;
retErrno = PROJ_ERR_COORD_TRANSFM;
}
}
}
proj_context_errno_set(P->ctx, retErrno);
return retErrno;
}
/*************************************************************************************/
size_t proj_trans_generic (
PJ *P,
PJ_DIRECTION direction,
double *x, size_t sx, size_t nx,
double *y, size_t sy, size_t ny,
double *z, size_t sz, size_t nz,
double *t, size_t st, size_t nt
) {
/**************************************************************************************
Transform a series of coordinates, where the individual coordinate dimension
may be represented by an array that is either
1. fully populated
2. a null pointer and/or a length of zero, which will be treated as a
fully populated array of zeroes
3. of length one, i.e. a constant, which will be treated as a fully
populated array of that constant value
The strides, sx, sy, sz, st, represent the step length, in bytes, between
consecutive elements of the corresponding array. This makes it possible for
proj_transform to handle transformation of a large class of application
specific data structures, without necessarily understanding the data structure
format, as in:
typedef struct {double x, y; int quality_level; char surveyor_name[134];} XYQS;
XYQS survey[345];
double height = 23.45;
PJ *P = {...};
size_t stride = sizeof (XYQS);
...
proj_transform (
P, PJ_INV, sizeof(XYQS),
&(survey[0].x), stride, 345, (* We have 345 eastings *)
&(survey[0].y), stride, 345, (* ...and 345 northings. *)
&height, 1, (* The height is the constant 23.45 m *)
0, 0 (* and the time is the constant 0.00 s *)
);
This is similar to the inner workings of the pj_transform function, but the
stride functionality has been generalized to work for any size of basic unit,
not just a fixed number of doubles.
In most cases, the stride will be identical for x, y,z, and t, since they will
typically be either individual arrays (stride = sizeof(double)), or strided
views into an array of application specific data structures (stride = sizeof (...)).
But in order to support cases where x, y, z, and t come from heterogeneous
sources, individual strides, sx, sy, sz, st, are used.
Caveat: Since proj_transform does its work *in place*, this means that even the
supposedly constants (i.e. length 1 arrays) will return from the call in altered
state. Hence, remember to reinitialize between repeated calls.
Return value: Number of transformations completed.
**************************************************************************************/
PJ_COORD coord = {{0,0,0,0}};
size_t i, nmin;
double null_broadcast = 0;
double invalid_time = HUGE_VAL;
if (nullptr==P)
return 0;
if (P->inverted)
direction = opposite_direction(direction);
/* ignore lengths of null arrays */
if (nullptr==x) nx = 0;
if (nullptr==y) ny = 0;
if (nullptr==z) nz = 0;
if (nullptr==t) nt = 0;
/* and make the nullities point to some real world memory for broadcasting nulls */
if (0==nx) x = &null_broadcast;
if (0==ny) y = &null_broadcast;
if (0==nz) z = &null_broadcast;
if (0==nt) t = &invalid_time;
/* nothing to do? */
if (0==nx+ny+nz+nt)
return 0;
/* arrays of length 1 are constants, which we broadcast along the longer arrays */
/* so we need to find the length of the shortest non-unity array to figure out */
/* how many coordinate pairs we must transform */
nmin = (nx > 1)? nx: (ny > 1)? ny: (nz > 1)? nz: (nt > 1)? nt: 1;
if ((nx > 1) && (nx < nmin)) nmin = nx;
if ((ny > 1) && (ny < nmin)) nmin = ny;
if ((nz > 1) && (nz < nmin)) nmin = nz;
if ((nt > 1) && (nt < nmin)) nmin = nt;
/* Check validity of direction flag */
switch (direction) {
case PJ_FWD:
case PJ_INV:
break;
case PJ_IDENT:
return nmin;
}
/* Arrays of length==0 are broadcast as the constant 0 */
/* Arrays of length==1 are broadcast as their single value */
/* Arrays of length >1 are iterated over (for the first nmin values) */
/* The slightly convolved incremental indexing is used due */
/* to the stride, which may be any size supported by the platform */
for (i = 0; i < nmin; i++) {
coord.xyzt.x = *x;
coord.xyzt.y = *y;
coord.xyzt.z = *z;
coord.xyzt.t = *t;
coord = proj_trans(P, direction, coord);
/* in all full length cases, we overwrite the input with the output, */
/* and step on to the next element. */
/* The casts are somewhat funky, but they compile down to no-ops and */
/* they tell compilers and static analyzers that we know what we do */
if (nx > 1) {
*x = coord.xyzt.x;
x = (double *) ((void *) ( ((char *) x) + sx));
}
if (ny > 1) {
*y = coord.xyzt.y;
y = (double *) ((void *) ( ((char *) y) + sy));
}
if (nz > 1) {
*z = coord.xyzt.z;
z = (double *) ((void *) ( ((char *) z) + sz));
}
if (nt > 1) {
*t = coord.xyzt.t;
t = (double *) ((void *) ( ((char *) t) + st));
}
}
/* Last time around, we update the length 1 cases with their transformed alter egos */
if (nx==1)
*x = coord.xyzt.x;
if (ny==1)
*y = coord.xyzt.y;
if (nz==1)
*z = coord.xyzt.z;
if (nt==1)
*t = coord.xyzt.t;
return i;
}
/*************************************************************************************/
PJ_COORD pj_geocentric_latitude (const PJ *P, PJ_DIRECTION direction, PJ_COORD coord) {
/**************************************************************************************
Convert geographical latitude to geocentric (or the other way round if
direction = PJ_INV)
The conversion involves a call to the tangent function, which goes through the
roof at the poles, so very close (the last centimeter) to the poles no
conversion takes place and the input latitude is copied directly to the output.
Fortunately, the geocentric latitude converges to the geographical at the
poles, so the difference is negligible.
For the spherical case, the geographical latitude equals the geocentric, and
consequently, the input is copied directly to the output.
**************************************************************************************/
const double limit = M_HALFPI - 1e-9;
PJ_COORD res = coord;
if ((coord.lp.phi > limit) || (coord.lp.phi < -limit) || (P->es==0))
return res;
if (direction==PJ_FWD)
res.lp.phi = atan (P->one_es * tan (coord.lp.phi) );
else
res.lp.phi = atan (P->rone_es * tan (coord.lp.phi) );
return res;
}
double proj_torad (double angle_in_degrees) { return PJ_TORAD (angle_in_degrees);}
double proj_todeg (double angle_in_radians) { return PJ_TODEG (angle_in_radians);}
double proj_dmstor(const char *is, char **rs) {
return dmstor(is, rs);
}
char* proj_rtodms(char *s, double r, int pos, int neg) {
return rtodms(s, r, pos, neg);
}
/*************************************************************************************/
static PJ* skip_prep_fin(PJ *P) {
/**************************************************************************************
Skip prepare and finalize function for the various "helper operations" added to P when
in cs2cs compatibility mode.
**************************************************************************************/
P->skip_fwd_prepare = 1;
P->skip_fwd_finalize = 1;
P->skip_inv_prepare = 1;
P->skip_inv_finalize = 1;
return P;
}
/*************************************************************************************/
static int cs2cs_emulation_setup (PJ *P) {
/**************************************************************************************
If any cs2cs style modifiers are given (axis=..., towgs84=..., ) create the 4D API
equivalent operations, so the preparation and finalization steps in the pj_inv/pj_fwd
invocators can emulate the behavior of pj_transform and the cs2cs app.
Returns 1 on success, 0 on failure
**************************************************************************************/
PJ *Q;
paralist *p;
int do_cart = 0;
if (nullptr==P)
return 0;
/* Don't recurse when calling proj_create (which calls us back) */
if (pj_param_exists (P->params, "break_cs2cs_recursion"))
return 1;
/* Swap axes? */
p = pj_param_exists (P->params, "axis");
const bool disable_grid_presence_check = pj_param_exists (
P->params, "disable_grid_presence_check") != nullptr;
/* Don't axisswap if data are already in "enu" order */
if (p && (0!=strcmp ("enu", p->param))) {
char *def = static_cast<char*>(malloc (100+strlen(P->axis)));
if (nullptr==def)
return 0;
sprintf (def, "break_cs2cs_recursion proj=axisswap axis=%s", P->axis);
Q = pj_create_internal (P->ctx, def);
free (def);
if (nullptr==Q)
return 0;
P->axisswap = skip_prep_fin(Q);
}
/* Geoid grid(s) given? */
p = pj_param_exists (P->params, "geoidgrids");
if (!disable_grid_presence_check && p && strlen (p->param) > strlen ("geoidgrids=")) {
char *gridnames = p->param + strlen ("geoidgrids=");
char *def = static_cast<char*>(malloc (100+2*strlen(gridnames)));
if (nullptr==def)
return 0;
sprintf (def, "break_cs2cs_recursion proj=vgridshift grids=%s",
pj_double_quote_string_param_if_needed(gridnames).c_str());
Q = pj_create_internal (P->ctx, def);
free (def);
if (nullptr==Q)
return 0;
P->vgridshift = skip_prep_fin(Q);
}
/* Datum shift grid(s) given? */
p = pj_param_exists (P->params, "nadgrids");
if (!disable_grid_presence_check && p && strlen (p->param) > strlen ("nadgrids=")) {
char *gridnames = p->param + strlen ("nadgrids=");
char *def = static_cast<char*>(malloc (100+2*strlen(gridnames)));
if (nullptr==def)
return 0;
sprintf (def, "break_cs2cs_recursion proj=hgridshift grids=%s",
pj_double_quote_string_param_if_needed(gridnames).c_str());
Q = pj_create_internal (P->ctx, def);
free (def);
if (nullptr==Q)
return 0;
P->hgridshift = skip_prep_fin(Q);
}
/* We ignore helmert if we have grid shift */
p = P->hgridshift ? nullptr : pj_param_exists (P->params, "towgs84");
while (p) {
char *def;
char *s = p->param;
double *d = P->datum_params;
size_t n = strlen (s);
/* We ignore null helmert shifts (common in auto-translated resource files, e.g. epsg) */
if (0==d[0] && 0==d[1] && 0==d[2] && 0==d[3] && 0==d[4] && 0==d[5] && 0==d[6]) {
/* If the current ellipsoid is not WGS84, then make sure the */
/* change in ellipsoid is still done. */
if (!(fabs(P->a_orig - 6378137.0) < 1e-8 && fabs(P->es_orig - 0.0066943799901413) < 1e-15)) {
do_cart = 1;
}
break;
}
if (n <= 8) /* 8==strlen ("towgs84=") */
return 0;
def = static_cast<char*>(malloc (100+n));
if (nullptr==def)
return 0;
sprintf (def, "break_cs2cs_recursion proj=helmert exact %s convention=position_vector", s);
Q = pj_create_internal (P->ctx, def);
free(def);
if (nullptr==Q)
return 0;
pj_inherit_ellipsoid_def (P, Q);
P->helmert = skip_prep_fin (Q);
break;
}
/* We also need cartesian/geographical transformations if we are working in */
/* geocentric/cartesian space or we need to do a Helmert transform. */
if (P->is_geocent || P->helmert || do_cart) {
char def[150];
sprintf (def, "break_cs2cs_recursion proj=cart a=%40.20g es=%40.20g", P->a_orig, P->es_orig);
{
/* In case the current locale does not use dot but comma as decimal */
/* separator, replace it with dot, so that proj_atof() behaves */
/* correctly. */
/* TODO later: use C++ ostringstream with imbue(std::locale::classic()) */
/* to be locale unaware */
char* next_pos;
for (next_pos = def; (next_pos = strchr (next_pos, ',')) != nullptr; next_pos++) {
*next_pos = '.';
}
}
Q = pj_create_internal (P->ctx, def);
if (nullptr==Q)
return 0;
P->cart = skip_prep_fin (Q);
if (!P->is_geocent) {
sprintf (def, "break_cs2cs_recursion proj=cart ellps=WGS84");
Q = pj_create_internal (P->ctx, def);
if (nullptr==Q)
return 0;
P->cart_wgs84 = skip_prep_fin (Q);
}
}
return 1;
}
/*************************************************************************************/
PJ *pj_create_internal (PJ_CONTEXT *ctx, const char *definition) {
/*************************************************************************************/
/**************************************************************************************
Create a new PJ object in the context ctx, using the given definition. If ctx==0,
the default context is used, if definition==0, or invalid, a null-pointer is
returned. The definition may use '+' as argument start indicator, as in
"+proj=utm +zone=32", or leave it out, as in "proj=utm zone=32".
It may even use free formatting "proj = utm; zone =32 ellps= GRS80".
Note that the semicolon separator is allowed, but not required.
**************************************************************************************/
char *args, **argv;
size_t argc, n;
if (nullptr==ctx)
ctx = pj_get_default_ctx ();
/* Make a copy that we can manipulate */
n = strlen (definition);
args = (char *) malloc (n + 1);
if (nullptr==args) {
proj_context_errno_set(ctx, PROJ_ERR_OTHER /*ENOMEM*/);
return nullptr;
}
strcpy (args, definition);
argc = pj_trim_argc (args);
if (argc==0) {
free (args);
proj_context_errno_set(ctx, PROJ_ERR_INVALID_OP_MISSING_ARG);
return nullptr;
}
argv = pj_trim_argv (argc, args);
if (!argv) {
free(args);
proj_context_errno_set(ctx, PROJ_ERR_OTHER /*ENOMEM*/);
return nullptr;
}
PJ* P = pj_create_argv_internal (ctx, (int) argc, argv);
free (argv);
free (args);
return P;
}
/*************************************************************************************/
PJ *proj_create_argv (PJ_CONTEXT *ctx, int argc, char **argv) {
/**************************************************************************************
Create a new PJ object in the context ctx, using the given definition argument
array argv. If ctx==0, the default context is used, if definition==0, or invalid,
a null-pointer is returned. The definition arguments may use '+' as argument start
indicator, as in {"+proj=utm", "+zone=32"}, or leave it out, as in {"proj=utm",
"zone=32"}.
**************************************************************************************/
PJ *P;
const char *c;
if (nullptr==ctx)
ctx = pj_get_default_ctx ();
if (nullptr==argv) {
proj_context_errno_set(ctx, PROJ_ERR_INVALID_OP_MISSING_ARG);
return nullptr;
}
/* We assume that free format is used, and build a full proj_create compatible string */
c = pj_make_args (argc, argv);
if (nullptr==c) {
proj_context_errno_set(ctx, PROJ_ERR_INVALID_OP /* ENOMEM */);
return nullptr;
}
P = proj_create (ctx, c);
free ((char *) c);
return P;
}
/*************************************************************************************/
PJ *pj_create_argv_internal (PJ_CONTEXT *ctx, int argc, char **argv) {
/**************************************************************************************
For use by pipeline init function.
**************************************************************************************/
if (nullptr==ctx)
ctx = pj_get_default_ctx ();
if (nullptr==argv) {
proj_context_errno_set(ctx, PROJ_ERR_INVALID_OP_MISSING_ARG);
return nullptr;
}
/* ...and let pj_init_ctx do the hard work */
/* New interface: forbid init=epsg:XXXX syntax by default */
const int allow_init_epsg = proj_context_get_use_proj4_init_rules(ctx, FALSE);
PJ* P = pj_init_ctx_with_allow_init_epsg (ctx, argc, argv, allow_init_epsg);
/* Support cs2cs-style modifiers */
int ret = cs2cs_emulation_setup (P);
if (0==ret)
return proj_destroy (P);
return P;
}
/** Create an area of use */
PJ_AREA * proj_area_create(void) {
return static_cast<PJ_AREA*>(calloc(1, sizeof(PJ_AREA)));
}
/** Assign a bounding box to an area of use. */
void proj_area_set_bbox(PJ_AREA *area,
double west_lon_degree,
double south_lat_degree,
double east_lon_degree,
double north_lat_degree) {
area->bbox_set = TRUE;
area->west_lon_degree = west_lon_degree;
area->south_lat_degree = south_lat_degree;
area->east_lon_degree = east_lon_degree;
area->north_lat_degree = north_lat_degree;
}
/** Free an area of use */
void proj_area_destroy(PJ_AREA* area) {
free(area);
}
/************************************************************************/
/* proj_context_use_proj4_init_rules() */
/************************************************************************/
void proj_context_use_proj4_init_rules(PJ_CONTEXT *ctx, int enable) {
if( ctx == nullptr ) {
ctx = pj_get_default_ctx();
}
ctx->use_proj4_init_rules = enable;
}
/************************************************************************/
/* EQUAL() */
/************************************************************************/
static int EQUAL(const char* a, const char* b) {
#ifdef _MSC_VER
return _stricmp(a, b) == 0;
#else
return strcasecmp(a, b) == 0;
#endif
}
/************************************************************************/
/* proj_context_get_use_proj4_init_rules() */
/************************************************************************/
int proj_context_get_use_proj4_init_rules(PJ_CONTEXT *ctx, int from_legacy_code_path) {
const char* val = getenv("PROJ_USE_PROJ4_INIT_RULES");
if( ctx == nullptr ) {
ctx = pj_get_default_ctx();
}
if( val ) {
if( EQUAL(val, "yes") || EQUAL(val, "on") || EQUAL(val, "true") ) {
return TRUE;
}
if( EQUAL(val, "no") || EQUAL(val, "off") || EQUAL(val, "false") ) {
return FALSE;
}
pj_log(ctx, PJ_LOG_ERROR, "Invalid value for PROJ_USE_PROJ4_INIT_RULES");
}
if( ctx->use_proj4_init_rules >= 0 ) {
return ctx->use_proj4_init_rules;
}
return from_legacy_code_path;
}
/** Adds a " +type=crs" suffix to a PROJ string (if it is a PROJ string) */
std::string pj_add_type_crs_if_needed(const std::string& str)
{
std::string ret(str);
if( (starts_with(str, "proj=") ||
starts_with(str, "+proj=") ||
starts_with(str, "+init=") ||
starts_with(str, "+title=")) &&
str.find("type=crs") == std::string::npos )
{
ret += " +type=crs";
}
return ret;
}
// ---------------------------------------------------------------------------
static double simple_min(const double* data, const int arr_len) {
double min_value = data[0];
for( int iii = 1; iii < arr_len; iii++ ) {
if (data[iii] < min_value)
min_value = data[iii];
}
return min_value;
}
// ---------------------------------------------------------------------------
static double simple_max(const double* data, const int arr_len) {
double max_value = data[0];
for( int iii = 1; iii < arr_len; iii++ ) {
if ((data[iii] > max_value || max_value == HUGE_VAL) && data[iii] != HUGE_VAL)
max_value = data[iii];
}
return max_value;
}
// ---------------------------------------------------------------------------
static int _find_previous_index(const int iii, const double* data, const int arr_len) {
// find index of nearest valid previous value if exists
int prev_iii = iii - 1;
if (prev_iii == -1) // handle wraparound
prev_iii = arr_len - 1;
while (data[prev_iii] == HUGE_VAL && prev_iii != iii) {
prev_iii --;
if (prev_iii == -1) // handle wraparound
prev_iii = arr_len - 1;
}
return prev_iii;
}
// ---------------------------------------------------------------------------
/******************************************************************************
Handles the case when longitude values cross the antimeridian
when calculating the minimum.
Note: The data array must be in a linear ring.
Note: This requires a densified ring with at least 2 additional
points per edge to correctly handle global extents.
If only 1 additional point:
| |
|RL--x0--|RL--
| |
-180 180|-180
If they are evenly spaced and it crosses the antimeridian:
x0 - L = 180
R - x0 = -180
For example:
Let R = -179.9, x0 = 0.1, L = -179.89
x0 - L = 0.1 - -179.9 = 180
R - x0 = -179.89 - 0.1 ~= -180
This is the same in the case when it didn't cross the antimeridian.
If you have 2 additional points:
| |
|RL--x0--x1--|RL--
| |
-180 180|-180
If they are evenly spaced and it crosses the antimeridian:
x0 - L = 120
x1 - x0 = 120
R - x1 = -240
For example:
Let R = -179.9, x0 = -59.9, x1 = 60.1 L = -179.89
x0 - L = 59.9 - -179.9 = 120
x1 - x0 = 60.1 - 59.9 = 120
R - x1 = -179.89 - 60.1 ~= -240
However, if they are evenly spaced and it didn't cross the antimeridian:
x0 - L = 120
x1 - x0 = 120
R - x1 = 120
From this, we have a delta that is guaranteed to be significantly
large enough to tell the difference reguarless of the direction
the antimeridian was crossed.
However, even though the spacing was even in the source projection, it isn't
guaranteed in the target geographic projection. So, instead of 240, 200 is used
as it significantly larger than 120 to be sure that the antimeridian was crossed
but smalller than 240 to account for possible irregularities in distances
when re-projecting. Also, 200 ensures latitudes are ignored for axis order handling.
******************************************************************************/
static double antimeridian_min(const double* data, const int arr_len) {
double positive_min = HUGE_VAL;
double min_value = HUGE_VAL;
int crossed_meridian_count = 0;
bool positive_meridian = false;
for( int iii = 0; iii < arr_len; iii++ ) {
if (data[iii] == HUGE_VAL)
continue;
int prev_iii = _find_previous_index(iii, data, arr_len);
// check if crossed meridian
double delta = data[prev_iii] - data[iii];
// 180 -> -180
if (delta >= 200 && delta != HUGE_VAL) {
if (crossed_meridian_count == 0)
positive_min = min_value;
crossed_meridian_count ++;
positive_meridian = false;
// -180 -> 180
} else if (delta <= -200 && delta != HUGE_VAL) {
if (crossed_meridian_count == 0)
positive_min = data[iii];
crossed_meridian_count ++;
positive_meridian = true;
}
// positive meridian side min
if (positive_meridian && data[iii] < positive_min)
positive_min = data[iii];
// track general min value
if (data[iii] < min_value)
min_value = data[iii];
}
if (crossed_meridian_count == 2)
return positive_min;
else if (crossed_meridian_count == 4)
// bounds extends beyond -180/180
return -180;
return min_value;
}
// ---------------------------------------------------------------------------
// Handles the case when longitude values cross the antimeridian
// when calculating the minimum.
// Note: The data array must be in a linear ring.
// Note: This requires a densified ring with at least 2 additional
// points per edge to correctly handle global extents.
// See antimeridian_min docstring for reasoning.
static double antimeridian_max(const double* data, const int arr_len) {
double negative_max = -HUGE_VAL;
double max_value = -HUGE_VAL;
bool negative_meridian = false;
int crossed_meridian_count = 0;
for( int iii = 0; iii < arr_len; iii++ ) {
if (data[iii] == HUGE_VAL)
continue;
int prev_iii = _find_previous_index(iii, data, arr_len);
// check if crossed meridian
double delta = data[prev_iii] - data[iii];
// 180 -> -180
if (delta >= 200 && delta != HUGE_VAL) {
if (crossed_meridian_count == 0)
negative_max = data[iii];
crossed_meridian_count ++;
negative_meridian = true;
// -180 -> 180
} else if (delta <= -200 && delta != HUGE_VAL){
if (crossed_meridian_count == 0)
negative_max = max_value;
negative_meridian = false;
crossed_meridian_count++;
}
// negative meridian side max
if (negative_meridian
&& (data[iii] > negative_max || negative_max == HUGE_VAL)
&& data[iii] != HUGE_VAL
)
negative_max = data[iii];
// track general max value
if ((data[iii] > max_value || max_value == HUGE_VAL) && data[iii] != HUGE_VAL)
max_value = data[iii];
}
if (crossed_meridian_count == 2)
return negative_max;
else if (crossed_meridian_count == 4)
// bounds extends beyond -180/180
return 180;
return max_value;
}
// ---------------------------------------------------------------------------
// Check if the original projected bounds contains
// the north pole.
// This assumes that the destination CRS is geographic.
static bool contains_north_pole(
PJ* projobj,
PJ_DIRECTION pj_direction,
const double xmin,
const double ymin,
const double xmax,
const double ymax,
bool lon_lat_order
) {
double pole_y = 90;
double pole_x = 0;
if (!lon_lat_order) {
pole_y = 0;
pole_x = 90;
}
proj_trans_generic(
projobj,
opposite_direction(pj_direction),
&pole_x, sizeof(double), 1,
&pole_y, sizeof(double), 1,
nullptr, sizeof(double), 0,
nullptr, sizeof(double), 0
);
if (xmin < pole_x && pole_x < xmax && ymax > pole_y && pole_y > ymin)
return true;
return false;
}
// ---------------------------------------------------------------------------
// Check if the original projected bounds contains
// the south pole.
// This assumes that the destination CRS is geographic.
static bool contains_south_pole(
PJ* projobj,
PJ_DIRECTION pj_direction,
const double xmin,
const double ymin,
const double xmax,
const double ymax,
bool lon_lat_order
) {
double pole_y = -90;
double pole_x = 0;
if (!lon_lat_order) {
pole_y = 0;
pole_x = -90;
}
proj_trans_generic(
projobj,
opposite_direction(pj_direction),
&pole_x, sizeof(double), 1,
&pole_y, sizeof(double), 1,
nullptr, sizeof(double), 0,
nullptr, sizeof(double), 0
);
if (xmin < pole_x && pole_x < xmax && ymax > pole_y && pole_y > ymin)
return true;
return false;
}
// ---------------------------------------------------------------------------
// Check if the target CRS of the transformation
// has the longitude latitude axis order.
// This assumes that the destination CRS is geographic.
static int target_crs_lon_lat_order(
PJ_CONTEXT* transformer_ctx,
PJ* transformer_pj,
PJ_DIRECTION pj_direction
) {
PJ* target_crs = nullptr;
if (pj_direction == PJ_FWD)
target_crs = proj_get_target_crs(transformer_ctx, transformer_pj);
else if (pj_direction == PJ_INV)
target_crs = proj_get_source_crs(transformer_ctx, transformer_pj);
if (target_crs == nullptr) {
proj_context_log_debug(transformer_ctx, "Unable to retrieve target CRS");
return -1;
}
PJ* coord_system_pj = proj_crs_get_coordinate_system(
transformer_ctx,
target_crs
);
proj_destroy(target_crs);
if (coord_system_pj == nullptr) {
proj_context_log_debug(transformer_ctx, "Unable to get target CRS coordinate system.");
return -1;
}
const char* abbrev = nullptr;
int success = proj_cs_get_axis_info(
transformer_ctx,
coord_system_pj,
0,
nullptr,
&abbrev,
nullptr,
nullptr,
nullptr,
nullptr,
nullptr
);
proj_destroy(coord_system_pj);
if (success != 1)
return -1;
return strcmp(abbrev, "lon") == 0 || strcmp(abbrev, "Lon") == 0;
}
// ---------------------------------------------------------------------------
/** \brief Transform boundary,
*
* Transform boundary densifying the edges to account for nonlinear
* transformations along these edges and extracting the outermost bounds.
*
* If the destination CRS is geographic, the first axis is longitude,
* and xmax < xmin then the bounds crossed the antimeridian.
* In this scenario there are two polygons, one on each side of the antimeridian.
* The first polygon should be constructed with (xmin, ymin, 180, ymax)
* and the second with (-180, ymin, xmax, ymax).
*
* If the destination CRS is geographic, the first axis is latitude,
* and ymax < ymin then the bounds crossed the antimeridian.
* In this scenario there are two polygons, one on each side of the antimeridian.
* The first polygon should be constructed with (ymin, xmin, ymax, 180)
* and the second with (ymin, -180, ymax, xmax).
*
* @param context The PJ_CONTEXT object.
* @param P The PJ object representing the transformation.
* @param direction The direction of the transformation.
* @param xmin Minimum bounding coordinate of the first axis in source CRS
* (target CRS if direction is inverse).
* @param ymin Minimum bounding coordinate of the second axis in source CRS.
* (target CRS if direction is inverse).
* @param xmax Maximum bounding coordinate of the first axis in source CRS.
* (target CRS if direction is inverse).
* @param ymax Maximum bounding coordinate of the second axis in source CRS.
* (target CRS if direction is inverse).
* @param out_xmin Minimum bounding coordinate of the first axis in target CRS
* (source CRS if direction is inverse).
* @param out_ymin Minimum bounding coordinate of the second axis in target CRS.
* (source CRS if direction is inverse).
* @param out_xmax Maximum bounding coordinate of the first axis in target CRS.
* (source CRS if direction is inverse).
* @param out_ymax Maximum bounding coordinate of the second axis in target CRS.
* (source CRS if direction is inverse).
* @param densify_pts Recommended to use 21. This is the number of points
* to use to densify the bounding polygon in the transformation.
* @return an integer. 1 if successful. 0 if failures encountered.
* @since 8.2
*/
int proj_trans_bounds(PJ_CONTEXT* context,
PJ *P,
PJ_DIRECTION direction,
const double xmin,
const double ymin,
const double xmax,
const double ymax,
double* out_xmin,
double* out_ymin,
double* out_xmax,
double* out_ymax,
const int densify_pts
) {
*out_xmin = HUGE_VAL;
*out_ymin = HUGE_VAL;
*out_xmax = HUGE_VAL;
*out_ymax = HUGE_VAL;
if (P == nullptr) {
proj_log_error(P, _("NULL P object not allowed."));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
if (densify_pts < 0 || densify_pts > 10000) {
proj_log_error(P, _("densify_pts must be between 0-10000."));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
PJ_PROJ_INFO pj_info = proj_pj_info(P);
if (pj_info.id == nullptr) {
proj_log_error(P, _("NULL transformation not allowed,"));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
if (strcmp(pj_info.id, "noop") == 0 || direction == PJ_IDENT) {
*out_xmin = xmin;
*out_xmax = xmax;
*out_ymin = ymin;
*out_ymax = ymax;
return true;
}
bool degree_output = proj_degree_output(P, direction) != 0;
bool degree_input = proj_degree_input(P, direction) != 0;
if (degree_output && densify_pts < 2) {
proj_log_error(P, _("densify_pts must be at least 2 if the output is geograpic."));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
int side_pts = densify_pts + 1; // add one because we are densifying
const int boundary_len = side_pts * 4;
std::vector<double> x_boundary_array;
std::vector<double> y_boundary_array;
try
{
x_boundary_array.resize(boundary_len);
y_boundary_array.resize(boundary_len);
}
catch( const std::exception & e ) // memory allocation failure
{
proj_log_error(P, e.what());
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
double delta_x = 0;
double delta_y = 0;
bool north_pole_in_bounds = false;
bool south_pole_in_bounds = false;
bool input_lon_lat_order = false;
bool output_lon_lat_order = false;
if (degree_input) {
int in_order_lon_lat = target_crs_lon_lat_order(
context, P, opposite_direction(direction)
);
if (in_order_lon_lat == -1)
return false;
input_lon_lat_order = in_order_lon_lat != 0;
}
if (degree_output) {
int out_order_lon_lat = target_crs_lon_lat_order(context, P, direction);
if (out_order_lon_lat == -1)
return false;
output_lon_lat_order = out_order_lon_lat != 0;
north_pole_in_bounds = contains_north_pole(
P,
direction,
xmin,
ymin,
xmax,
ymax,
output_lon_lat_order
);
south_pole_in_bounds = contains_south_pole(
P,
direction,
xmin,
ymin,
xmax,
ymax,
output_lon_lat_order
);
}
if (degree_input && xmax < xmin) {
if (!input_lon_lat_order) {
proj_log_error(P, _("latitude max < latitude min."));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
// handle antimeridian
delta_x = (xmax - xmin + 360.0) / side_pts;
} else {
delta_x = (xmax - xmin) / side_pts;
}
if (degree_input && ymax < ymin) {
if (input_lon_lat_order) {
proj_log_error(P, _("latitude max < latitude min."));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return false;
}
// handle antimeridian
delta_y = (ymax - ymin + 360.0) / side_pts;
} else {
delta_y = (ymax - ymin) / side_pts;
}
// build densified bounding box
// Note: must be a linear ring for antimeridian logic
for( int iii = 0; iii < side_pts; iii++ )
{
// xmin boundary
y_boundary_array[iii] = ymax - iii * delta_y;
x_boundary_array[iii] = xmin;
// ymin boundary
y_boundary_array[iii + side_pts] = ymin;
x_boundary_array[iii + side_pts] = xmin + iii * delta_x;
// xmax boundary
y_boundary_array[iii + side_pts * 2] = ymin + iii * delta_y;
x_boundary_array[iii + side_pts * 2] = xmax;
// ymax boundary
y_boundary_array[iii + side_pts * 3] = ymax;
x_boundary_array[iii + side_pts * 3] = xmax - iii * delta_x;
}
proj_trans_generic (
P,
direction,
&x_boundary_array[0], sizeof(double), boundary_len,
&y_boundary_array[0], sizeof(double), boundary_len,
nullptr, 0, 0,
nullptr, 0, 0
);
if (!degree_output) {
*out_xmin = simple_min(&x_boundary_array[0], boundary_len);
*out_xmax = simple_max(&x_boundary_array[0], boundary_len);
*out_ymin = simple_min(&y_boundary_array[0], boundary_len);
*out_ymax = simple_max(&y_boundary_array[0], boundary_len);
} else if (north_pole_in_bounds && output_lon_lat_order) {
*out_xmin = -180;
*out_ymin = simple_min(&y_boundary_array[0], boundary_len);
*out_xmax = 180;
*out_ymax = 90;
} else if (north_pole_in_bounds) {
*out_xmin = simple_min(&x_boundary_array[0], boundary_len);
*out_ymin = -180;
*out_xmax = 90;
*out_ymax = 180;
} else if (south_pole_in_bounds && output_lon_lat_order) {
*out_xmin = -180;
*out_ymin = -90;
*out_xmax = 180;
*out_ymax = simple_max(&y_boundary_array[0], boundary_len);
} else if (south_pole_in_bounds) {
*out_xmin = -90;
*out_ymin = -180;
*out_xmax = simple_max(&x_boundary_array[0], boundary_len);
*out_ymax = 180;
} else if (output_lon_lat_order) {
*out_xmin = antimeridian_min(&x_boundary_array[0], boundary_len);
*out_xmax = antimeridian_max(&x_boundary_array[0], boundary_len);
*out_ymin = simple_min(&y_boundary_array[0], boundary_len);
*out_ymax = simple_max(&y_boundary_array[0], boundary_len);
} else {
*out_xmin = simple_min(&x_boundary_array[0], boundary_len);
*out_xmax = simple_max(&x_boundary_array[0], boundary_len);
*out_ymin = antimeridian_min(&y_boundary_array[0], boundary_len);
*out_ymax = antimeridian_max(&y_boundary_array[0], boundary_len);
}
return true;
}
/*****************************************************************************/
static void reproject_bbox(PJ* pjGeogToCrs,
double west_lon, double south_lat,
double east_lon, double north_lat,
double& minx,
double& miny,
double& maxx,
double& maxy) {
/*****************************************************************************/
minx = -std::numeric_limits<double>::max();
miny = -std::numeric_limits<double>::max();
maxx = std::numeric_limits<double>::max();
maxy = std::numeric_limits<double>::max();
if( !(west_lon == -180.0 && east_lon == 180.0 &&
south_lat == -90.0 && north_lat == 90.0) )
{
minx = -minx;
miny = -miny;
maxx = -maxx;
maxy = -maxy;
constexpr int N_STEPS = 20;
constexpr int N_STEPS_P1 = N_STEPS+1;
constexpr int XY_SIZE = N_STEPS_P1 * 4;
std::vector<double> x(XY_SIZE);
std::vector<double> y(XY_SIZE);
const double step_lon = (east_lon - west_lon) / N_STEPS;
const double step_lat = (north_lat - south_lat) / N_STEPS;
for( int j = 0; j <= N_STEPS; j++ )
{
x[j] = west_lon + j * step_lon;
y[j] = south_lat;
x[N_STEPS_P1+j] = x[j];
y[N_STEPS_P1+j] = north_lat;
x[N_STEPS_P1*2+j] = west_lon;
y[N_STEPS_P1*2+j] = south_lat + j * step_lat;
x[N_STEPS_P1*3+j] = east_lon;
y[N_STEPS_P1*3+j] = y[N_STEPS_P1*2+j];
}
proj_trans_generic (
pjGeogToCrs, PJ_FWD,
&x[0], sizeof(double), XY_SIZE,
&y[0], sizeof(double), XY_SIZE,
nullptr, 0, 0,
nullptr, 0, 0);
for( int j = 0; j < XY_SIZE; j++ )
{
if( x[j] != HUGE_VAL && y[j] != HUGE_VAL )
{
minx = std::min(minx, x[j]);
miny = std::min(miny, y[j]);
maxx = std::max(maxx, x[j]);
maxy = std::max(maxy, y[j]);
}
}
}
}
/*****************************************************************************/
static PJ* add_coord_op_to_list(
int idxInOriginalList,
PJ* op,
double west_lon, double south_lat,
double east_lon, double north_lat,
PJ* pjGeogToSrc,
PJ* pjGeogToDst,
bool isOffshore,
std::vector<PJCoordOperation>& altCoordOps) {
/*****************************************************************************/
double minxSrc;
double minySrc;
double maxxSrc;
double maxySrc;
double minxDst;
double minyDst;
double maxxDst;
double maxyDst;
reproject_bbox(pjGeogToSrc, west_lon, south_lat, east_lon, north_lat,
minxSrc, minySrc, maxxSrc, maxySrc);
reproject_bbox(pjGeogToDst, west_lon, south_lat, east_lon, north_lat,
minxDst, minyDst, maxxDst, maxyDst);
if( minxSrc <= maxxSrc && minxDst <= maxxDst )
{
const char* c_name = proj_get_name(op);
std::string name(c_name ? c_name : "");
const double accuracy = proj_coordoperation_get_accuracy(op->ctx, op);
altCoordOps.emplace_back(idxInOriginalList,
minxSrc, minySrc, maxxSrc, maxySrc,
minxDst, minyDst, maxxDst, maxyDst,
op, name, accuracy, isOffshore);
op = nullptr;
}
return op;
}
/*****************************************************************************/
static PJ* create_operation_to_geog_crs(PJ_CONTEXT* ctx, const PJ* crs) {
/*****************************************************************************/
// Create a geographic 2D long-lat degrees CRS that is related to the
// CRS
auto geodetic_crs = proj_crs_get_geodetic_crs(ctx, crs);
if( !geodetic_crs ) {
proj_context_log_debug(ctx, "Cannot find geodetic CRS matching CRS");
return nullptr;
}
auto geodetic_crs_type = proj_get_type(geodetic_crs);
if( geodetic_crs_type == PJ_TYPE_GEOCENTRIC_CRS ||
geodetic_crs_type == PJ_TYPE_GEOGRAPHIC_2D_CRS ||
geodetic_crs_type == PJ_TYPE_GEOGRAPHIC_3D_CRS )
{
auto datum = proj_crs_get_datum_forced(ctx, geodetic_crs);
assert( datum );
auto cs = proj_create_ellipsoidal_2D_cs(
ctx, PJ_ELLPS2D_LONGITUDE_LATITUDE, nullptr, 0);
auto ellps = proj_get_ellipsoid(ctx, datum);
proj_destroy(datum);
double semi_major_metre = 0;
double inv_flattening = 0;
proj_ellipsoid_get_parameters(ctx, ellps, &semi_major_metre,
nullptr, nullptr, &inv_flattening);
// It is critical to set the prime meridian to 0
auto temp = proj_create_geographic_crs(
ctx, "unnamed crs", "unnamed datum",
proj_get_name(ellps),
semi_major_metre, inv_flattening,
"Reference prime meridian", 0, nullptr, 0,
cs);
proj_destroy(ellps);
proj_destroy(cs);
proj_destroy(geodetic_crs);
geodetic_crs = temp;
geodetic_crs_type = proj_get_type(geodetic_crs);
}
if( geodetic_crs_type != PJ_TYPE_GEOGRAPHIC_2D_CRS )
{
// Shouldn't happen
proj_context_log_debug(ctx, "Cannot find geographic CRS matching CRS");
proj_destroy(geodetic_crs);
return nullptr;
}
// Create the transformation from this geographic 2D CRS to the source CRS
auto operation_ctx = proj_create_operation_factory_context(ctx, nullptr);
proj_operation_factory_context_set_spatial_criterion(
ctx, operation_ctx, PROJ_SPATIAL_CRITERION_PARTIAL_INTERSECTION);
proj_operation_factory_context_set_grid_availability_use(
ctx, operation_ctx, PROJ_GRID_AVAILABILITY_DISCARD_OPERATION_IF_MISSING_GRID);
auto target_crs_2D = proj_crs_demote_to_2D(ctx, nullptr, crs);
auto op_list_to_geodetic = proj_create_operations(
ctx, geodetic_crs, target_crs_2D, operation_ctx);
proj_destroy(target_crs_2D);
proj_operation_factory_context_destroy(operation_ctx);
proj_destroy(geodetic_crs);
const int nOpCount = op_list_to_geodetic == nullptr ? 0 :
proj_list_get_count(op_list_to_geodetic);
if( nOpCount == 0 )
{
proj_context_log_debug(ctx, "Cannot compute transformation from geographic CRS to CRS");
proj_list_destroy(op_list_to_geodetic);
return nullptr;
}
PJ* opGeogToCrs = nullptr;
// Use in priority operations *without* grids
for(int i = 0; i < nOpCount; i++ )
{
auto op = proj_list_get(ctx, op_list_to_geodetic, i);
assert(op);
if( proj_coordoperation_get_grid_used_count(ctx, op) == 0 )
{
opGeogToCrs = op;
break;
}
proj_destroy(op);
}
if( opGeogToCrs == nullptr )
{
opGeogToCrs = proj_list_get(ctx, op_list_to_geodetic, 0);
assert(opGeogToCrs);
}
proj_list_destroy(op_list_to_geodetic);
return opGeogToCrs;
}
/*****************************************************************************/
PJ *proj_create_crs_to_crs (PJ_CONTEXT *ctx, const char *source_crs, const char *target_crs, PJ_AREA *area) {
/******************************************************************************
Create a transformation pipeline between two known coordinate reference
systems.
See docs/source/development/reference/functions.rst
******************************************************************************/
if( !ctx ) {
ctx = pj_get_default_ctx();
}
PJ* src;
PJ* dst;
try
{
std::string source_crs_modified(pj_add_type_crs_if_needed(source_crs));
std::string target_crs_modified(pj_add_type_crs_if_needed(target_crs));
src = proj_create(ctx, source_crs_modified.c_str());
if( !src ) {
proj_context_log_debug(ctx, "Cannot instantiate source_crs");
return nullptr;
}
dst = proj_create(ctx, target_crs_modified.c_str());
if( !dst ) {
proj_context_log_debug(ctx, "Cannot instantiate target_crs");
proj_destroy(src);
return nullptr;
}
}
catch( const std::exception& )
{
return nullptr;
}
auto ret = proj_create_crs_to_crs_from_pj(ctx, src, dst, area, nullptr);
proj_destroy(src);
proj_destroy(dst);
return ret;
}
/*****************************************************************************/
std::vector<PJCoordOperation> pj_create_prepared_operations(PJ_CONTEXT *ctx,
const PJ *source_crs,
const PJ *target_crs,
PJ_OBJ_LIST* op_list)
/*****************************************************************************/
{
auto pjGeogToSrc = create_operation_to_geog_crs(ctx, source_crs);
if( !pjGeogToSrc )
{
proj_context_log_debug(ctx,
"Cannot create transformation from geographic CRS of source CRS to source CRS");
return {};
}
auto pjGeogToDst = create_operation_to_geog_crs(ctx, target_crs);
if( !pjGeogToDst )
{
proj_context_log_debug(ctx,
"Cannot create transformation from geographic CRS of target CRS to target CRS");
proj_destroy(pjGeogToSrc);
return {};
}
try
{
std::vector<PJCoordOperation> preparedOpList;
// Iterate over source->target candidate transformations and reproject
// their long-lat bounding box into the source CRS.
const auto op_count = proj_list_get_count(op_list);
for( int i = 0; i < op_count; i++ )
{
auto op = proj_list_get(ctx, op_list, i);
assert(op);
double west_lon = 0.0;
double south_lat = 0.0;
double east_lon = 0.0;
double north_lat = 0.0;
const char* areaName = nullptr;
if( proj_get_area_of_use(ctx, op,
&west_lon, &south_lat, &east_lon, &north_lat, &areaName) )
{
const bool isOffshore =
areaName && strstr(areaName, "- offshore");
if( west_lon <= east_lon )
{
op = add_coord_op_to_list(i, op,
west_lon, south_lat, east_lon, north_lat,
pjGeogToSrc, pjGeogToDst, isOffshore,
preparedOpList);
}
else
{
auto op_clone = proj_clone(ctx, op);
op = add_coord_op_to_list(i, op,
west_lon, south_lat, 180, north_lat,
pjGeogToSrc, pjGeogToDst, isOffshore,
preparedOpList);
op_clone = add_coord_op_to_list(i, op_clone,
-180, south_lat, east_lon, north_lat,
pjGeogToSrc, pjGeogToDst, isOffshore,
preparedOpList);
proj_destroy(op_clone);
}
}
proj_destroy(op);
}
proj_destroy(pjGeogToSrc);
proj_destroy(pjGeogToDst);
return preparedOpList;
}
catch( const std::exception& )
{
proj_destroy(pjGeogToSrc);
proj_destroy(pjGeogToDst);
return {};
}
}
// ---------------------------------------------------------------------------
//! @cond Doxygen_Suppress
static const char *getOptionValue(const char *option,
const char *keyWithEqual) noexcept {
if (ci_starts_with(option, keyWithEqual)) {
return option + strlen(keyWithEqual);
}
return nullptr;
}
//! @endcond
/*****************************************************************************/
PJ *proj_create_crs_to_crs_from_pj (PJ_CONTEXT *ctx, const PJ *source_crs, const PJ *target_crs, PJ_AREA *area, const char* const * options) {
/******************************************************************************
Create a transformation pipeline between two known coordinate reference
systems.
See docs/source/development/reference/functions.rst
******************************************************************************/
if( !ctx ) {
ctx = pj_get_default_ctx();
}
const char* authority = nullptr;
double accuracy = -1;
bool allowBallparkTransformations = true;
for (auto iter = options; iter && iter[0]; ++iter) {
const char *value;
if ((value = getOptionValue(*iter, "AUTHORITY="))) {
authority = value;
} else if ((value = getOptionValue(*iter, "ACCURACY="))) {
accuracy = pj_atof(value);
} else if ((value = getOptionValue(*iter, "ALLOW_BALLPARK="))) {
if( ci_equal(value, "yes") )
allowBallparkTransformations = true;
else if( ci_equal(value, "no") )
allowBallparkTransformations = false;
else {
ctx->logger(ctx->logger_app_data, PJ_LOG_ERROR,
"Invalid value for ALLOW_BALLPARK option.");
return nullptr;
}
} else {
std::string msg("Unknown option :");
msg += *iter;
ctx->logger(ctx->logger_app_data, PJ_LOG_ERROR, msg.c_str());
return nullptr;
}
}
auto operation_ctx = proj_create_operation_factory_context(ctx, authority);
if( !operation_ctx ) {
return nullptr;
}
proj_operation_factory_context_set_allow_ballpark_transformations(
ctx, operation_ctx, allowBallparkTransformations);
if( accuracy >= 0 ) {
proj_operation_factory_context_set_desired_accuracy(ctx, operation_ctx,
accuracy);
}
if( area && area->bbox_set ) {
proj_operation_factory_context_set_area_of_interest(
ctx,
operation_ctx,
area->west_lon_degree,
area->south_lat_degree,
area->east_lon_degree,
area->north_lat_degree);
}
proj_operation_factory_context_set_spatial_criterion(
ctx, operation_ctx, PROJ_SPATIAL_CRITERION_PARTIAL_INTERSECTION);
proj_operation_factory_context_set_grid_availability_use(
ctx, operation_ctx,
proj_context_is_network_enabled(ctx) ?
PROJ_GRID_AVAILABILITY_KNOWN_AVAILABLE:
PROJ_GRID_AVAILABILITY_DISCARD_OPERATION_IF_MISSING_GRID);
auto op_list = proj_create_operations(ctx, source_crs, target_crs, operation_ctx);
proj_operation_factory_context_destroy(operation_ctx);
if( !op_list ) {
return nullptr;
}
auto op_count = proj_list_get_count(op_list);
if( op_count == 0 ) {
proj_list_destroy(op_list);
proj_context_log_debug(ctx, "No operation found matching criteria");
return nullptr;
}
PJ* P = proj_list_get(ctx, op_list, 0);
assert(P);
if( P == nullptr || op_count == 1 || (area && area->bbox_set) ||
proj_get_type(source_crs) == PJ_TYPE_GEOCENTRIC_CRS ||
proj_get_type(target_crs) == PJ_TYPE_GEOCENTRIC_CRS ) {
proj_list_destroy(op_list);
return P;
}
auto preparedOpList = pj_create_prepared_operations(ctx, source_crs, target_crs,
op_list);
proj_list_destroy(op_list);
if( preparedOpList.empty() )
{
proj_destroy(P);
return nullptr;
}
// If there's finally juste a single result, return it directly
if( preparedOpList.size() == 1 )
{
auto retP = preparedOpList[0].pj;
preparedOpList[0].pj = nullptr;
proj_destroy(P);
return retP;
}
P->alternativeCoordinateOperations = std::move(preparedOpList);
// The returned P is rather dummy
P->descr = "Set of coordinate operations";
P->iso_obj = nullptr;
P->fwd = nullptr;
P->inv = nullptr;
P->fwd3d = nullptr;
P->inv3d = nullptr;
P->fwd4d = nullptr;
P->inv4d = nullptr;
return P;
}
/*****************************************************************************/
int proj_errno (const PJ *P) {
/******************************************************************************
Read an error level from the context of a PJ.
******************************************************************************/
return proj_context_errno (pj_get_ctx ((PJ *) P));
}
/*****************************************************************************/
int proj_context_errno (PJ_CONTEXT *ctx) {
/******************************************************************************
Read an error directly from a context, without going through a PJ
belonging to that context.
******************************************************************************/
if (nullptr==ctx)
ctx = pj_get_default_ctx();
return ctx->last_errno;
}
/*****************************************************************************/
int proj_errno_set (const PJ *P, int err) {
/******************************************************************************
Set context-errno, bubble it up to the thread local errno, return err
******************************************************************************/
/* Use proj_errno_reset to explicitly clear the error status */
if (0==err)
return 0;
/* For P==0 err goes to the default context */
proj_context_errno_set (pj_get_ctx ((PJ *) P), err);
errno = err;
return err;
}
/*****************************************************************************/
int proj_errno_restore (const PJ *P, int err) {
/******************************************************************************
Use proj_errno_restore when the current function succeeds, but the
error flag was set on entry, and stored/reset using proj_errno_reset
in order to monitor for new errors.
See usage example under proj_errno_reset ()
******************************************************************************/
if (0==err)
return 0;
proj_errno_set (P, err);
return 0;
}
/*****************************************************************************/
int proj_errno_reset (const PJ *P) {
/******************************************************************************
Clears errno in the context and thread local levels
through the low level pj_ctx interface.
Returns the previous value of the errno, for convenient reset/restore
operations:
int foo (PJ *P) {
// errno may be set on entry, but we need to reset it to be able to
// check for errors from "do_something_with_P(P)"
int last_errno = proj_errno_reset (P);
// local failure
if (0==P)
return proj_errno_set (P, 42);
// call to function that may fail
do_something_with_P (P);
// failure in do_something_with_P? - keep latest error status
if (proj_errno(P))
return proj_errno (P);
// success - restore previous error status, return 0
return proj_errno_restore (P, last_errno);
}
******************************************************************************/
int last_errno;
last_errno = proj_errno (P);
proj_context_errno_set (pj_get_ctx ((PJ *) P), 0);
errno = 0;
return last_errno;
}
/* Create a new context based on the default context */
PJ_CONTEXT *proj_context_create (void) {
return new (std::nothrow) pj_ctx(*pj_get_default_ctx());
}
PJ_CONTEXT *proj_context_destroy (PJ_CONTEXT *ctx) {
if (nullptr==ctx)
return nullptr;
/* Trying to free the default context is a no-op (since it is statically allocated) */
if (pj_get_default_ctx ()==ctx)
return nullptr;
delete ctx;
return nullptr;
}
/*****************************************************************************/
static char *path_append (char *buf, const char *app, size_t *buf_size) {
/******************************************************************************
Helper for proj_info() below. Append app to buf, separated by a
semicolon. Also handle allocation of longer buffer if needed.
Returns buffer and adjusts *buf_size through provided pointer arg.
******************************************************************************/
char *p;
size_t len, applen = 0, buflen = 0;
#ifdef _WIN32
const char *delim = ";";
#else
const char *delim = ":";
#endif
/* Nothing to do? */
if (nullptr == app)
return buf;
applen = strlen (app);
if (0 == applen)
return buf;
/* Start checking whether buf is long enough */
if (nullptr != buf)
buflen = strlen (buf);
len = buflen+applen+strlen (delim) + 1;
/* "pj_realloc", so to speak */
if (*buf_size < len) {
p = static_cast<char*>(calloc (2 * len, sizeof (char)));
if (nullptr==p) {
free (buf);
return nullptr;
}
*buf_size = 2 * len;
if (buf != nullptr)
strcpy (p, buf);
free (buf);
buf = p;
}
assert(buf);
/* Only append a delimiter if something's already there */
if (0 != buflen)
strcat (buf, delim);
strcat (buf, app);
return buf;
}
static const char *empty = {""};
static char version[64] = {""};
static PJ_INFO info = {0, 0, 0, nullptr, nullptr, nullptr, nullptr, 0};
/*****************************************************************************/
PJ_INFO proj_info (void) {
/******************************************************************************
Basic info about the current instance of the PROJ.4 library.
Returns PJ_INFO struct.
******************************************************************************/
size_t buf_size = 0;
char *buf = nullptr;
pj_acquire_lock ();
info.major = PROJ_VERSION_MAJOR;
info.minor = PROJ_VERSION_MINOR;
info.patch = PROJ_VERSION_PATCH;
/* This is a controlled environment, so no risk of sprintf buffer
overflow. A normal version string is xx.yy.zz which is 8 characters
long and there is room for 64 bytes in the version string. */
sprintf (version, "%d.%d.%d", info.major, info.minor, info.patch);
info.version = version;
info.release = pj_get_release ();
/* build search path string */
auto ctx = pj_get_default_ctx();
if (ctx->search_paths.empty()) {
const auto searchpaths = pj_get_default_searchpaths(ctx);
for( const auto& path: searchpaths ) {
buf = path_append(buf, path.c_str(), &buf_size);
}
} else {
for (const auto &path : ctx->search_paths) {
buf = path_append(buf, path.c_str(), &buf_size);
}
}
free(const_cast<char*>(info.searchpath));
info.searchpath = buf ? buf : empty;
info.paths = ctx->c_compat_paths;
info.path_count = static_cast<int>(ctx->search_paths.size());
pj_release_lock ();
return info;
}
/*****************************************************************************/
PJ_PROJ_INFO proj_pj_info(PJ *P) {
/******************************************************************************
Basic info about a particular instance of a projection object.
Returns PJ_PROJ_INFO struct.
******************************************************************************/
PJ_PROJ_INFO pjinfo;
char *def;
memset(&pjinfo, 0, sizeof(PJ_PROJ_INFO));
pjinfo.accuracy = -1.0;
if (nullptr==P)
return pjinfo;
/* coordinate operation description */
if( P->iCurCoordOp >= 0 ) {
P = P->alternativeCoordinateOperations[P->iCurCoordOp].pj;
} else if( !P->alternativeCoordinateOperations.empty() ) {
pjinfo.id = "unknown";
pjinfo.description = "unavailable until proj_trans is called";
pjinfo.definition = "unavailable until proj_trans is called";
return pjinfo;
}
/* projection id */
if (pj_param(P->ctx, P->params, "tproj").i)
pjinfo.id = pj_param(P->ctx, P->params, "sproj").s;
if( P->iso_obj ) {
pjinfo.description = P->iso_obj->nameStr().c_str();
} else {
pjinfo.description = P->descr;
}
// accuracy
if( P->iso_obj ) {
auto conv = dynamic_cast<const NS_PROJ::operation::Conversion*>(P->iso_obj.get());
if( conv ) {
pjinfo.accuracy = 0.0;
} else {
auto op = dynamic_cast<const NS_PROJ::operation::CoordinateOperation*>(P->iso_obj.get());
if( op ) {
const auto& accuracies = op->coordinateOperationAccuracies();
if( !accuracies.empty() ) {
try {
pjinfo.accuracy = std::stod(accuracies[0]->value());
} catch ( const std::exception& ) {}
}
}
}
}
/* projection definition */
if (P->def_full)
def = P->def_full;
else
def = pj_get_def(P, 0); /* pj_get_def takes a non-const PJ pointer */
if (nullptr==def)
pjinfo.definition = empty;
else
pjinfo.definition = pj_shrink (def);
/* Make proj_destroy clean this up eventually */
P->def_full = def;
pjinfo.has_inverse = pj_has_inverse(P);
return pjinfo;
}
/*****************************************************************************/
PJ_GRID_INFO proj_grid_info(const char *gridname) {
/******************************************************************************
Information about a named datum grid.
Returns PJ_GRID_INFO struct.
******************************************************************************/
PJ_GRID_INFO grinfo;
/*PJ_CONTEXT *ctx = proj_context_create(); */
PJ_CONTEXT *ctx = pj_get_default_ctx();
memset(&grinfo, 0, sizeof(PJ_GRID_INFO));
const auto fillGridInfo = [&grinfo, ctx, gridname]
(const NS_PROJ::Grid& grid, const std::string& format)
{
const auto& extent = grid.extentAndRes();
/* name of grid */
strncpy (grinfo.gridname, gridname, sizeof(grinfo.gridname) - 1);
/* full path of grid */
pj_find_file(ctx, gridname, grinfo.filename, sizeof(grinfo.filename) - 1);
/* grid format */
strncpy (grinfo.format, format.c_str(), sizeof(grinfo.format) - 1);
/* grid size */
grinfo.n_lon = grid.width();
grinfo.n_lat = grid.height();
/* cell size */
grinfo.cs_lon = extent.resX;
grinfo.cs_lat = extent.resY;
/* bounds of grid */
grinfo.lowerleft.lam = extent.west;
grinfo.lowerleft.phi = extent.south;
grinfo.upperright.lam = extent.east;
grinfo.upperright.phi = extent.north;
};
{
const auto gridSet = NS_PROJ::VerticalShiftGridSet::open(ctx, gridname);
if( gridSet )
{
const auto& grids = gridSet->grids();
if( !grids.empty() )
{
const auto& grid = grids.front();
fillGridInfo(*grid, gridSet->format());
return grinfo;
}
}
}
{
const auto gridSet = NS_PROJ::HorizontalShiftGridSet::open(ctx, gridname);
if( gridSet )
{
const auto& grids = gridSet->grids();
if( !grids.empty() )
{
const auto& grid = grids.front();
fillGridInfo(*grid, gridSet->format());
return grinfo;
}
}
}
strcpy(grinfo.format, "missing");
return grinfo;
}
/*****************************************************************************/
PJ_INIT_INFO proj_init_info(const char *initname){
/******************************************************************************
Information about a named init file.
Maximum length of initname is 64.
Returns PJ_INIT_INFO struct.
If the init file is not found all members of the return struct are set
to the empty string.
If the init file is found, but the metadata is missing, the value is
set to "Unknown".
******************************************************************************/
int file_found;
char param[80], key[74];
paralist *start, *next;
PJ_INIT_INFO ininfo;
PJ_CONTEXT *ctx = pj_get_default_ctx();
memset(&ininfo, 0, sizeof(PJ_INIT_INFO));
file_found = pj_find_file(ctx, initname, ininfo.filename, sizeof(ininfo.filename));
if (!file_found || strlen(initname) > 64) {
if( strcmp(initname, "epsg") == 0 || strcmp(initname, "EPSG") == 0 ) {
const char* val;
proj_context_errno_set( ctx, 0 );
strncpy (ininfo.name, initname, sizeof(ininfo.name) - 1);
strcpy(ininfo.origin, "EPSG");
val = proj_context_get_database_metadata(ctx, "EPSG.VERSION");
if( val ) {
strncpy(ininfo.version, val, sizeof(ininfo.version) - 1);
}
val = proj_context_get_database_metadata(ctx, "EPSG.DATE");
if( val ) {
strncpy(ininfo.lastupdate, val, sizeof(ininfo.lastupdate) - 1);
}
return ininfo;
}
if( strcmp(initname, "IGNF") == 0 ) {
const char* val;
proj_context_errno_set( ctx, 0 );
strncpy (ininfo.name, initname, sizeof(ininfo.name) - 1);
strcpy(ininfo.origin, "IGNF");
val = proj_context_get_database_metadata(ctx, "IGNF.VERSION");
if( val ) {
strncpy(ininfo.version, val, sizeof(ininfo.version) - 1);
}
val = proj_context_get_database_metadata(ctx, "IGNF.DATE");
if( val ) {
strncpy(ininfo.lastupdate, val, sizeof(ininfo.lastupdate) - 1);
}
return ininfo;
}
return ininfo;
}
/* The initial memset (0) makes strncpy safe here */
strncpy (ininfo.name, initname, sizeof(ininfo.name) - 1);
strcpy(ininfo.origin, "Unknown");
strcpy(ininfo.version, "Unknown");
strcpy(ininfo.lastupdate, "Unknown");
strncpy (key, initname, 64); /* make room for ":metadata\0" at the end */
key[64] = 0;
memcpy(key + strlen(key), ":metadata", 9 + 1);
strcpy(param, "+init=");
/* The +strlen(param) avoids a cppcheck false positive warning */
strncat(param + strlen(param), key, sizeof(param)-1-strlen(param));
start = pj_mkparam(param);
pj_expand_init(ctx, start);
if (pj_param(ctx, start, "tversion").i)
strncpy(ininfo.version, pj_param(ctx, start, "sversion").s, sizeof(ininfo.version) - 1);
if (pj_param(ctx, start, "torigin").i)
strncpy(ininfo.origin, pj_param(ctx, start, "sorigin").s, sizeof(ininfo.origin) - 1);
if (pj_param(ctx, start, "tlastupdate").i)
strncpy(ininfo.lastupdate, pj_param(ctx, start, "slastupdate").s, sizeof(ininfo.lastupdate) - 1);
for ( ; start; start = next) {
next = start->next;
free(start);
}
return ininfo;
}
/*****************************************************************************/
PJ_FACTORS proj_factors(PJ *P, PJ_COORD lp) {
/******************************************************************************
Cartographic characteristics at point lp.
Characteristics include meridian, parallel and areal scales, angular
distortion, meridian/parallel, meridian convergence and scale error.
returns PJ_FACTORS. If unsuccessful, error number is set and the
struct returned contains NULL data.
******************************************************************************/
PJ_FACTORS factors = {0,0,0, 0,0,0, 0,0, 0,0,0,0};
struct FACTORS f;
if (nullptr==P)
return factors;
const auto type = proj_get_type(P);
if( type == PJ_TYPE_PROJECTED_CRS )
{
// If it is a projected CRS, then compute the factors on the conversion
// associated to it. We need to start from a temporary geographic CRS
// using the same datum as the one of the projected CRS, and with
// input coordinates being in longitude, latitude order in radian,
// to be consistent with the expectations of the lp input parameter.
auto ctx = P->ctx;
auto geodetic_crs = proj_get_source_crs(ctx, P);
assert(geodetic_crs);
auto datum = proj_crs_get_datum(ctx, geodetic_crs);
auto datum_ensemble = proj_crs_get_datum_ensemble(ctx, geodetic_crs);
auto cs = proj_create_ellipsoidal_2D_cs(
ctx, PJ_ELLPS2D_LONGITUDE_LATITUDE, "Radian", 1.0);
auto temp = proj_create_geographic_crs_from_datum(
ctx, "unnamed crs", datum ? datum : datum_ensemble,
cs);
proj_destroy(datum);
proj_destroy(datum_ensemble);
proj_destroy(cs);
proj_destroy(geodetic_crs);
auto newOp = proj_create_crs_to_crs_from_pj(ctx, temp, P, nullptr, nullptr);
proj_destroy(temp);
assert(newOp);
auto ret = proj_factors(newOp, lp);
proj_destroy(newOp);
return ret;
}
if( type != PJ_TYPE_CONVERSION &&
type != PJ_TYPE_TRANSFORMATION &&
type != PJ_TYPE_CONCATENATED_OPERATION &&
type != PJ_TYPE_OTHER_COORDINATE_OPERATION )
{
proj_log_error(P, _("Invalid type for P object"));
proj_errno_set (P, PROJ_ERR_INVALID_OP_ILLEGAL_ARG_VALUE);
return factors;
}
if (pj_factors(lp.lp, P, 0.0, &f))
return factors;
factors.meridional_scale = f.h;
factors.parallel_scale = f.k;
factors.areal_scale = f.s;
factors.angular_distortion = f.omega;
factors.meridian_parallel_angle = f.thetap;
factors.meridian_convergence = f.conv;
factors.tissot_semimajor = f.a;
factors.tissot_semiminor = f.b;
/* Raw derivatives, for completeness's sake */
factors.dx_dlam = f.der.x_l;
factors.dx_dphi = f.der.x_p;
factors.dy_dlam = f.der.y_l;
factors.dy_dphi = f.der.y_p;
return factors;
}
|