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
|
/*
* Copyright (c) 1999-2009 Delft University of Technology, The Netherlands
*
* This file is part of Doris, the Delft o-o radar interferometric software.
*
* Doris program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* Doris is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*
*/
/****************************************************************
* $Source: /users/kampes/DEVELOP/DORIS/doris/src/RCS/referencephase.cc,v $ *
* $Revision: 3.22 $ *
* $Date: 2006/05/18 11:09:20 $ *
* $Author: kampes $ *
* *
* -computation flat earth correction. *
* -computation radarcoding dem + interpolation to (1,1) grid *
****************************************************************/
#include "matrixbk.hh"
#include "orbitbk.hh"
#include "slcimage.hh" // my slc image class
#include "productinfo.hh" // my 'products' class
#include "constants.hh"
#include "referencephase.hh" // proto types
#include "ioroutines.hh" // error etc.
#include "utilities.hh" // isodd; ones(), etc.
#include "coregistration.hh" // distribute points
#include "exceptions.hh" // my exceptions class
#include <iomanip> // setw only for test..
#include <cstdlib> // system
#include <algorithm> // max
#ifdef WIN32
// Jia defined this.
// Bert Kampes, 24-Aug-2005
#include "winsock2.h"
#else
#include <netinet/in.h> // ntohl byteorder x86-HP unix
#endif
// Using A Two-Dimensional Quality Mesh Generator and Delaunay Triangulator
// from Jonathan Richard Shewchuk
// Some definition for triangulate call
#define VOID int
#define REAL double
#define ANSI_DECLARATORS
extern "C" {
#include "triangle.h"
}
/****************************************************************
* flatearth *
* *
* Compute polynomial model for 'flat earth' correction. *
* fie(l,p) = sumj=0:d sumk=0:d Ajk l^j p^k (NOT bert 8sept99) *
* precise orbits are used to compute delta range for Npoints *
* after which the polynomial model is fitted (LS). *
* *
* input: *
* - inputoptions *
* - info structs *
* - platform data points *
* output: *
* - void (result to file "scratchresflat") *
* - coefficients normalized wrt. original window of master *
* *
* Bert Kampes, 09-Mar-1999 *
* Bert Kampes, 26-Oct-1999 normalization of coeff., *
* dump to logfile: var(unknowns) == diag(inv(AtA)) *
****************************************************************/
void flatearth(
const input_comprefpha &comprefphainput,
const input_ell &ellips,
const slcimage &master,
const slcimage &slave,
const productinfo &interferogram,
orbit &masterorbit,
orbit &slaveorbit)
{
TRACE_FUNCTION("flatearth (BK 26-Oct-1999)")
const int32 MAXITER = 10;
const real8 CRITERPOS = 1e-6;
const real8 CRITERTIM = 1e-10;
char dummyline[2*ONE27];
INFO << "FLATEARTH: MAXITER: " << MAXITER << "; "
<< "CRITERPOS: " << CRITERPOS << " m; "
<< "CRITERTIM: " << CRITERTIM << " s";
INFO.print();
// ______ Normalization factors for polynomial ______
const real8 minL = master.originalwindow.linelo;
const real8 maxL = master.originalwindow.linehi;
const real8 minP = master.originalwindow.pixlo;
const real8 maxP = master.originalwindow.pixhi;
INFO << "flatearth: polynomial normalized by factors: "
<< minL << " " << maxL << " " << minP << " " << maxP
<< " to [-2,2]";
INFO.print();
// ______Handling of input______
const real8 m_minpi4cdivlam = (-4.0*PI*SOL)/master.wavelength;
const real8 s_minpi4cdivlam = (-4.0*PI*SOL)/slave.wavelength;
DEBUG << "master wavelength = " << master.wavelength;
DEBUG.print();
DEBUG << "slave wavelength = " << slave.wavelength;
DEBUG.print();
const int32 DEGREE = comprefphainput.degree;
const int32 Nunk = Ncoeffs(DEGREE); // Number of unknowns
bool pointsrandom = true;
if (specified(comprefphainput.ifpositions))
pointsrandom = false; // only use those points
// ______ Distribute points wel distributed over win ______
// ______ or read from ascii file ______
// ______(i,0): line, (i,1): pixel, (i,2) flagfromdisk______
matrix<uint> Position;
const uint Npoints = comprefphainput.Npoints;
register int32 i,j,k,index;
if (pointsrandom) // no filename specified
{
Position = distributepoints(Npoints,interferogram.win);
}
else // read from file
{
Position.resize(Npoints,3);
//ifstream ifpos(comprefphainput.ifpositions, ios::in);
ifstream ifpos;
openfstream(ifpos,comprefphainput.ifpositions);
bk_assert(ifpos,comprefphainput.ifpositions,__FILE__,__LINE__);
uint ll,pp;
for (i=0; i<Npoints; ++i)
{
ifpos >> ll >> pp;
Position(i,0) = uint(ll);
Position(i,1) = uint(pp);
Position(i,2) = uint(1); // flag from file
ifpos.getline(dummyline,2*ONE27,'\n'); // goto next line.
}
ifpos.close();
// ______ Check last point ivm. EOL after last position in file ______
if (Position(Npoints-1,0) == Position(Npoints-2,0) &&
Position(Npoints-1,1) == Position(Npoints-2,1))
{
Position(Npoints-1,0) = uint(.5*(minL + maxL) + 27); // random
Position(Npoints-1,1) = uint(.5*(minP + maxP) + 37); // random
WARNING << "refpha: there should be no EOL after last point in file: "
<< comprefphainput.ifpositions;
WARNING.print();
}
// ______ Check if points are in overlap ______
// ______ no check for uniqueness of points ______
}
matrix<real8> y(Npoints,1); // observation
matrix<real8> y_h2ph(Npoints,1); // observation, h2ph factors, added by FvL
matrix<real8> A(Npoints,Nunk); // designmatrix
// ______Check redundancy______
if (Npoints < Nunk)
{
PRINT_ERROR("flatearth: Number of points is smaller than parameters solved for.")
throw(input_error);
}
// ======Compute delta r for all points======
for (i=0; i<Npoints; ++i)
{
const real8 line = Position(i,0);
const real8 pixel = Position(i,1);
// ______ Compute azimuth/range time of this pixel______
//const real8 m_trange = pix2tr(pixel,master.t_range1,master.rsr2x);
const real8 m_trange = master.pix2tr(pixel);
const real8 m_tazi = master.line2ta(line); // added by FvL
// ______ Compute xyz of this point P from position in image ______
cn P; // point, returned by lp2xyz
lp2xyz(line,pixel,ellips,master,masterorbit,
P,MAXITER,CRITERPOS);
// ______ Compute xyz for slave satelite from P ______
real8 s_tazi; // returned
real8 s_trange; // returned
xyz2t(s_tazi,s_trange,slave,
slaveorbit,
P,MAXITER,CRITERTIM);
// ______Compute delta range ~= phase______
// ______ real8 dr = dist(m_possat,pospoint) - dist(s_possat,pospoint);
// ______ real8 phase = -pi4*(dr/LAMBDA);
// ______ dr == M-S want if no flatearth M-S - flatearth = M-S-(M-S)=0
// ______ phase == -4pi*dr/lambda == 4pi*(S-M)/lambda
// BK: 24-9: actually defined as: phi = +pi4/lambda * (r1-r2) ???
// real8 phase = pi4*((dist(s_possat,pospoint)-dist(m_possat,pospoint))/LAMBDA);
//y(i,0) = pi4*((dist(s_possat,pospoint)-dist(m_possat,pospoint))/LAMBDA);
//y(i,0) = pi4divlam*(s_possat.dist(pospoint)-m_possat.dist(pospoint));
//y(i,0) = minpi4cdivlam * (m_trange - s_trange);
y(i,0) = m_minpi4cdivlam*m_trange - s_minpi4cdivlam*s_trange;
DEBUG << "l=" << line << " p=" << pixel
<< " t1=" << m_trange << " t2=" << s_trange
<< " fe=" << y(i,0) << " [rad]";
DEBUG.print();
// ____________________________________________________________________________________
// _____________ Vector with h2ph factors for random number of points by FvL __________
//_____________________________________________________________________________________
cn Psat_master = masterorbit.getxyz(m_tazi);
cn Psat_slave = slaveorbit.getxyz(s_tazi);
real8 B = Psat_master.dist(Psat_slave); // abs. value
// const real8 Bpar = P.dist(M) - P.dist(S); // sign ok
real8 Bpar = SOL*(m_trange-s_trange); // sign ok
// ______ if (MP>SP) then S is to the right of slant line, then B perp is positive.
cn r1 = Psat_master.min(P);
cn r2 = Psat_slave.min(P);
// real8 theta = Psat_master.angle(r1); // look angle
real8 theta = P.angle(r1); // incidence angle
real8 theta_slave = P.angle(r2); // incidence angle slave
real8 Bperp = (theta > theta_slave ) ? // sign ok
sqrt(sqr(B)-sqr(Bpar)) :
-sqrt(sqr(B)-sqr(Bpar)) ;
y_h2ph(i,0) = Bperp/(m_trange*SOL*sin(theta));
// ____________________________________________________________________________________
// _____________ End added part by FvL ________________________________________________
//_____________________________________________________________________________________
// ______Set up system of equations______
// ______Order unknowns: A00 A10 A01 A20 A11 A02 A30 A21 A12 A03 for degree=3______
// ______ normalize data [-2,2] ______
real8 posL = normalize(line,minL,maxL);
real8 posP = normalize(pixel,minP,maxP);
index = 0;
for (j=0; j<=DEGREE; j++)
{
for (k=0; k<=j; k++)
{
A(i,index) = pow(posL,real8(j-k)) * pow(posP,real8(k));
index++;
}
}
}
// ======Compute polynomial for these phases (LS)======
// ______Compute Normalmatrix, rghthandside______
matrix<real8> N = matTxmat(A,A);
matrix<real8> rhs = matTxmat(A,y);
matrix<real8> rhs_h2ph = matTxmat(A,y_h2ph); // Added by FvL, same A matrix can be used
// ______Compute solution______
matrix<real8> Qx_hat = N;
choles(Qx_hat); // Cholesky factorisation normalmatrix
solvechol(Qx_hat,rhs); // Estimate of unknowns in rhs
solvechol(Qx_hat,rhs_h2ph); // Estimate of unknowns in rhs_h2ph, added by FvL
invertchol(Qx_hat); // Covariance matrix
// ______Test inverse______
for (i=0; i<Qx_hat.lines(); i++)
for (j=0; j<i; j++)
Qx_hat(j,i) = Qx_hat(i,j);// repair Qx_hat
const real8 maxdev = max(abs(N*Qx_hat-eye(real8(Qx_hat.lines()))));
INFO << "flatearth: max(abs(N*inv(N)-I)) = " << maxdev;
INFO.print();
if (maxdev > .01)
{
ERROR << "Deviation too large. Decrease degree or number of points?";
PRINT_ERROR(ERROR.get_str())
throw(some_error);
}
else if (maxdev > .001)
{
WARNING << "Deviation quite large. Decrease degree or number of points?";
WARNING.print();
}
else
{
INFO.print("Deviation is OK.");
}
// ______Some other stuff, scale is ok______
// matrix<real8> Qy_hat = A * (matxmatT(Qx_hat,A));
matrix<real8> y_hat = A * rhs;
matrix<real8> y_hat_h2ph = A * rhs_h2ph; // added by FvL
matrix<real8> e_hat = y - y_hat;
matrix<real8> e_hat_h2ph = y_h2ph - y_hat_h2ph; // added by FvL
// ______Overall model test (variance factor)______
// ... ?
// ______ Wrap offset ______
// BK 30/9/99 do not do this, later absolute ref. phase is used.
// in s2h rodriguez for example.
// it does not change anything for compinterfero etc.
// rhs(0,0) = remainder(rhs(0,0),2*PI);
// ______Write results to file______
ofstream scratchlogfile("scratchlogflat", ios::out | ios::trunc);
bk_assert(scratchlogfile,"flatearth: scratchlogflat",__FILE__,__LINE__);
scratchlogfile << "\n\n*******************************************************************"
<< "\n* FLATEARTH: "
//<< "\n*_Start_" << processcontrol[pr_i_comprefpha]
<< "\n*******************************************************************"
<< "\nDegree_flat:\t" << DEGREE
<< "\nEstimated coefficients:\n"
<< "\nx_hat \tstd:\n";
for (i=0; i<Nunk; i++)
scratchlogfile
<< setiosflags(ios::fixed)
<< setiosflags(ios::showpoint)
<< setiosflags(ios::right)
<< setw(8) << setprecision(4)
<< rhs(i,0) << " \t" << sqrt(Qx_hat(i,i)) << endl;
// ___________________ added by FvL _________________________________________________________
scratchlogfile << "\n" << "\nDegree_h2ph:\t" << DEGREE
<< "\nEstimated coefficients:\n"
<< "\nx_hat \tstd:\n";
for (i=0; i<Nunk; i++)
scratchlogfile
<< setiosflags(ios::fixed)
<< setiosflags(ios::showpoint)
<< setiosflags(ios::right)
<< setw(8) << setprecision(4)
<< rhs_h2ph(i,0) << " \t" << sqrt(Qx_hat(i,i)) << endl;
// ___________________ end added by FvL _________________________________________________________
scratchlogfile << "\nCovariance matrix estimated parameters:"
<< "\n---------------------------------------\n";
for (i=0; i<Nunk; i++)
{
for (j=0; j<Nunk; j++)
{
scratchlogfile
<< setiosflags(ios::fixed)
<< setiosflags(ios::showpoint)
<< setiosflags(ios::right)
<< setw(8) << setprecision(4)
<< Qx_hat(i,j) << " ";
}
scratchlogfile << endl;
}
scratchlogfile << "\nMaximum deviation N*inv(N):"
<< setiosflags(ios::scientific)
<< maxdev
<< "\nSome more info for each observation:"
<< "\nline \tpixel \tobs \t\tobs_hat \t\t err_hat\n";
for (i=0; i<Npoints; i++)
scratchlogfile << Position(i,0) << "\t" << Position(i,1) << "\t"
<< y(i,0) << "\t" << y_hat(i,0) << "\t"
<< setiosflags(ios::fixed)
<< setiosflags(ios::showpoint)
<< setiosflags(ios::right)
<< setw(8) << setprecision(4)
<< e_hat(i,0) << "\n";
scratchlogfile << "\nMaximum absolute error: \t\t\t"
<< max(abs(e_hat))
<< "\n*******************************************************************\n";
// scratchlogfile.close(); // commented out for coordinates dumped below
ofstream scratchresfile("scratchresflat", ios::out | ios::trunc);
bk_assert(scratchresfile,"flatearth: scratchresflat",__FILE__,__LINE__);
scratchresfile.setf(ios::scientific, ios::floatfield);
scratchresfile.setf(ios::right, ios::adjustfield);
scratchresfile.precision(8);
scratchresfile.width(18);
scratchresfile << "\n\n*******************************************************************"
//<< "\n*_Start_flat_earth"
<< "\n*_Start_" << processcontrol[pr_i_comprefpha]
<< "\n*******************************************************************"
<< "\nDegree_flat:\t" << DEGREE
<< "\nEstimated_coefficients_flatearth:\n";
int32 coeffL = 0;
int32 coeffP = 0;
for (i=0; i<Nunk; i++)
{
if (rhs(i,0) < 0.)
scratchresfile << rhs(i,0);
else
scratchresfile << " " << rhs(i,0);
// ______ Add coefficient number behind value ______
scratchresfile << " \t" << coeffL << " " << coeffP << "\n";
coeffL--;
coeffP++;
if (coeffL == -1)
{
coeffL = coeffP;
coeffP = 0;
}
}
//_________ added by FvL _______________________________________
scratchresfile << "\n" << "\nDegree_h2ph:\t" << DEGREE
<< "\nEstimated_coefficients_h2ph:\n";
coeffL = 0;
coeffP = 0;
for (i=0; i<Nunk; i++)
{
if (rhs_h2ph(i,0) < 0.)
scratchresfile << rhs_h2ph(i,0);
else
scratchresfile << " " << rhs_h2ph(i,0);
// ______ Add coefficient number behind value ______
scratchresfile << " \t" << coeffL << " " << coeffP << "\n";
coeffL--;
coeffP++;
if (coeffL == -1)
{
coeffL = coeffP;
coeffP = 0;
}
}
//_________ end added by FvL _______________________________________
scratchresfile << "*******************************************************************"
<< "\n* End_" << processcontrol[pr_i_comprefpha] << "_NORMAL"
<< "\n*******************************************************************\n";
// ====== Compute coordinates of corners of interferogram here ======
// ______ (though better place this somewhere else ....
real8 phi,lambda,height; // rad, returned
lp2ell(interferogram.win.linelo,interferogram.win.pixlo,
ellips,master,masterorbit,
phi,lambda,height,MAXITER,CRITERPOS);
INFO << "Coordinates of corner interferogram: "
<< interferogram.win.linelo << ", " << interferogram.win.pixlo
<< " = " << rad2deg(phi) << ", " << rad2deg(lambda);
INFO.print();
scratchlogfile << "\n\n********************************************"
<< "\n* [Lat_Long] coordinates of crop [START] *"
<< "\n********************************************\n";
scratchlogfile << "\nCoords_of_ifg_corner [l,p] : [phi,lam]: "
<< interferogram.win.linelo << " , " << interferogram.win.pixlo
<< " = " << rad2deg(phi) << " , " << rad2deg(lambda);
lp2ell(interferogram.win.linehi,interferogram.win.pixlo,
ellips,master,masterorbit,
phi,lambda,height,MAXITER,CRITERPOS);
INFO << "\nCoordinates of corner interferogram: "
<< interferogram.win.linehi << ", " << interferogram.win.pixlo
<< " = " << rad2deg(phi) << ", " << rad2deg(lambda);
INFO.print();
scratchlogfile << "\nCoords_of_ifg_corner [l,p] : [phi,lam]: "
<< interferogram.win.linehi << " , " << interferogram.win.pixlo
<< " = " << rad2deg(phi) << " , " << rad2deg(lambda);
lp2ell(interferogram.win.linelo,interferogram.win.pixhi,
ellips,master,masterorbit,
phi,lambda,height,MAXITER,CRITERPOS);
INFO << "\nCoordinates of corner interferogram: "
<< interferogram.win.linelo << ", " << interferogram.win.pixhi
<< " = " << rad2deg(phi) << ", " << rad2deg(lambda);
INFO.print();
scratchlogfile << "\nCoordinates of ifgs [low,hi] corner [l,p] : [phi,lam]: "
<< interferogram.win.linelo << " , " << interferogram.win.pixhi
<< " = " << rad2deg(phi) << " , " << rad2deg(lambda);
lp2ell(interferogram.win.linehi,interferogram.win.pixhi,
ellips,master,masterorbit,
phi,lambda,height,MAXITER,CRITERPOS);
INFO << "\nCoordinates of corner interferogram: "
<< interferogram.win.linehi << ", " << interferogram.win.pixhi
<< " = " << rad2deg(phi) << ", " << rad2deg(lambda);
INFO.print();
scratchlogfile << "\nCoords_of_ifg_corner [l,p] : [phi,lam]: "
<< interferogram.win.linehi << " , " << interferogram.win.pixhi
<< " = " << rad2deg(phi) << " , " << rad2deg(lambda);
scratchlogfile << "\n\n******************************************"
<< "\n* [Lat_Long] coordinates of crop [STOP] *"
<< "\n******************************************\n";
// ______Tidy up______
scratchresfile.close(); // close res file
scratchlogfile.close(); // close log.out file
} // END flatearth
/****************************************************************
* demassist *
* *
* Coregistration based on DEM (SRTM) *
* DEM on equiangular grid (lat/lon) assumed *
* DEM seems stored from North to South *
* *
* Freek van Leijen, Liu Guang, Mahmut Arikan, 21-Sep-2007 *
****************************************************************/
void demassist(
const input_gen &generalinput,
const input_ell &ellips,
const input_demassist &demassistinput,
const slcimage &master,
const slcimage &slave,
orbit &masterorbit,
orbit &slaveorbit)
{
TRACE_FUNCTION("demassist (FvL 21-SEP-2007)")
const string STEP="DAC: ";
const int32 MAXITER = 10;
const real8 CRITERPOS = 1e-6;
const real8 CRITERTIM = 1e-10;
const real8 lat0file = demassistinput.demlatleftupper;// first pix on disk w02090
const real8 lon0file = demassistinput.demlonleftupper;// first pix on disk
const real8 DEMdeltalat = demassistinput.demdeltalat; // in radians
const real8 DEMdeltalon = demassistinput.demdeltalon; // in radians
const int32 numberoflonpixels = demassistinput.demcols; // NCOLS on file
const int32 numberoflatpixels = demassistinput.demrows; // NROWS on file
const real8 NODATA = demassistinput.demnodata; // (BK 4 may 2001)
const bool outputdemi = specified(demassistinput.fodemi);// if spec. then output
const bool outputrefdemhei = specified(demassistinput.forefdemhei);
////const real8 m_min4picdivlam = (-4.0*PI*SOL)/master.wavelength;
////const real8 s_min4picdivlam = (-4.0*PI*SOL)/slave.wavelength;
const real8 latNfile = lat0file-DEMdeltalat*(numberoflatpixels-1); // upper=max. lat value
const real8 lonNfile = lon0file+DEMdeltalon*(numberoflonpixels-1); // left=min. lon value
// ______ Extra info ______
INFO << "DEM input: w/e/s/n: \t"
<< rad2deg(lon0file) << "/" << rad2deg(lonNfile) << "/"
<< rad2deg(latNfile) << "/" << rad2deg(lat0file);
INFO.print();
// ______ Get corners of master (approx) to select DEM ______
// ______ in radians (if height were zero)______
real8 extralat = (1.5*DEMdeltalat + deg2rad(4.0/25.0));
real8 extralong = (1.5*DEMdeltalon + deg2rad(4.0/25.0));
real8 phimin;
real8 phimax;
real8 lambdamin;
real8 lambdamax;
int32 indexphi0DEM;
int32 indexphiNDEM;
int32 indexlambda0DEM;
int32 indexlambdaNDEM;
getcorners(master.currentwindow.linelo,master.currentwindow.linehi,
master.currentwindow.pixlo,master.currentwindow.pixhi,
extralat,extralong,lat0file,lon0file,
DEMdeltalat,DEMdeltalon,numberoflatpixels,numberoflonpixels,
ellips,master,masterorbit,phimin,phimax,lambdamin,lambdamax,
indexphi0DEM,indexphiNDEM,indexlambda0DEM,indexlambdaNDEM);
// ______ Extra info ______
INFO << "DEM input required: w/e/s/n: \t"
<< rad2deg(lambdamin) << "/" << rad2deg(lambdamax) << "/"
<< rad2deg(phimin) << "/" << rad2deg(phimax);
INFO.print();
INFO << "For window (l0,lN,p0,pN): \t"
<< master.currentwindow.linelo << " "
<< master.currentwindow.linehi << " "
<< master.currentwindow.pixlo << " "
<< master.currentwindow.pixhi;
INFO.print();
// ______ Check corners of DEM ______
// check if DEM is appropriate for master crop
// DEM should at least partially cover master crop
// note: phi is [90:-90]
if (phimax <= latNfile)// DEM is more north than master
{
ERROR << "master crop outside DEM: most South latitude: " << rad2deg(latNfile)
<< " [deg]; master crop requires: " << rad2deg(phimax)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
//throw(some_error);
}
// DEM is more south than master crop
if (phimin >= lat0file)// largest latitude at first line of file
{
ERROR << "master crop outside DEM: most North latitude: " << rad2deg(lat0file)
<< " [deg]; master crop requires: " << rad2deg(phimax)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
//throw(some_error);
}
if (lambdamax <= lon0file)
{
ERROR << "master crop outside DEM: most West longitude: " << rad2deg(lon0file)
<< " [deg]; master crop window requires: " << rad2deg(lambdamax)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
//throw(some_error);
}
if (lambdamin >= lonNfile)
{
ERROR << "master crop outside DEM: most East longitude: " << rad2deg(lonNfile)
<< " [deg]; master crop window requires: " << rad2deg(lambdamin)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
//throw(some_error);
}
//===================================================================
//============ First loop: radarcode DEM ============================
//============ (DEM geometry) ============================
//===================================================================
int32 numvalid = 0;// number of good values, not NODATA in buffer
int32 numNODATA = 0;// number of NODATA values in buffer
real8 meancroppedDEM = 0.0;// to detect byte order problems, formats
real8 min_input_dem = 100000.0;// stats
real8 max_input_dem = -100000.0;// stats
// ______ Compute buffer size radarcoding DEM______
const real8 BUFFERMEMSIZE = generalinput.memory;// Bytes
int32 NcolsDEM = indexlambdaNDEM-indexlambda0DEM+1;
int32 NrowsDEM = indexphiNDEM-indexphi0DEM+1;
const uint32 NcolsDEMlog = NcolsDEM; // since NcolsDEM updated after getcorners() call.
const uint32 NrowsDEMlog = NrowsDEM;
const real8 Nrows_possible_DEM = BUFFERMEMSIZE / (5*8*NcolsDEM);
int32 bufferlines = int32(ceil(Nrows_possible_DEM)); // [MA] checked ok. Sinces SLC is not multilooked, see comprefdem for solution
if (bufferlines>NrowsDEM) bufferlines=NrowsDEM;
int32 numfullbuffers = NrowsDEM / bufferlines;
int32 restlines = NrowsDEM % bufferlines;
int32 extrabuffer = (restlines == 0) ? 0 : 1;
// ______ Extra info ______
INFO << "DEM output total pixels: " << NcolsDEM;
INFO.print();
INFO << "DEM output total lines : " << NrowsDEM;
INFO.print();
INFO << "Radar coding of DEM in: " << numfullbuffers << " buffers of "
<< bufferlines << " lines and " << extrabuffer << " extra buffer of "
<< restlines << " lines.";
INFO.print();
// ______ Open (temporary) output files ______
// DEM heights
ofstream demofile;
openfstream(demofile,demassistinput.fodem,generalinput.overwrit); // dem_crop radarcoded
bk_assert(demofile,demassistinput.fodem,__FILE__,__LINE__);
// master line coordinates of DEM
ofstream masterdemlineoutfile("dac_m_demline.temp", ios::out | ios::trunc);
bk_assert(masterdemlineoutfile,"dac_m_demline.temp",__FILE__,__LINE__);
// master pixel coordinates of DEM
ofstream masterdempixeloutfile("dac_m_dempixel.temp", ios::out | ios::trunc);
bk_assert(masterdempixeloutfile,"dac_m_dempixel.temp",__FILE__,__LINE__);
// delta line coordinates of DEM ( slave-master )
ofstream deltademlineoutfile("dac_delta_demline.temp", ios::out | ios::trunc);
bk_assert(deltademlineoutfile,"dac_delta_demline.temp",__FILE__,__LINE__);
// delta pixel coordinates of DEM
ofstream deltadempixeloutfile("dac_delta_dempixel.temp", ios::out | ios::trunc);
bk_assert(deltadempixeloutfile,"dac_delta_dempixel.temp",__FILE__,__LINE__);
// ______ DEM loop per buffer ______
register int32 j,i;// DEM index grid counter, register j first to ensure allocation
for (register int32 buffer=0; buffer<numfullbuffers+extrabuffer; ++buffer)
{
// Determine indices for buffer
const int32 indexphi0BUFFER = indexphi0DEM+buffer*bufferlines;
const int32 blines = (buffer == numfullbuffers) ? restlines : bufferlines;
const int32 indexphiNBUFFER = indexphi0BUFFER+(blines-1);
matrix<real4> DEM(blines,NcolsDEM);
// ______ Extra info ______
PROGRESS << STEP << "Buffer# [l0:lN, p0:pN]: " << buffer+1 << " ["
<< indexphi0BUFFER << ": " << indexphiNBUFFER << ", "
<< indexlambda0DEM << ": " << indexlambdaNDEM << "]";
PROGRESS.print();
// ______ lat/lon for first pixel in matrix read from file ______
// ______ upper is max. latitude, left is min. longitude ______
const real8 upperleftphi = lat0file-indexphi0BUFFER*DEMdeltalat;
const real8 upperleftlambda = lon0file+indexlambda0DEM*DEMdeltalon;
window zerooffset (0,0,0,0);
window winfromfile (indexphi0BUFFER,indexphiNBUFFER,
indexlambda0DEM,indexlambdaNDEM);
// ______ Read in grdfile of DEM in matrix R4 (raw data, no header) _______
// ______ added formats (BK 4-May-2001) ______
PROGRESS << STEP << "Reading crop of DEM for buffer: " << buffer+1;
PROGRESS.print();
DEBUG.print("Reading input DEM into real4 matrix (buffer).");
switch (demassistinput.iformatflag)
{
// ______ Read as short BE, then convert to host order ______
case FORMATI2_BIGENDIAN:
{
matrix<int16> DEMi2(blines,NcolsDEM);
readfile(DEMi2,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = real4(ntohs(DEMi2(iii,jjj)));// cast to real4
DEMi2.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: SHORT SIGNED INT.");
break;
}
case FORMATI2:
{
matrix<int16> DEMi2(blines,NcolsDEM);
readfile(DEMi2,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = DEMi2(iii,jjj);// cast to real4
DEMi2.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: SHORT SIGNED INT.");
break;
}
case FORMATR4:
readfile(DEM,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
INFO.print("Read crop of input DEM: format: REAL4.");
break;
case FORMATR8:
{
matrix<real8> DEMr8(blines,NcolsDEM);
readfile(DEMr8,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = DEMr8(iii,jjj);// cast to real4
DEMr8.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: REAL8.");
break;
}
default:
PRINT_ERROR("totally impossible, checked input.")
//throw(unhandled_case_error);
}
// ----- Loop over DEM for stats ------------------------
real8 min_dem_buffer = 100000.0;
real8 max_dem_buffer = -100000.0;
for (i=0; i<DEM.lines(); ++i)
{
// ----- Loop over oversampled matrix in x ------
for (j=0; j<DEM.pixels(); ++j)
{
if(DEM(i,j)!=NODATA)
{
numvalid++;
meancroppedDEM += DEM(i,j);// divide by numvalid later
if (DEM(i,j)<min_dem_buffer) min_dem_buffer=DEM(i,j);//buffer
if (DEM(i,j)>max_dem_buffer) max_dem_buffer=DEM(i,j);// stats
}
else
{
numNODATA++;
}
}//loop dem for stats
}//loop dem for stats
min_input_dem = min(min_input_dem,min_dem_buffer);//global stats
max_input_dem = max(max_input_dem,max_dem_buffer);//global stats
// ====== Radarcoding DEM ==============================
// ______ DEM contains values from leftupper with ______
// ______ spacing (DEMdeltalat,DEMdeltalon) ______
// ______ Transform DEM to l,p,refphase ______
PROGRESS.print("Converting DEM to radar system for this buffer.");
const int32 NpointsDEM = DEM.size();
const int32 NpixelsDEM = DEM.pixels();
// ______ Extra info ______
INFO << "Number of points in DEM: "
<< NpointsDEM;
INFO.print();
matrix<real8> masterDEMline(DEM.lines(),DEM.pixels());
matrix<real8> masterDEMpixel(DEM.lines(),DEM.pixels());
matrix<real8> deltaDEMline(DEM.lines(),DEM.pixels());
matrix<real8> deltaDEMpixel(DEM.lines(),DEM.pixels());
// --- Loop DEM ---
real8 phi,lambda,height,m_l,m_p,s_l,s_p;
phi = upperleftphi;
for (i=0; i<DEM.lines(); ++i)
{
if ((i%100)==0)
{
// ______ Extra info ______
PROGRESS << STEP << "Radarcoding buffer: " << buffer+1 << "of" << numfullbuffers + extrabuffer << " DEM line: " << i << " ("
<< floor(.5+(100.*real8(i)/real8(DEM.lines())))
<< "%)";
PROGRESS.print();
}
lambda = upperleftlambda;
for (j=0; j<DEM.pixels(); ++j)
{
height = DEM(i,j);
ell2lp(m_l,m_p,ellips,master,masterorbit,phi,lambda,height,MAXITER,CRITERTIM);
ell2lp(s_l,s_p,ellips,slave,slaveorbit,phi,lambda,height,MAXITER,CRITERTIM);
masterDEMline(i,j) = m_l;
masterDEMpixel(i,j) = m_p;
deltaDEMline(i,j) = s_l-m_l;
deltaDEMpixel(i,j) = s_p-m_p;
lambda += DEMdeltalon;
} // loop DEM pixels
// ______ update latitude of next line ______
phi -= DEMdeltalat; // upper left is max. value
} // loop DEM lines
// Write results to output files
PROGRESS << STEP << "Writing radar coded DEM to file, buffer: " << buffer+1 << " of " << numfullbuffers + extrabuffer ;
PROGRESS.print();
demofile << DEM;
masterdemlineoutfile << masterDEMline;
masterdempixeloutfile << masterDEMpixel;
deltademlineoutfile << deltaDEMline;
deltadempixeloutfile << deltaDEMpixel;
masterDEMline.resize(1,1); //deallocate
masterDEMpixel.resize(1,1); //deallocate
deltaDEMline.resize(1,1); //deallocate
deltaDEMpixel.resize(1,1); //deallocate
DEM.resize(1,1); //deallocate
} // buffer loop
demofile.close();
masterdemlineoutfile.close();
masterdempixeloutfile.close();
deltademlineoutfile.close();
deltadempixeloutfile.close();
//===================================================================
//============ End first loop: radarcode DEM ========================
//============ (DEM geometry) ============================
//===================================================================
//===================================================================
//============ Second loop: interpolation =============
//============ (radar geometry) =============
//===================================================================
INFO << STEP << "Start interpolation...";
INFO.print();
// ______ Line/pixel of first point in original master coordinates ______
const int32 Nlinesml = master.currentwindow.lines();
const int32 Npixelsml = master.currentwindow.pixels();
const real8 veryfirstline = real8(master.currentwindow.linelo);
const real8 verylastline = real8(master.currentwindow.linehi);
const real8 firstpixel = real8(master.currentwindow.pixlo);
const real8 lastpixel = real8(master.currentwindow.pixhi);
//Determine range-azimuth spacing ratio, needed for proper triangulation
cn P1, P2 , P3, P4;
lp2xyz(veryfirstline,firstpixel,ellips,master,masterorbit,
P1,MAXITER,CRITERPOS);
lp2xyz(veryfirstline,lastpixel,ellips,master,masterorbit,
P2,MAXITER,CRITERPOS);
lp2xyz(verylastline,firstpixel,ellips,master,masterorbit,
P3,MAXITER,CRITERPOS);
lp2xyz(verylastline,lastpixel,ellips,master,masterorbit,
P4,MAXITER,CRITERPOS);
const real8 r_spacing = ( (P1.min(P2)).norm() + (P3.min(P4)).norm() ) / 2 /(lastpixel - firstpixel) ;
const real8 az_spacing = ( (P1.min(P3)).norm() + (P2.min(P4)).norm() ) /2 /(verylastline - veryfirstline);
const real8 r_az_ratio = r_spacing/az_spacing;
INFO << "Master azimuth spacing: " << az_spacing;
INFO.print();
INFO << "Master range spacing: " << r_spacing;
INFO.print();
INFO << "Range-azimuth spacing ratio: " << r_az_ratio;
INFO.print();
// ______ Compute buffer size interpolation______
const real8 Nlinesml_possible = BUFFERMEMSIZE / (6*8*Npixelsml);
bufferlines = int32(ceil(Nlinesml_possible));
if (bufferlines > Nlinesml) bufferlines=Nlinesml;
numfullbuffers = Nlinesml / bufferlines;
restlines = Nlinesml % bufferlines; // the number of lines in extra buffer
extrabuffer = (restlines == 0) ? 0 : 1;
// ______ Extra info ______
INFO << "Interpolation in: " << numfullbuffers << " buffers of "
<< bufferlines << " lines and " << extrabuffer << " extra buffer of "
<< restlines << " lines.";
INFO.print();
// ______ Open output files ______
ofstream deltalineofile("dac_delta_line.raw", ios::out | ios::trunc);
bk_assert(deltalineofile,"dac_delta_line.raw",__FILE__,__LINE__);
ofstream deltapixelofile("dac_delta_pixel.raw", ios::out | ios::trunc);
bk_assert(deltapixelofile,"dac_delta_pixel.raw",__FILE__,__LINE__);
// if request for height in radar coordinates l,p
ofstream refdemheiofile;
if (outputrefdemhei==true)
{
openfstream(refdemheiofile,demassistinput.forefdemhei,generalinput.overwrit);
bk_assert(refdemheiofile,demassistinput.forefdemhei,__FILE__,__LINE__);
}
// ______ interpolation loop per buffer ______
for (register int32 buffer = 0; buffer < numfullbuffers + extrabuffer; ++buffer)
{
// Determine indices for buffer
const int32 blines = (buffer == numfullbuffers) ? restlines : bufferlines;
const real8 firstline_buffer = veryfirstline+buffer*bufferlines;
const real8 lastline_buffer = firstline_buffer+blines-1;
// ______ Extra info ______
PROGRESS << STEP << "Interpolation buffer: " << buffer+1 << "of" << numfullbuffers + extrabuffer << " [l0:lN, p0:pN]: " << " ["
<< firstline_buffer << ": " << lastline_buffer << ", "
<< firstpixel << ": " << lastpixel << "]";
PROGRESS.print();
// Get corners of buffer
real8 phimin_az;
real8 phimax_az;
real8 lambdamin_az;
real8 lambdamax_az;
getcorners(firstline_buffer,lastline_buffer,
firstpixel,lastpixel,
extralat,extralong,phimax,lambdamin,
DEMdeltalat,DEMdeltalon,NrowsDEM,NcolsDEM,
ellips,master,masterorbit,phimin_az,phimax_az,lambdamin_az,lambdamax,
indexphi0DEM,indexphiNDEM,indexlambda0DEM,indexlambdaNDEM);
window zerooffset (0,0,0,0);
window winfromfile (indexphi0DEM,indexphiNDEM,
indexlambda0DEM,indexlambdaNDEM);
const int32 NrowsDEM_buffer = indexphiNDEM-indexphi0DEM+1;
const int32 NcolsDEM_buffer = indexlambdaNDEM-indexlambda0DEM+1;
PROGRESS << STEP << "Reading input for interpolation buffer: " << buffer+1 << "of" << numfullbuffers + extrabuffer;
PROGRESS.print();
// read x,y
matrix<real8> DEMline_buffer(NrowsDEM_buffer,NcolsDEM_buffer);
matrix<real8> DEMpixel_buffer(NrowsDEM_buffer,NcolsDEM_buffer);
readfile(DEMline_buffer,"dac_m_demline.temp",NrowsDEM,winfromfile,zerooffset);
readfile(DEMpixel_buffer,"dac_m_dempixel.temp",NrowsDEM,winfromfile,zerooffset);
// read z (multiple, number can easily be increased, e.g. simulated intensity)
int32 Nz = 2; //number of z
matrix<real8> input_buffer(NrowsDEM_buffer *Nz ,NcolsDEM_buffer);
matrix<real8> temp_input_buffer(NrowsDEM_buffer,NcolsDEM_buffer);
if (outputrefdemhei==true)
{
Nz += 1;
input_buffer.resize(NrowsDEM_buffer *Nz ,NcolsDEM_buffer);
}
readfile(temp_input_buffer,"dac_delta_demline.temp",NrowsDEM,winfromfile,zerooffset);
input_buffer.setdata(0, 0, temp_input_buffer);
readfile(temp_input_buffer,"dac_delta_dempixel.temp",NrowsDEM,winfromfile,zerooffset);
input_buffer.setdata(NrowsDEM_buffer, 0, temp_input_buffer);
Nz = 2;
if (outputrefdemhei==true)
{
Nz += 1;
/// i would like to use real4, test later on
matrix<real4> dem_input(NrowsDEM_buffer,NcolsDEM_buffer);
readfile(dem_input,demassistinput.fodem,NrowsDEM,winfromfile,zerooffset);
for (register int32 i =0 ; i < NrowsDEM_buffer ; i ++)
for(register int32 j = 0; j < NcolsDEM_buffer; j++)
temp_input_buffer(i,j) = real8(dem_input(i,j));
input_buffer.setdata(NrowsDEM_buffer * (Nz-1), 0, temp_input_buffer);
}
// initialize output array
Nz = 2;
matrix<real8> output_buffer(blines * Nz, Npixelsml);
if (outputrefdemhei==true)
{
Nz += 1;
output_buffer.resize(blines * Nz, Npixelsml);
}
// interpolation
griddatalinear(DEMline_buffer,DEMpixel_buffer,input_buffer,
firstline_buffer,lastline_buffer,firstpixel,lastpixel,
1,1,r_az_ratio,0,NODATA,output_buffer);
deltalineofile << output_buffer(window(0, blines - 1, 0, Npixelsml -1 ));
deltapixelofile << output_buffer(window(blines , 2 * blines - 1, 0, Npixelsml -1 ));
Nz = 2;
if (outputrefdemhei==true)
{
Nz += 1;
refdemheiofile << output_buffer(window((Nz-1) * blines,Nz * blines - 1, 0, Npixelsml -1 ));
}
DEMline_buffer.resize(1,1);
DEMpixel_buffer.resize(1,1);
input_buffer.resize(1,1);
temp_input_buffer.resize(1,1);
output_buffer.resize(1,1);
} // end loop azimuth direction
INFO << "Closing output files";
INFO.print();
deltalineofile.close();
deltapixelofile.close();
if (outputrefdemhei==true) // Radarcoded DEM
refdemheiofile.close();
//===================================================================
//============ End second loop: interpolation =============
//============ (radar geometry) =============
//===================================================================
//===================================================================
//============ Determine inverse transformation =============
//============ (slave corners only, needed for overlap) =============
//===================================================================
real8 line, pixel;
real8 deltaline_slave00,deltapixel_slave00,
deltaline_slave0N,deltapixel_slave0N,
deltaline_slaveN0,deltapixel_slaveN0,
deltaline_slaveNN,deltapixel_slaveNN;
real8 phimin_az,phimax_az,lambdamin_az,lambdamax_az;
for (register int16 corner = 0 ; corner < 4 ; corner ++)
{
PROGRESS << "Radarcoding slave corner: " << corner+1;
PROGRESS.print();
switch (corner)
{
case 0:
{
line=slave.currentwindow.linelo;
pixel=slave.currentwindow.pixlo;
break;
}
case 1:
{
line=slave.currentwindow.linelo;
pixel=slave.currentwindow.pixhi;
break;
}
case 2:
{
line=slave.currentwindow.linehi;
pixel=slave.currentwindow.pixlo;
break;
}
case 3:
{
line=slave.currentwindow.linehi;
pixel=slave.currentwindow.pixhi;
break;
}
default:
PRINT_ERROR("totally impossible, checked input.");
}
//use getcorners with line,line,pixel,pixel for single point
getcorners(line,line,pixel,pixel,
extralat,extralong,lat0file,lon0file,
DEMdeltalat,DEMdeltalon,numberoflatpixels,numberoflonpixels,
ellips,slave,slaveorbit,
phimin_az,phimax_az,lambdamin_az,lambdamax,
indexphi0DEM,indexphiNDEM,indexlambda0DEM,indexlambdaNDEM);
NcolsDEM = indexlambdaNDEM-indexlambda0DEM+1;
NrowsDEM = indexphiNDEM-indexphi0DEM+1;
const real8 upperleftphi = lat0file-indexphi0DEM*DEMdeltalat;
const real8 upperleftlambda = lon0file+indexlambda0DEM*DEMdeltalon;
window zerooffset (0,0,0,0);
window winfromfile (indexphi0DEM,indexphiNDEM,
indexlambda0DEM,indexlambdaNDEM);
// ______ Read in DEM in matrix R4 (raw data, no header) _______
matrix<real4> DEM(NrowsDEM,NcolsDEM);
switch (demassistinput.iformatflag)
{
// ______ Read as short BE, then convert to host order ______
case FORMATI2_BIGENDIAN:
{
matrix<int16> DEMi2(NrowsDEM,NcolsDEM);
readfile(DEMi2,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = real4(ntohs(DEMi2(iii,jjj)));// cast to real4
DEMi2.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: SHORT SIGNED INT.");
break;
}
case FORMATI2:
{
matrix<int16> DEMi2(NrowsDEM,NcolsDEM);
readfile(DEMi2,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = DEMi2(iii,jjj);// cast to real4
DEMi2.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: SHORT SIGNED INT.");
break;
}
case FORMATR4:
readfile(DEM,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
INFO.print("Read crop of input DEM: format: REAL4.");
break;
case FORMATR8:
{
matrix<real8> DEMr8(NrowsDEM,NcolsDEM);
readfile(DEMr8,demassistinput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = DEMr8(iii,jjj);// cast to real4
DEMr8.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: REAL8.");
break;
}
default:
PRINT_ERROR("totally impossible, checked input.");
//throw(unhandled_case_error);
}
// radarcode dem
matrix<real8> slaveDEMline(DEM.lines(),DEM.pixels());
matrix<real8> slaveDEMpixel(DEM.lines(),DEM.pixels());
matrix<real8> deltaDEMline(DEM.lines(),DEM.pixels());
matrix<real8> deltaDEMpixel(DEM.lines(),DEM.pixels());
// --- Loop DEM ---
real8 phi,lambda,height,m_l,m_p,s_l,s_p;
phi = upperleftphi;
for (i=0; i<DEM.lines(); ++i)
{
lambda = upperleftlambda;
for (j=0; j<DEM.pixels(); ++j)
{
height = DEM(i,j);
ell2lp(m_l,m_p,ellips,master,masterorbit,phi,lambda,height,MAXITER,CRITERTIM);
ell2lp(s_l,s_p,ellips,slave,slaveorbit,phi,lambda,height,MAXITER,CRITERTIM);
slaveDEMline(i,j) = s_l;
slaveDEMpixel(i,j) = s_p;
deltaDEMline(i,j) = m_l-s_l;
deltaDEMpixel(i,j) = m_p-s_p;
lambda += DEMdeltalon;
} // loop DEM pixels
// ______ update latitude of next line ______
phi -= DEMdeltalat; // upper left is max. value
} // loop DEM lines
// interpolate to slave corner
matrix<real8> input_buffer(DEM.lines()*2,DEM.pixels());
input_buffer.setdata(0, 0, deltaDEMline);
input_buffer.setdata(DEM.lines(), 0, deltaDEMpixel);
matrix<real8> output_buffer(2,1);
griddatalinear(slaveDEMline,slaveDEMpixel,input_buffer,
line,line,pixel,pixel,
1,1,r_az_ratio,0,NODATA,output_buffer);
switch (corner)
{
case 0:
{
deltaline_slave00 = output_buffer(0,0);
deltapixel_slave00 = output_buffer(1,0);
INFO << "Deltaline_slave00: " << deltaline_slave00;
INFO.print();
INFO << "Deltapixel_slave00: " << deltapixel_slave00;
INFO.print();
break;
}
case 1:
{
deltaline_slave0N = output_buffer(0,0);
deltapixel_slave0N = output_buffer(1,0);
INFO << "Deltaline_slave0N: " << deltaline_slave0N;
INFO.print();
INFO << "Deltapixel_slave0N: " << deltapixel_slave0N;
INFO.print();
break;
}
case 2:
{
deltaline_slaveN0 = output_buffer(0,0);
deltapixel_slaveN0 = output_buffer(1,0);
INFO << "Deltaline_slaveN0: " << deltaline_slaveN0;
INFO.print();
INFO << "Deltapixel_slaveN0: " << deltapixel_slaveN0;
INFO.print();
break;
}
case 3:
{
deltaline_slaveNN = output_buffer(0,0);
deltapixel_slaveNN = output_buffer(1,0);
INFO << "Deltaline_slaveNN: " << deltaline_slaveNN;
INFO.print();
INFO << "Deltapixel_slaveNN: " << deltapixel_slaveNN;
INFO.print();
break;
}
default:
PRINT_ERROR("totally impossible, checked input.");
}
}
//===================================================================
//============ End determine inverse transformation =============
//============ (slave corners only, needed for overlap) =============
//===================================================================
// ====== Write output information ======
char croppeddemi[4*ONE27];
strcpy(croppeddemi,"NO output requested");
if (outputdemi) strcpy(croppeddemi,demassistinput.fodemi);
INFO << "Min. value of input DEM covering master: " << min_input_dem;
INFO.print();
INFO << "Max. value of input DEM covering master: " << max_input_dem;
INFO.print();
ofstream scratchlogfile("scratchlogdemassist", ios::out | ios::trunc);
bk_assert(scratchlogfile,"demassist: scratchlogdemassist",__FILE__,__LINE__);
scratchlogfile
<< "\n*******************************************************************"
<< "\n* " << processcontrol[pr_i_demassist]
<< "\n*******************************************************************"
<< "\n1) DEM source file: \t" << demassistinput.firefdem
<< "\nFormat: \t";
switch (demassistinput.iformatflag)
{
case FORMATI2:
{
scratchlogfile << "SHORT SIGNED INTEGER (HOST ENDIANNESS)";
break;
}
case FORMATI2_BIGENDIAN:
{
scratchlogfile << "SHORT SIGNED INTEGER, BIG ENDIAN";
break;
}
case FORMATR4:
{
scratchlogfile << "REAL4 SIGNED FLOAT";
break;
}
case FORMATR8:
{
scratchlogfile << "REAL8 SIGNED DOUBLE";
break;
}
default:
{
scratchlogfile << "UNKNOWN? IMPOSSIBLE...";
break;
}
}
scratchlogfile
<< "\nByte order: \t" << "check it yourself..."
<< "\nNumber of lines: \t" << numberoflatpixels
<< "\nNumber of pixels: \t" << numberoflonpixels
<< "\nResolution latitude: \t" << rad2deg(DEMdeltalat) << " [deg]"
<< "\nResolution longitude: \t" << rad2deg(DEMdeltalon) << " [deg]"
<< "\nMost West point in input DEM: \t" << rad2deg(lon0file)
<< "\nMost East point in input DEM: \t" << rad2deg(lonNfile)
<< "\nMost South point in input DEM: \t" << rad2deg(latNfile)
<< "\nMost North point in input DEM: \t" << rad2deg(lat0file)
<< "\nMin. value of input DEM covering master: " << min_input_dem
<< "\nMax. value of input DEM covering master: " << max_input_dem
<< "\n2) Output file cropped DEM: \t" << demassistinput.fodem //[FVL
<< "\nFormat: \t" << "REAL4"
<< "\nByte order: \t" << "(same as host)"
<< "\nNumber of lines : \t" << NrowsDEMlog
<< "\nNumber of pixels : \t" << NcolsDEMlog
<< "\nDEM extend w/e/s/n : \t" << rad2deg(lambdamin) << "/"
<< rad2deg(lambdamax) << "/" << rad2deg(phimin) << "/" << rad2deg(phimax)
// << "\nMean value: \t" << meancroppedDEM
<< "\n3) Output file interpolated crop DEM: \t" << croppeddemi
<< "\nFormat: \t" << "REAL4"
<< "\nByte order: \t" << "(same as host)"
<< "\nNumber of lines (multilooked): \t" << Nlinesml
<< "\nNumber of pixels (multilooked): \t" << Npixelsml
<< "\nDeltaline_slave00_dem: \t" << deltaline_slave00
<< "\nDeltapixel_slave00_dem: \t" << deltapixel_slave00
<< "\nDeltaline_slave0N_dem: \t" << deltaline_slave0N
<< "\nDeltapixel_slave0N_dem: \t" << deltapixel_slave0N
<< "\nDeltaline_slaveN0_dem: \t" << deltaline_slaveN0
<< "\nDeltapixel_slaveN0_dem: \t" << deltapixel_slaveN0
<< "\nDeltaline_slaveNN_dem: \t" << deltaline_slaveNN
<< "\nDeltapixel_slaveNN_dem: \t" << deltapixel_slaveNN
<< "\n*******************************************************************\n\n";
scratchlogfile.close();
ofstream scratchresfile("scratchresdemassist", ios::out | ios::trunc);
bk_assert(scratchresfile,"demassist: scratchresdemassist",__FILE__,__LINE__);
scratchresfile
<< "\n\n*******************************************************************"
<< "\n*_Start_" << processcontrol[pr_i_demassist]
<< "\n*******************************************************************";
scratchresfile
<< "\nDEM source file: \t" << demassistinput.firefdem
<< "\nMin. of input DEM: \t" << min_input_dem
<< "\nMax. of input DEM: \t" << max_input_dem
<< "\nFirst_line (w.r.t. original_master): \t"
<< master.currentwindow.linelo
<< "\nLast_line (w.r.t. original_master): \t"
<< master.currentwindow.linehi
<< "\nFirst_pixel (w.r.t. original_master): \t"
<< master.currentwindow.pixlo
<< "\nLast_pixel (w.r.t. original_master): \t"
<< master.currentwindow.pixhi
<< "\nNumber of lines: \t" << Nlinesml
<< "\nNumber of pixels: \t" << Npixelsml
<< "\nDeltaline_slave00_dem: \t" << deltaline_slave00
<< "\nDeltapixel_slave00_dem: \t" << deltapixel_slave00
<< "\nDeltaline_slave0N_dem: \t" << deltaline_slave0N
<< "\nDeltapixel_slave0N_dem: \t" << deltapixel_slave0N
<< "\nDeltaline_slaveN0_dem: \t" << deltaline_slaveN0
<< "\nDeltapixel_slaveN0_dem: \t" << deltapixel_slaveN0
<< "\nDeltaline_slaveNN_dem: \t" << deltaline_slaveNN
<< "\nDeltapixel_slaveNN_dem: \t" << deltapixel_slaveNN
<< "\n*******************************************************************"
<< "\n* End_" << processcontrol[pr_i_demassist] << "_NORMAL"
<< "\n*******************************************************************\n";
scratchresfile.close();
} // END demassist
/****************************************************************
* radarcodedem (a.k.a comprefdem) *
* *
* Compute reference phase based on DEM (SRTM) *
* DEM on equiangular grid (lat/lon) assumed *
* DEM seems stored from North to South *
* *
* Freek van Leijen, 26-Sep-2007 *
****************************************************************/
void radarcodedem(
const input_gen &generalinput,
const input_ell &ellips,
const input_comprefdem &refdeminput,
const slcimage &master,
const slcimage &slave,
const productinfo &interferogram,
orbit &masterorbit,
orbit &slaveorbit)
{
TRACE_FUNCTION("radarcodedem (FvL 26-Sep-2007)")
const string STEP="CRD: ";
const int32 MAXITER = 10;
const real8 CRITERPOS = 1e-6;
const real8 CRITERTIM = 1e-10;
const real8 lat0file = refdeminput.demlatleftupper; // first pix on disk w02090
const real8 lon0file = refdeminput.demlonleftupper; // first pix on disk
const real8 DEMdeltalat = refdeminput.demdeltalat; // in radians
const real8 DEMdeltalon = refdeminput.demdeltalon; // in radians
const int32 numberoflonpixels = refdeminput.demcols; // NCOLS on file
const int32 numberoflatpixels = refdeminput.demrows; // NROWS on file
const real8 NODATA = refdeminput.demnodata; // (BK 4 may 2001)
bool onlyrefphasetopo = !refdeminput.includerefpha; // true: phase DEM w.r.t. ellipsoid
const bool outputdemi = specified(refdeminput.fodemi); // if spec. then output
const bool outputh2ph = specified(refdeminput.foh2ph); // if spec. then output, added by FvL
const bool outputrefdemhei = specified(refdeminput.forefdemhei);
// _____ start added by MA _____
bool mlookedIFG = false; // true: ifg is multilooked
int32 mlL = interferogram.multilookL; // initialize multilookfactor
int32 mlP = interferogram.multilookP;
const int32 &ifgmlL = interferogram.multilookL; // multilookfactor of interferogram
const int32 &ifgmlP = interferogram.multilookP; // multilookfactor of interferogram
if ( ifgmlL != 1 || ifgmlP != 1 ) // [MA] additional entry for Coherence comptation using refdem.
{ // always do computation without multilooking
mlL = 1; // set multilookfactor for interpolation
mlP = 1; // set multilookfactor for interpolation
mlookedIFG = true; // dealing with mlooked ifg.
}
// _____ end added by MA _____
const real8 m_min4picdivlam = (-4.0*PI*SOL)/master.wavelength;
const real8 s_min4picdivlam = (-4.0*PI*SOL)/slave.wavelength;
DEBUG << "master wavelength = " << master.wavelength;
DEBUG.print();
DEBUG << "slave wavelength = " << slave.wavelength;
DEBUG.print();
const real8 latNfile = lat0file-DEMdeltalat*(numberoflatpixels-1); // upper=max. lat value
const real8 lonNfile = lon0file+DEMdeltalon*(numberoflonpixels-1); // left=min. lon value
// ______ Extra info ______
INFO << "DEM input: w/e/s/n: \t"
<< rad2deg(lon0file) << "/" << rad2deg(lonNfile) << "/"
<< rad2deg(latNfile) << "/" << rad2deg(lat0file);
INFO.print();
// ______ Get corners of interferogram (approx) to select DEM ______
// ______ in radians (if height were zero)______
real8 extralat = (1.5*DEMdeltalat + deg2rad(4.0/25.0));
real8 extralong = (1.5*DEMdeltalon + deg2rad(4.0/25.0));
real8 phimin;
real8 phimax;
real8 lambdamin;
real8 lambdamax;
int32 indexphi0DEM;
int32 indexphiNDEM;
int32 indexlambda0DEM;
int32 indexlambdaNDEM;
const uint &ifglinelo = interferogram.win.linelo; // [MA] win no-mlooked master coords
const uint &ifglinehi = interferogram.win.linehi;
const uint &ifgpixlo = interferogram.win.pixlo;
const uint &ifgpixhi = interferogram.win.pixhi;
getcorners(ifglinelo,ifglinehi,
ifgpixlo,ifgpixhi,
extralat,extralong,lat0file,lon0file,
DEMdeltalat,DEMdeltalon,numberoflatpixels,numberoflonpixels,
ellips,master,masterorbit,phimin,phimax,lambdamin,lambdamax,
indexphi0DEM,indexphiNDEM,indexlambda0DEM,indexlambdaNDEM);
// ______ Extra info ______
INFO << "DEM input required: w/e/s/n: \t"
<< rad2deg(lambdamin) << "/" << rad2deg(lambdamax) << "/"
<< rad2deg(phimin) << "/" << rad2deg(phimax);
INFO.print();
INFO << "For window (l0,lN,p0,pN): \t"
<< ifglinelo << " "
<< ifglinehi << " "
<< ifgpixlo << " "
<< ifgpixhi;
INFO.print();
// ______ Check corners of DEM ______
// check if DEM is appropriate for interferogram
// DEM should at least partially cover IFG
// note: phi is [90:-90]
if (phimax <= latNfile)// DEM is more north than IFG
{
ERROR << "IFG outside DEM: most South latitude: " << rad2deg(latNfile)
<< " [deg]; IFG requires: " << rad2deg(phimax)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
throw(some_error);
}
// DEM is more south than IFG
if (phimin >= lat0file)// largest latitude at first line of file
{
ERROR << "IFG outside DEM: most North latitude: " << rad2deg(lat0file)
<< " [deg]; IFG requires: " << rad2deg(phimax)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
throw(some_error);
}
if (lambdamax <= lon0file)
{
ERROR << "IFG outside DEM: most West longitude: " << rad2deg(lon0file)
<< " [deg]; IFG window requires: " << rad2deg(lambdamax)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
throw(some_error);
}
if (lambdamin >= lonNfile)
{
ERROR << "IFG outside DEM: most East longitude: " << rad2deg(lonNfile)
<< " [deg]; IFG window requires: " << rad2deg(lambdamin)
<< " [deg]";
PRINT_ERROR(ERROR.get_str())
throw(some_error);
}
//===================================================================
//============ First loop: radarcode DEM ============================
//============ (DEM geometry) ============================
//===================================================================
int32 numvalid = 0;// number of good values, not NODATA in buffer
int32 numNODATA = 0;// number of NODATA values in buffer
real8 meancroppedDEM = 0.0;// to detect byte order problems, formats
real8 min_input_dem = 100000.0;// stats
real8 max_input_dem = -100000.0;// stats
// ______ Compute buffer size radarcoding DEM______
const real8 BUFFERMEMSIZE = generalinput.memory;// Bytes
const int32 NcolsDEM = indexlambdaNDEM-indexlambda0DEM+1;
const int32 NrowsDEM = indexphiNDEM-indexphi0DEM+1;
const real8 Nrows_possible_DEM = BUFFERMEMSIZE / (5*8*NcolsDEM);
int32 bufferlines = int32(ceil(Nrows_possible_DEM));
INFO << "Possible max. buffer lines: " << bufferlines << " for " << BUFFERMEMSIZE << " memory size.";
INFO.print();
if (bufferlines>NrowsDEM) bufferlines=NrowsDEM;
int32 numfullbuffers = NrowsDEM / bufferlines;
int32 restlines = NrowsDEM % bufferlines;
int32 extrabuffer = (restlines == 0) ? 0 : 1;
// ______ Extra info ______
INFO << "DEM output total pixels: " << NcolsDEM;
INFO.print();
INFO << "DEM output total lines : " << NrowsDEM;
INFO.print();
INFO << "Radar coding of DEM in: " << numfullbuffers << " buffers of "
<< bufferlines << " lines and " << extrabuffer << " extra buffer of "
<< restlines << " lines.";
INFO.print();
// ______ Open (temporary) output files ______
// DEM heights
ofstream demofile;
openfstream(demofile,refdeminput.fodem,generalinput.overwrit);
bk_assert(demofile,refdeminput.fodem,__FILE__,__LINE__);
// master line coordinates of DEM
ofstream masterdemlineoutfile("crd_m_demline.temp", ios::out | ios::trunc);
bk_assert(masterdemlineoutfile,"crd_m_demline.temp",__FILE__,__LINE__);
// master pixel coordinates of DEM
ofstream masterdempixeloutfile("crd_m_dempixel.temp", ios::out | ios::trunc);
bk_assert(masterdempixeloutfile,"crd_m_dempixel.temp",__FILE__,__LINE__);
// ref phase in DEM geometry
ofstream demrefphaseoutfile("crd_dem_refphase.temp", ios::out | ios::trunc);
bk_assert(demrefphaseoutfile,"crd_dem_refphase.temp",__FILE__,__LINE__);
// h2ph factor in DEM geometry
ofstream demh2phoutfile;
if (outputh2ph==true)
{
openfstream(demh2phoutfile,"crd_dem_h2ph.temp",generalinput.overwrit);
bk_assert(demh2phoutfile,"crd_dem_h2ph.temp",__FILE__,__LINE__);
}
// ______ DEM loop per buffer ______
register int32 j,i;// DEM index grid counter, register j first to ensure allocation
for (register int32 buffer=0; buffer<numfullbuffers+extrabuffer; ++buffer)
{
// Determine indices for buffer
const int32 indexphi0BUFFER = indexphi0DEM+buffer*bufferlines;
const int32 blines = (buffer == numfullbuffers) ? restlines : bufferlines;
const int32 indexphiNBUFFER = indexphi0BUFFER+(blines-1);
matrix<real4> DEM(blines,NcolsDEM);
// ______ Extra info ______
PROGRESS << STEP << "Buffer# [l0:lN, p0:pN]: " << buffer+1 << " ["
<< indexphi0BUFFER << ": " << indexphiNBUFFER << ", "
<< indexlambda0DEM << ": " << indexlambdaNDEM << "]";
PROGRESS.print();
// ______ lat/lon for first pixel in matrix read from file ______
// ______ upper is max. latitude, left is min. longitude ______
const real8 upperleftphi = lat0file-indexphi0BUFFER*DEMdeltalat;
const real8 upperleftlambda = lon0file+indexlambda0DEM*DEMdeltalon;
window zerooffset (0,0,0,0);
window winfromfile (indexphi0BUFFER,indexphiNBUFFER,
indexlambda0DEM,indexlambdaNDEM);
// ______ Read in grdfile of DEM in matrix R4 (raw data, no header) _______
// ______ added formats (BK 4-May-2001) ______
PROGRESS << STEP << "Reading crop of DEM for buffer: " << buffer+1;
PROGRESS.print();
DEBUG.print("Reading input DEM into real4 matrix (buffer).");
switch (refdeminput.iformatflag)
{
// ______ Read as short BE, then convert to host order ______
case FORMATI2_BIGENDIAN:
{
matrix<int16> DEMi2(blines,NcolsDEM);
readfile(DEMi2,refdeminput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = real4(ntohs(DEMi2(iii,jjj)));// cast to real4
DEMi2.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: SHORT SIGNED INT.");
break;
}
case FORMATI2:
{
matrix<int16> DEMi2(blines,NcolsDEM);
readfile(DEMi2,refdeminput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = DEMi2(iii,jjj);// cast to real4
DEMi2.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: SHORT SIGNED INT.");
break;
}
case FORMATR4:
readfile(DEM,refdeminput.firefdem,numberoflatpixels,winfromfile,zerooffset);
INFO.print("Read crop of input DEM: format: REAL4.");
break;
case FORMATR8:
{
matrix<real8> DEMr8(blines,NcolsDEM);
readfile(DEMr8,refdeminput.firefdem,numberoflatpixels,winfromfile,zerooffset);
for (int32 iii=0; iii<DEM.lines(); ++iii)
for (int32 jjj=0; jjj<DEM.pixels(); ++jjj)
DEM(iii,jjj) = DEMr8(iii,jjj);// cast to real4
DEMr8.resize(1,1);// dealloc...
INFO.print("Read crop of input DEM: format: REAL8.");
break;
}
default:
PRINT_ERROR("totally impossible, checked input.")
//throw(unhandled_case_error);
}
// ----- Loop over DEM for stats ------------------------
real8 min_dem_buffer = 100000.0;
real8 max_dem_buffer = -100000.0;
for (i=0; i<DEM.lines(); ++i)
{
// ----- Loop over oversampled matrix in x ------
for (j=0; j<DEM.pixels(); ++j)
{
if(DEM(i,j)!=NODATA)
{
numvalid++;
meancroppedDEM += DEM(i,j);// divide by numvalid later
if (DEM(i,j)<min_dem_buffer) min_dem_buffer=DEM(i,j);//buffer
if (DEM(i,j)>max_dem_buffer) max_dem_buffer=DEM(i,j);// stats
}
else
{
numNODATA++;
}
}//loop dem for stats
}//loop dem for stats
min_input_dem = min(min_input_dem,min_dem_buffer);//global stats
max_input_dem = max(max_input_dem,max_dem_buffer);//global stats
// ====== Radarcoding DEM ==============================
// ______ DEM contains values from leftupper with ______
// ______ spacing (DEMdeltalat,DEMdeltalon) ______
// ______ Transform DEM to l,p,refphase ______
PROGRESS.print("Converting DEM to radar system for this buffer.");
const int32 NpointsDEM = DEM.size();
const int32 NpixelsDEM = DEM.pixels();
// ______ Extra info ______
INFO << "Number of points in DEM: "
<< NpointsDEM;
INFO.print();
matrix<real8> masterDEMline(DEM.lines(),DEM.pixels());
matrix<real8> masterDEMpixel(DEM.lines(),DEM.pixels());
matrix<real8> ref_phase_array(DEM.lines(),DEM.pixels());
matrix<real8> h2ph_array(DEM.lines(),DEM.pixels());
// --- Loop DEM ---
cn P;
real8 phi,lambda,height,l,p,ref_phase;
phi = upperleftphi;
for (i=0; i<DEM.lines(); ++i)
{
if ((i%100)==0)
{
// ______ Extra info ______
PROGRESS << STEP << "Radarcoding buffer: " << buffer+1 << "of" << numfullbuffers + extrabuffer << " DEM line: " << i << " ("
<< floor(.5+(100.*real8(i)/real8(DEM.lines())))
<< "%)";
PROGRESS.print();
}
lambda = upperleftlambda;
for (j=0; j<DEM.pixels(); ++j)
{
height = DEM(i,j);
ell2lp(l,p,ellips,master,masterorbit,phi,lambda,height,MAXITER,CRITERTIM);
masterDEMline(i,j) = l;
masterDEMpixel(i,j) = p;
P = ellips.ell2xyz(phi,lambda,height); // returns P(x,y,z)
real8 t_range_master;
real8 t_azi_master;
xyz2t(t_azi_master,t_range_master,
master,masterorbit,
P,MAXITER,CRITERTIM);
real8 t_azi_slave;
real8 t_range_slave;
xyz2t(t_azi_slave,t_range_slave, // t returned
slave,slaveorbit,
P,MAXITER,CRITERTIM);
if (outputh2ph==true)
{
// compute h2ph factor
cn Psat_master = masterorbit.getxyz(t_azi_master);
cn Psat_slave = slaveorbit.getxyz(t_azi_slave);
real8 B = Psat_master.dist(Psat_slave);
real8 Bpar = SOL*(t_range_master-t_range_slave);
cn r1 = Psat_master.min(P);
cn r2 = Psat_slave.min(P);
// real8 theta = Psat_master.angle(r1); // look angle
real8 theta = P.angle(r1); // incidence angle
real8 theta_slave = P.angle(r2); // incidence angle slave
real8 Bperp = (theta > theta_slave ) ? // sign ok
sqrt(sqr(B)-sqr(Bpar)) :
-sqrt(sqr(B)-sqr(Bpar)) ;
h2ph_array(i,j) = Bperp/(t_range_master*SOL*sin(theta));
}
if (onlyrefphasetopo) // do not include flat earth phase
{
lp2xyz(l,p,ellips, // h==0
master, masterorbit,
P,MAXITER,CRITERPOS); // P returned
real8 t_range_flatearth,t_azi_dummy;
xyz2t(t_azi_dummy,t_range_flatearth,
slave,slaveorbit,
P,MAXITER,CRITERTIM); // P on h=0
ref_phase = s_min4picdivlam*t_range_flatearth-
s_min4picdivlam*t_range_slave;
}
else // include flatearth, ref.pha = phi_topo+phi_flatearth
{
ref_phase =
m_min4picdivlam*master.pix2tr(p)-
s_min4picdivlam*t_range_slave;
}
ref_phase_array(i,j) = ref_phase;
lambda += DEMdeltalon;
} // loop DEM pixels
// ______ update latitude of next line ______
phi -= DEMdeltalat; // upper left is max. value
} // loop DEM lines
// Write results to output files
PROGRESS << STEP << "Writing radar coded DEM to file, buffer: " << buffer+1 << " of " << numfullbuffers + extrabuffer ;
PROGRESS.print();
demofile << DEM;
masterdemlineoutfile << masterDEMline;
masterdempixeloutfile << masterDEMpixel;
demrefphaseoutfile << ref_phase_array;
if (outputh2ph==true)
demh2phoutfile << h2ph_array;
masterDEMline.resize(1,1); //deallocate
masterDEMpixel.resize(1,1); //deallocate
DEM.resize(1,1); //deallocate
ref_phase_array(1,1); //deallocate
h2ph_array(1,1); //deallocate
} // buffer loop
demofile.close();
masterdemlineoutfile.close();
masterdempixeloutfile.close();
demrefphaseoutfile.close();
if (outputh2ph==true)
demh2phoutfile.close();
//===================================================================
//============ End first loop: radarcode DEM ========================
//============ (DEM geometry) ============================
//===================================================================
//===================================================================
//============ Second loop: interpolation =============
//============ (radar geometry) =============
//===================================================================
INFO << STEP << "Start interpolation...";
INFO.print();
// ______ Line/pixel of first point in original master coordinates ______
// ______ maybe this should be changed to be x+(ml/2) ?? but depends on
// ______ definition of range_to_first_bin is to center or start..
// Bert Kampes, 08-Apr-2005: chose center by adding ml/2
const int32 Nlinesml = interferogram.win.lines() / mlL; // ifg lines when mlL = 1 (no multilooking)
const int32 Npixelsml = interferogram.win.pixels() / mlP;
const int32 ifgNlinesml = interferogram.win.lines() / ifgmlL; // for the result file when mlL != 1
const int32 ifgNpixelsml = interferogram.win.pixels() / ifgmlP;
const real8 offset = 0;
//cerr << "xNFO: linesnoml: " << Nlinesml << " pixnoml: " << Npixelsml << endl;
//cerr << "xNFO: ifglinesml: " << ifgNlinesml << " ifgpixml: " << ifgNpixelsml << endl;
const real8 veryfirstline = real8(ifglinelo) + (real8(mlL)-1.0)/2.0;
const real8 verylastline = veryfirstline + real8((Nlinesml-1)*mlL);
const real8 firstpixel = real8(ifgpixlo) + (real8(mlP)-1.0)/2.0;
const real8 lastpixel = firstpixel + real8((Npixelsml-1)*mlP);
//cerr << "xNFO: vl0:vlN " << veryfirstline << ":" << verylastline << " p1:pN " << firstpixel << ":" << lastpixel << endl;
//Determine range-azimuth spacing ratio, needed for proper triangulation
cn P1, P2 , P3, P4;
lp2xyz(veryfirstline,firstpixel,ellips,master,masterorbit,
P1,MAXITER,CRITERPOS);
lp2xyz(veryfirstline,lastpixel,ellips,master,masterorbit,
P2,MAXITER,CRITERPOS);
lp2xyz(verylastline,firstpixel,ellips,master,masterorbit,
P3,MAXITER,CRITERPOS);
lp2xyz(verylastline,lastpixel,ellips,master,masterorbit,
P4,MAXITER,CRITERPOS);
const real8 r_spacing = ( (P1.min(P2)).norm() + (P3.min(P4)).norm() ) / 2 /(lastpixel - firstpixel) ;
const real8 az_spacing = ( (P1.min(P3)).norm() + (P2.min(P4)).norm() ) /2 /(verylastline - veryfirstline);
const real8 r_az_ratio = r_spacing/az_spacing;
INFO << "Interferogram azimuth spacing: " << az_spacing;
INFO.print();
INFO << "Interferogram range spacing: " << r_spacing;
INFO.print();
INFO << "Range-azimuth spacing ratio: " << r_az_ratio;
INFO.print();
// ______ Compute buffer size interpolation______
const real8 Nlinesml_possible = BUFFERMEMSIZE / (6*8*Npixelsml);
bufferlines = int32(ceil(Nlinesml_possible)); // initialized
if ( mlookedIFG == true) // if ifg is multilooked by a factor
{
bufferlines = int32( floor(Nlinesml_possible/ifgmlL) * ifgmlL ); // [HB] Hermann provided the fix:
// Successive bufferlines must have a size which is multiple of multilooking factor
// unless data can fit completely to an initial single buffer.
// Extra buffer will scale correctly and rounding due to multilooking
// will yield correct number lines for the output file.
// Ex: floor(2097/25)*2 + floor(1578/25) = 229 lines (wrong) (bufferlines/mlL)*Nbuffers+extrabufferlines/mlL
} // Ex: floor(2075/25)*25*2 + floor(1622/25) = 230 lines (correct)
else // no-multilooking
{
bufferlines = int32( floor(Nlinesml_possible) ); // [MA] instead of ceil, prefered floor to use less mem
}
if (bufferlines > Nlinesml) bufferlines=Nlinesml; // if bufferlines > Nlines then shrink bufferlines to Nlines, no extra buffer requested.
INFO << "Possible max. buffer lines: " << bufferlines << " for " << BUFFERMEMSIZE << " memory size.";
INFO.print();
numfullbuffers = Nlinesml / bufferlines;
restlines = Nlinesml % bufferlines; // the number of lines in extra buffer
extrabuffer = (restlines == 0) ? 0 : 1;
// ______ Extra info ______
INFO << "Interpolation in: " << numfullbuffers << " buffers of "
<< bufferlines << " lines and " << extrabuffer << " extra buffer of "
<< restlines << " lines.";
INFO.print();
// ______ Open output files ______
ofstream refdemofile; // refdem phase
openfstream(refdemofile,refdeminput.forefdem,generalinput.overwrit);
bk_assert(refdemofile,refdeminput.forefdem,__FILE__,__LINE__);
// _____ start added by MA _____
ofstream refdemofilenoML; // [MA] refdem phase no-multilooked
{ // local scope practice
string fname = string(refdeminput.forefdem) + ".noML"; // new name as m_s_refdemphase.raw.noML
if ( mlookedIFG == true) // if ifg is multilooked by a factor
{
openfstream(refdemofilenoML,fname.c_str(),generalinput.overwrit);
bk_assert(refdemofilenoML,fname.c_str(),__FILE__,__LINE__);
}
else // no-multilooking
{
if(!remove(fname.c_str())) // when success report removed.
{
WARNING << "Removed existing " << fname << "file.";
WARNING.print();
}
}
}
// _____ end added by MA _____
// if request for height in radar coordinates l,p
ofstream refdemheiofile;// Radarcoded DEM (Z.Perski)
if (outputrefdemhei==true)
{
openfstream(refdemheiofile,refdeminput.forefdemhei,generalinput.overwrit);
bk_assert(refdemheiofile,refdeminput.forefdemhei,__FILE__,__LINE__);
}
// if request for h2ph in radar coordinates l,p
ofstream h2phofile;
if (outputh2ph==true)
{
openfstream(h2phofile,refdeminput.foh2ph,generalinput.overwrit);
bk_assert(h2phofile,refdeminput.foh2ph,__FILE__,__LINE__);
}
// ______ interpolation loop per buffer ______
for (register int32 buffer = 0; buffer < numfullbuffers + extrabuffer; ++buffer)
{
// Determine indices for buffer
const int32 blines = (buffer == numfullbuffers) ? restlines : bufferlines;
const real8 firstline_buffer = veryfirstline+buffer*bufferlines*mlL;
const real8 lastline_buffer = firstline_buffer+(blines-1)*mlL;
// ______ Extra info ______
PROGRESS << STEP << "Interpolation buffer: " << buffer+1 << "of" << numfullbuffers + extrabuffer << " [l0:lN, p0:pN]: " << " ["
<< firstline_buffer << ": " << lastline_buffer << ", "
<< firstpixel << ": " << lastpixel << "]";
PROGRESS.print();
// Get corners of buffer
real8 phimin_az;
real8 phimax_az;
real8 lambdamin_az;
real8 lambdamax_az;
getcorners(firstline_buffer+offset,lastline_buffer+offset,
firstpixel+offset,lastpixel+offset,
extralat,extralong,phimax,lambdamin,
DEMdeltalat,DEMdeltalon,NrowsDEM,NcolsDEM,
ellips,master,masterorbit,phimin_az,phimax_az,lambdamin_az,lambdamax,
indexphi0DEM,indexphiNDEM,indexlambda0DEM,indexlambdaNDEM);
window zerooffset (0,0,0,0);
window winfromfile (indexphi0DEM,indexphiNDEM,
indexlambda0DEM,indexlambdaNDEM);
const int32 NrowsDEM_buffer = indexphiNDEM-indexphi0DEM+1;
const int32 NcolsDEM_buffer = indexlambdaNDEM-indexlambda0DEM+1;
PROGRESS << STEP << "Reading input for interpolation buffer: " << buffer+1 << "of" << numfullbuffers + extrabuffer;
PROGRESS.print();
// read x,y
matrix<real8> DEMline_buffer(NrowsDEM_buffer,NcolsDEM_buffer);
matrix<real8> DEMpixel_buffer(NrowsDEM_buffer,NcolsDEM_buffer);
readfile(DEMline_buffer,"crd_m_demline.temp",NrowsDEM,winfromfile,zerooffset);
readfile(DEMpixel_buffer,"crd_m_dempixel.temp",NrowsDEM,winfromfile,zerooffset);
// read z (multiple, number can easily be increased, e.g. simulated intensity)
int32 Nz = 1; //number of z
matrix<real8> input_buffer(NrowsDEM_buffer *Nz ,NcolsDEM_buffer);
matrix<real8> temp_input_buffer(NrowsDEM_buffer,NcolsDEM_buffer);
if (outputrefdemhei==true)
{
Nz += 1;
input_buffer.resize(NrowsDEM_buffer *Nz ,NcolsDEM_buffer);
}
if (outputh2ph==true)
{
Nz += 1;
input_buffer.resize(NrowsDEM_buffer *Nz ,NcolsDEM_buffer);
}
readfile(temp_input_buffer,"crd_dem_refphase.temp",NrowsDEM,winfromfile,zerooffset);
input_buffer.setdata(0, 0, temp_input_buffer);
Nz = 1;
if (outputrefdemhei==true)
{
Nz += 1;
/// i would like to use real4, test later on
matrix<real4> dem_input(NrowsDEM_buffer,NcolsDEM_buffer);
readfile(dem_input,refdeminput.fodem,NrowsDEM,winfromfile,zerooffset);
for (register int32 i =0 ; i < NrowsDEM_buffer ; i ++)
for(register int32 j = 0; j < NcolsDEM_buffer; j++)
temp_input_buffer(i,j) = real8(dem_input(i,j));
input_buffer.setdata(NrowsDEM_buffer * (Nz-1), 0, temp_input_buffer);
}
if (outputh2ph==true)
{
Nz += 1;
readfile(temp_input_buffer,"crd_dem_h2ph.temp",NrowsDEM,winfromfile,zerooffset);
input_buffer.setdata(NrowsDEM_buffer * (Nz-1) , 0, temp_input_buffer);
}
// initialize output array
Nz = 1;
matrix<real8> output_buffer(blines * Nz, Npixelsml);
if (outputrefdemhei==true)
{
Nz += 1;
output_buffer.resize(blines * Nz, Npixelsml);
}
if (outputh2ph==true)
{
Nz += 1;
output_buffer.resize(blines * Nz, Npixelsml);
}
// interpolation
griddatalinear(DEMline_buffer,DEMpixel_buffer,input_buffer,
firstline_buffer,lastline_buffer,firstpixel,lastpixel,
mlL,mlP,r_az_ratio,offset,NODATA,output_buffer);
//MA multilooking will start here.
//MA cast all output files to type real4
matrix<real4> output_layer(blines,Npixelsml);
Nz = 1; // matrix 3rd dimension counter, 1 --> PHASE
for (register int32 i =0 ; i < blines ; i++)
for(register int32 j = 0; j < Npixelsml; j++)
output_layer(i,j) = real4(output_buffer(i,j)); // real8 --> real4
// convert_type(output_buffer,output_layer); // TODO MA, should replace above 3 lines but this one doesn't work properly yet
//cerr << "refphase: blines: " << blines << " Npixelsml " << Npixelsml << endl;
// _____ start added by MA _____
if ( mlookedIFG == true) // [MA] if ifg is multilooked by a factor
{
refdemofile << multilook(output_layer, ifgmlL, ifgmlP); // multilook to output
refdemofilenoML << output_layer; // default output: PHASE
}
else
{
refdemofile << output_layer;
//refdemofile << output_buffer(window(0, blines - 1, 0, Npixelsml -1 ));
}
// _____ end added by MA _____
if (outputrefdemhei==true)
{
Nz += 1;
for (register int32 i = blines * (Nz-1) ; i < blines * Nz ; i++)
for(register int32 j = 0; j < Npixelsml; j++)
output_layer(i-(blines * (Nz-1)),j) = real4(output_buffer(i,j)); // real8 --> real4
//refdemheiofile << output_layer; // output reference dem heights
(mlookedIFG == true) ? refdemheiofile << multilook(output_layer, ifgmlL, ifgmlP) // [MA]
: refdemheiofile << output_layer ;
//refdemheiofile << output_buffer(window((Nz-1) * blines,Nz * blines - 1, 0, Npixelsml -1 ));
}
if (outputh2ph==true)
{
Nz += 1;
for (register int32 i = blines * (Nz-1) ; i < blines * Nz ; i++)
for(register int32 j = 0; j < Npixelsml; j++)
output_layer(i-(blines * (Nz-1)),j) = real4(output_buffer(i,j)); // real8 --> real4
//h2phofile << output_layer; // output h2ph matrix
(mlookedIFG == true) ? h2phofile << multilook(output_layer, ifgmlL, ifgmlP) // [MA]
: h2phofile << output_layer;
//h2phofile << output_buffer(window((Nz-1) * blines,Nz * blines - 1, 0, Npixelsml -1 ));
}
DEMline_buffer.resize(1,1); // deallocate
DEMpixel_buffer.resize(1,1);
input_buffer.resize(1,1);
temp_input_buffer.resize(1,1);
output_buffer.resize(1,1);
} // end loop azimuth direction
INFO << "Closing output files";
INFO.print();
refdemofile.close(); // has the same multilook as interferrogram
if (mlookedIFG==true)
refdemofilenoML.close(); // [MA] if interferogram is mlooked then this is
// generated as non-multilooked for coherence estimation.
if (outputrefdemhei==true) // For Zbigniew Perski
refdemheiofile.close();
if (outputh2ph==true)
h2phofile.close();
//===================================================================
//============ End second loop: interpolation =============
//============ (radar geometry) =============
//===================================================================
// === Clean up temporary outfiles === [MA]
if (remove("crd_m_demline.temp")) // remove files
WARNING.print("code 101: could not remove file: crd_m_demline.temp.");
if (remove("crd_m_dempixel.temp"))
WARNING.print("code 101: could not remove file: crd_m_dempixel.temp.");
if (remove("crd_dem_refphase.temp"))
WARNING.print("code 101: could not remove file: crd_dem_refphase.temp.");
if ( outputh2ph==true && remove("crd_dem_h2ph.temp")) // output available
WARNING.print("code 101: could not remove file: crd_dem_h2ph.temp.");
// === Clean up Done. ===
// ====== Write output information ======
char croppeddemi[4*ONE27];
strcpy(croppeddemi,"NO output requested");
if (outputdemi) strcpy(croppeddemi,refdeminput.fodemi);
INFO << "Min. value of input DEM covering interferogram: " << min_input_dem;
INFO.print();
INFO << "Max. value of input DEM covering interferogram: " << max_input_dem;
INFO.print();
ofstream scratchlogfile("scratchlogcomprefdem", ios::out | ios::trunc);
bk_assert(scratchlogfile,"comprefdem: scratchlogcomprefdem",__FILE__,__LINE__);
scratchlogfile
<< "\n*******************************************************************"
<< "\n* " << processcontrol[pr_i_comprefdem]
<< "\n*******************************************************************"
<< "\n1) DEM source file: \t" << refdeminput.firefdem
<< "\nFormat: \t";
switch (refdeminput.iformatflag)
{
case FORMATI2:
{
scratchlogfile << "SHORT SIGNED INTEGER (HOST ENDIANNESS)";
break;
}
case FORMATI2_BIGENDIAN:
{
scratchlogfile << "SHORT SIGNED INTEGER, BIG ENDIAN";
break;
}
case FORMATR4:
{
scratchlogfile << "REAL4 SIGNED FLOAT";
break;
}
case FORMATR8:
{
scratchlogfile << "REAL8 SIGNED DOUBLE";
break;
}
default:
{
scratchlogfile << "UNKNOWN? IMPOSSIBLE...";
break;
}
}
scratchlogfile
<< "\nByte order: \t" << "check it yourself..."
<< "\nNumber of lines: \t" << numberoflatpixels
<< "\nNumber of pixels: \t" << numberoflonpixels
<< "\nResolution latitude: \t" << rad2deg(DEMdeltalat) << " [deg]"
<< "\nResolution longitude: \t" << rad2deg(DEMdeltalon) << " [deg]"
<< "\nMost West point in input DEM: \t" << rad2deg(lon0file)
<< "\nMost East point in input DEM: \t" << rad2deg(lonNfile)
<< "\nMost South point in input DEM: \t" << rad2deg(latNfile)
<< "\nMost North point in input DEM: \t" << rad2deg(lat0file)
<< "\nMin. value of input DEM covering interferogram: " << min_input_dem
<< "\nMax. value of input DEM covering interferogram: " << max_input_dem
<< "\n2) Output file cropped DEM: \t" << refdeminput.fodem //[FVL
<< "\nFormat: \t" << "REAL4"
<< "\nByte order: \t" << "(same as host)"
<< "\nNumber of lines (multilooked): \t" << NcolsDEM
<< "\nNumber of pixels (multilooked): \t" << NrowsDEM
<< "\nDEM extend w/e/s/n : \t" << rad2deg(lambdamin) << "/"
<< rad2deg(lambdamax) << "/" << rad2deg(phimin) << "/" << rad2deg(phimax)
// << "\nMean value: \t" << meancroppedDEM
<< "\n3) Output file interpolated crop DEM: \t" << croppeddemi
<< "\nFormat: \t" << "REAL4"
<< "\nByte order: \t" << "(same as host)"
<< "\n4) Output file synthetic phase: \t" << refdeminput.forefdem
<< "\nFormat: \t" << "REAL4"
<< "\nByte order: \t" << "(same as host)"
<< "\nNumber of lines (multilooked): \t" << ifgNlinesml
<< "\nNumber of pixels (multilooked): \t" << ifgNpixelsml
// this is not correct, only stats per buffer...
// << "\n\n----- Other STATS -----"
// << "\nTotal points in cropped DEM: \t" << numpoints
// << "\nNumber of valid points in DEM: \t" << numvalid
// << " (" << 100*numvalid/numpoints << "%)"
// << "\nNumber of NODATA points in DEM: \t" << numNODATA
// << " (" << 100*numNODATA/numpoints << "%)"
// << "\nMean height in meters at valid points:\t" << meancroppedDEM
<< "\n*******************************************************************\n\n";
scratchlogfile.close();
ofstream scratchresfile("scratchrescomprefdem", ios::out | ios::trunc);
bk_assert(scratchresfile,"comprefdem: scratchrescomprefdem",__FILE__,__LINE__);
scratchresfile
<< "\n\n*******************************************************************"
<< "\n*_Start_" << processcontrol[pr_i_comprefdem]
<< "\n*******************************************************************";
if (onlyrefphasetopo==true) scratchresfile // [don]
<< "\nInclude_flatearth: \tNo";
else scratchresfile
<< "\nInclude_flatearth: \tYes";
scratchresfile
<< "\nDEM source file: \t" << refdeminput.firefdem
<< "\nMin. of input DEM: \t" << min_input_dem
<< "\nMax. of input DEM: \t" << max_input_dem
<< "\nData_output_file: \t" << refdeminput.forefdem
<< "\nData_output_format: \t" << "real4"
<< "\nFirst_line (w.r.t. original_master): \t"
<< interferogram.win.linelo
<< "\nLast_line (w.r.t. original_master): \t"
<< interferogram.win.linehi
<< "\nFirst_pixel (w.r.t. original_master): \t"
<< interferogram.win.pixlo
<< "\nLast_pixel (w.r.t. original_master): \t"
<< interferogram.win.pixhi
<< "\nMultilookfactor_azimuth_direction: \t" << ifgmlL
<< "\nMultilookfactor_range_direction: \t" << ifgmlP
<< "\nNumber of lines (multilooked): \t" << ifgNlinesml
<< "\nNumber of pixels (multilooked): \t" << ifgNpixelsml
<< "\n*******************************************************************"
<< "\n* End_" << processcontrol[pr_i_comprefdem] << "_NORMAL"
<< "\n*******************************************************************\n";
scratchresfile.close();
} // END radarcodedem
/****************************************************************
* getcorners *
* *
* Get corners of window (approx) to select DEM in radians (if *
* height were zero) *
* *
* Implementation: *
* 1) calculate phi, lambda of corners *
* 2) select the extreme values *
* 3) add extra overlap *
* 4) determine the indices in the file *
* *
* Freek van Leijen, 07-AUG-2006 *
* *
****************************************************************/
void getcorners(
const real8 &l0, // master.currentwindow.linelo
const real8 &lN, // master.currentwindow.
const real8 &p0, // master.currentwindow.
const real8 &pN, // master.currentwindow.
const real8 &extralat, //
const real8 &extralong, //
const real8 &lat0, // lat0file
const real8 &long0, // lon0file
const real8 &DEMdeltalat, //
const real8 &DEMdeltalong, //
const int32 &Nlatpixels, // numberoflatpixels
const int32 &Nlongpixels, // numberoflonpixels
const input_ell &ellips,
const slcimage &master,
orbit &masterorbit,
real8 &phimin,
real8 &phimax,
real8 &lambdamin,
real8 &lambdamax,
int32 &indexphi0DEM, // returned
int32 &indexphiNDEM, // returned
int32 &indexlambda0DEM, // returned
int32 &indexlambdaNDEM) // returned
{
TRACE_FUNCTION("getcorners (FvL 07-AUG-2006)")
DEBUG << "(getcorners) l0 :" << l0;
DEBUG.print();
DEBUG << "(getcorners) lN :" << lN;
DEBUG.print();
DEBUG << "(getcorners) p0 :" << p0;
DEBUG.print();
DEBUG << "(getcorners) pN :" << pN;
DEBUG.print();
DEBUG << "(getcorners) extralat :" << extralat;
DEBUG.print();
DEBUG << "(getcorners) extralong :" << extralong;
DEBUG.print();
DEBUG << "(getcorners) lat0 :" << lat0;
DEBUG.print();
DEBUG << "(getcorners) long0 :" << long0;
DEBUG.print();
DEBUG << "(getcorners) DEMdeltalat :" << DEMdeltalat;
DEBUG.print();
DEBUG << "(getcorners) DEMdeltalong :" << DEMdeltalong;
DEBUG.print();
DEBUG << "(getcorners) Total Nlatpixels :" << Nlatpixels;
DEBUG.print();
DEBUG << "(getcorners) Total Nlongpixels :" << Nlongpixels;
DEBUG.print();
real8 phi;
real8 lambda;
real8 height;
lp2ell(l0,p0,
ellips, master, masterorbit,
phi, lambda, height); // returned
real8 phil0p0 = phi;
real8 lambdal0p0 = lambda;
lp2ell(lN,p0,
ellips, master, masterorbit,
phi, lambda, height); // returned
real8 philNp0 = phi;
real8 lambdalNp0 = lambda;
lp2ell(lN,pN,
ellips, master, masterorbit,
phi, lambda, height); // returned
real8 philNpN = phi;
real8 lambdalNpN = lambda;
lp2ell(l0,pN,
ellips, master, masterorbit,
phi, lambda, height); // returned
real8 phil0pN = phi;
real8 lambdal0pN = lambda;
// ______ Select DEM values based on rectangle outside l,p border ______
phimin = min(min(min(phil0p0,philNp0),philNpN),phil0pN);
phimax = max(max(max(phil0p0,philNp0),philNpN),phil0pN);
lambdamin = min(min(min(lambdal0p0,lambdalNp0),lambdalNpN),lambdal0pN);
lambdamax = max(max(max(lambdal0p0,lambdalNp0),lambdalNpN),lambdal0pN);
// ______ a little bit extra at edges to be sure ______
// TODO this should change based on ascending or descending [MA]
phimin -= extralat;
phimax += extralat;
lambdamax += extralong;
lambdamin -= extralong;
DEBUG << "(getcorners) phimin :" << phimin;
DEBUG.print();
DEBUG << "(getcorners) phimax :" << phimax;
DEBUG.print();
DEBUG << "(getcorners) lambdamin :" << lambdamin;
DEBUG.print();
DEBUG << "(getcorners) lambdamax :" << lambdamax;
DEBUG.print();
// ______ Get indices of DEM needed ______
// ______ Index boundary: [0:numberofx-1] ______
indexphi0DEM = int32(floor((lat0-phimax)/DEMdeltalat));
if (indexphi0DEM < 0)
{
WARNING << "(getcorners) indexphi0DEM: " << indexphi0DEM; WARNING.print();
indexphi0DEM=0; // default start at first
WARNING.print("DEM does not cover entire interferogram.");
WARNING.print("input DEM should be extended to the North.");
}
indexphiNDEM = int32(ceil((lat0-phimin)/DEMdeltalat));
if (indexphiNDEM > Nlatpixels-1)
{
WARNING << "(getcorners) indexphiNDEM: " << indexphiNDEM; WARNING.print();
indexphiNDEM=Nlatpixels-1;
WARNING.print("DEM does not cover entire interferogram.");
WARNING.print("input DEM should be extended to the South.");
}
indexlambda0DEM = int32(floor((lambdamin-long0)/DEMdeltalong));
if (indexlambda0DEM < 0)
{
WARNING << "(getcorners) indexlambda0DEM: " << indexlambda0DEM; WARNING.print();
indexlambda0DEM=0; // default start at first
WARNING.print("DEM does not cover entire interferogram.");
WARNING.print("input DEM should be extended to the West.");
}
indexlambdaNDEM = int32(ceil((lambdamax-long0)/DEMdeltalong));
if (indexlambdaNDEM > Nlongpixels-1)
{
WARNING << "(getcorners) indexlambdaNDEM: " << indexlambdaNDEM; WARNING.print();
indexlambdaNDEM=Nlongpixels-1;
WARNING.print("DEM does not cover entire interferogram.");
WARNING.print("input DEM should be extended to the East.");
}
DEBUG << "(getcorners) indexphi0DEM :" << indexphi0DEM;
DEBUG.print();
DEBUG << "(getcorners) indexphiNDEM :" << indexphiNDEM;
DEBUG.print();
DEBUG << "(getcorners) indexlambda0DEM :" << indexlambda0DEM;
DEBUG.print();
DEBUG << "(getcorners) indexlambdaNDEM :" << indexlambdaNDEM;
DEBUG.print();
} // END getcorners
/****************************************************************
* griddatalinear (naming after Matlab function) *
* *
* Implementation after GMT function triangulate.c *
****************************************************************/
void griddatalinear(
const matrix<real8> &x_in,
const matrix<real8> &y_in,
const matrix<real8> &z_in,
const real8 &x_min,
const real8 &x_max,
const real8 &y_min,
const real8 &y_max,
const int32 &x_inc,
const int32 &y_inc,
const real8 &r_az_ratio,
const real8 &offset,
const real8 &NODATA,
matrix<real8> &grd
)
{
TRACE_FUNCTION("griddatalinear (LG&FvL 13-AUG-2006)")
INFO << "griddataLinear interpolation.";
INFO.print();
int32 i, j, k, ij, p, i_min, i_max, j_min, j_max;
int32 n, nx, ny, zLoops, zLoop, zBlockSize,indexFirstPoint,zInterpolateBlockSize;
real8 vx[4], vy[4], xkj, xlj, ykj, ylj, zj, zk, zl, zlj, zkj, xp, yp;
real8 f, *a = NULL , *b = NULL , *c = NULL ; // linear interpolation parameters
struct triangulateio In, Out, vorOut;
// Initialize variables
n = zBlockSize = x_in.size(); // block size of x and y coordination
/* How many groups of z value should be interpolated */
if (( z_in.size() % zBlockSize ) != 0)
{
INFO << "The input of the DEM buffer and z is not the same...";
INFO.print();
return;
}
else
zLoops = z_in.size()/x_in.size();
a = new real8[zLoops];
b = new real8[zLoops];
c = new real8[zLoops];
if (a == NULL || b == NULL || c == NULL)
{
ERROR << "Memory ERROR in source file: " << __FILE__
<< " at line: " << __LINE__;
PRINT_ERROR(ERROR.get_str());
throw(memory_error);
}
nx = grd.lines()/zLoops;
ny = grd.pixels();
zInterpolateBlockSize = grd.size()/zLoops;
/* Set everything to 0 and NULL */
memset ((void *)&In, 0, sizeof (struct triangulateio));
memset ((void *)&Out, 0, sizeof (struct triangulateio));
memset ((void *)&vorOut, 0, sizeof (struct triangulateio));
/* Allocate memory for input points */
In.numberofpoints = n ;
In.pointlist = new real8 [2 * n];
/* Copy x,y points to In structure array */
for (i = j = 0; i < n; i++)
{
In.pointlist[j++] = *(x_in[0] + i);
In.pointlist[j++] = (*(y_in[0] + i)) * r_az_ratio ;
// to eliminate the effect of difference in range and azimuth spacing;
}
/* Call Jonathan Shewchuk's triangulate algorithm. This is 64-bit safe since
* all the structures use 4-byte ints (longs are used internally). */
triangulate ("zIQB", &In, &Out, &vorOut);
//triangulate ("zIBV", &In, &Out, &vorOut); // [MA] for verbosing
int32 *link = Out.trianglelist; /* List of node numbers to return via link */
int32 np = Out.numberoftriangles;
for (k = ij = 0; k < np; k++)
{
DEBUG << "k of np, ij: " << k << " of " << np << ", :" << ij;
DEBUG.print();
//Store the Index of the first Point of this triangle.
indexFirstPoint = ij;
vx[0] = vx[3] = *(x_in[0] + link[ij]); vy[0] = vy[3] = *(y_in[0] + link[ij]); ij++;
vx[1] = *(x_in[0] + link[ij]); vy[1] = *(y_in[0]+link[ij]); ij++;
vx[2] = *(x_in[0] + link[ij]); vy[2] = *(y_in[0]+link[ij]); ij++;
if ( vx[0] == NODATA || vx[1] == NODATA || vx[2] == NODATA ) continue;
if ( vy[0] == NODATA || vy[1] == NODATA || vy[2] == NODATA ) continue;
/* Compute grid indices the current triangle may cover.*/
xp = min (min (vx[0], vx[1]), vx[2]); i_min = x_to_i (xp, x_min, x_inc, offset, nx);
//INFO << "xp: " << xp;
//INFO.print();
xp = max (max (vx[0], vx[1]), vx[2]); i_max = x_to_i (xp, x_min, x_inc, offset, nx);
//INFO << "xp: " << xp;
//INFO.print();
yp = min (min (vy[0], vy[1]), vy[2]); j_min = y_to_j (yp, y_min, y_inc, offset, ny);
//INFO << "yp: " << yp;
//INFO.print();
yp = max (max (vy[0], vy[1]), vy[2]); j_max = y_to_j (yp, y_min, y_inc, offset, ny);
//INFO << "yp: " << yp;
//INFO.print();
/* Adjustments for triangles outside -R region. */
/* Triangle to the left or right. */
if ((i_max < 0) || (i_min >= nx)) continue;
/* Triangle Above or below */
if ((j_max < 0) || (j_min >= ny)) continue;
/* Triangle covers boundary, left or right. */
if (i_min < 0) i_min = 0; if (i_max >= nx) i_max = nx - 1;
/* Triangle covers boundary, top or bottom. */
if (j_min < 0) j_min = 0; if (j_max >= ny) j_max = ny - 1;
// for (kk = 0; kk<npar;kk++) { //do for each parameter LIUG
// read zj, zk, zl (instead of above) LIUG
/* Find equation for the plane as z = ax + by + c */
xkj = vx[1] - vx[0]; ykj = vy[1] - vy[0];
xlj = vx[2] - vx[0]; ylj = vy[2] - vy[0];
f = 1.0 / (xkj * ylj - ykj * xlj);
for( zLoop = 0 ; zLoop < zLoops; zLoop++ )
{
zj = *(z_in[0] + zLoop * zBlockSize + link[indexFirstPoint]);
zk = *(z_in[0] + zLoop * zBlockSize + link[indexFirstPoint + 1]);
zl = *(z_in[0] + zLoop * zBlockSize + link[indexFirstPoint + 2]);
zkj = zk - zj; zlj = zl - zj;
a[zLoop] = -f * (ykj * zlj - zkj * ylj);
b[zLoop] = -f * (zkj * xlj - xkj * zlj);
c[zLoop] = -a[zLoop] * vx[1] - b[zLoop] * vy[1] + zk;
}
for (i = i_min; i <= i_max; i++)
{
xp = i_to_x (i, x_min, x_max, x_inc, offset, nx);
p = i * ny + j_min;
for (j = j_min; j <= j_max; j++, p++)
{
yp = j_to_y (j, y_min, y_max, y_inc, offset, ny);
if (!pointintriangle(vx, vy, xp, yp)) continue; /* Outside */
for( zLoop = 0 ; zLoop < zLoops; zLoop++ )
*(grd[0] + zLoop * zInterpolateBlockSize + p) = a[zLoop] * xp + b[zLoop] * yp + c[zLoop] ;
} //LIUG
}
}
if(a) delete[] a;
if(b) delete[] b;
if(c) delete[] c;
if(In.pointlist) delete[] In.pointlist; // only this use delete
if(Out.pointlist) free(Out.pointlist);
if(Out.trianglelist) free(Out.trianglelist);
} //END griddatalinear
/****************************************************************
* pointintriangle *
* *
* Liu, Guang and Freek van Leijen, 16-AUG-2006 *
* *
* *
****************************************************************/
int pointintriangle(double *xt,double * yt,double x ,double y)
{
int iRet0 = ((xt[2] - xt[0]) * (y - yt[0])) > ((x - xt[0]) * ( yt[2] - yt[0])) ? 1:-1;
int iRet1 = ((xt[0] - xt[1]) * (y - yt[1])) > ((x - xt[1]) * ( yt[0] - yt[1])) ? 1:-1;
int iRet2 = ((xt[1] - xt[2]) * (y - yt[2])) > ((x - xt[2]) * ( yt[1] - yt[2])) ? 1:-1;
if ((iRet0 >0 && iRet1 > 0 && iRet2 > 0 ) || (iRet0 <0 && iRet1 < 0 && iRet2 < 0 ))
return 1;
else
return 0;
} //END pointintriangle
|