41 #include "model/georef.h"
42 #include "model/cutil.h"
45 #define snprintf mysnprintf
49 static const long long lNaN = 0xfff8000000000000;
50 #define NAN (*(double *)&lNaN)
69 struct DATUM const gDatum[] = {
71 {
"Adindan", 5, -162, -12, 206},
72 {
"Afgooye", 15, -43, -163, 45},
73 {
"Ain el Abd 1970", 14, -150, -251, -2},
74 {
"Anna 1 Astro 1965", 2, -491, -22, 435},
75 {
"Arc 1950", 5, -143, -90, -294},
76 {
"Arc 1960", 5, -160, -8, -300},
77 {
"Ascension Island 58", 14, -207, 107, 52},
78 {
"Astro B4 Sorol Atoll", 14, 114, -116, -333},
79 {
"Astro Beacon E ", 14, 145, 75, -272},
80 {
"Astro DOS 71/4", 14, -320, 550, -494},
81 {
"Astronomic Stn 52", 14, 124, -234, -25},
82 {
"Australian Geod 66", 2, -133, -48, 148},
83 {
"Australian Geod 84", 2, -134, -48, 149},
84 {
"Bellevue (IGN)", 14, -127, -769, 472},
85 {
"Bermuda 1957", 4, -73, 213, 296},
86 {
"Bogota Observatory", 14, 307, 304, -318},
87 {
"Campo Inchauspe", 14, -148, 136, 90},
88 {
"Canton Astro 1966", 14, 298, -304, -375},
89 {
"Cape", 5, -136, -108, -292},
90 {
"Cape Canaveral", 4, -2, 150, 181},
91 {
"Carthage", 5, -263, 6, 431},
92 {
"CH-1903", 3, 674, 15, 405},
93 {
"Chatham 1971", 14, 175, -38, 113},
94 {
"Chua Astro", 14, -134, 229, -29},
95 {
"Corrego Alegre", 14, -206, 172, -6},
96 {
"Djakarta (Batavia)", 3, -377, 681, -50},
97 {
"DOS 1968", 14, 230, -199, -752},
98 {
"Easter Island 1967", 14, 211, 147, 111},
99 {
"European 1950", 14, -87, -98, -121},
100 {
"European 1979", 14, -86, -98, -119},
101 {
"Finland Hayford", 14, -78, -231, -97},
102 {
"Gandajika Base", 14, -133, -321, 50},
103 {
"Geodetic Datum 49", 14, 84, -22, 209},
104 {
"Guam 1963", 4, -100, -248, 259},
105 {
"GUX 1 Astro", 14, 252, -209, -751},
106 {
"Hermannskogel Datum", 3, 682, -203, 480},
107 {
"Hjorsey 1955", 14, -73, 46, -86},
108 {
"Hong Kong 1963", 14, -156, -271, -189},
109 {
"Indian Bangladesh", 6, 289, 734, 257},
110 {
"Indian Thailand", 6, 214, 836, 303},
111 {
"Ireland 1965", 1, 506, -122, 611},
112 {
"ISTS 073 Astro 69", 14, 208, -435, -229},
113 {
"Johnston Island", 14, 191, -77, -204},
114 {
"Kandawala", 6, -97, 787, 86},
115 {
"Kerguelen Island", 14, 145, -187, 103},
116 {
"Kertau 1948", 7, -11, 851, 5},
117 {
"L.C. 5 Astro", 4, 42, 124, 147},
118 {
"Liberia 1964", 5, -90, 40, 88},
119 {
"Luzon Mindanao", 4, -133, -79, -72},
120 {
"Luzon Philippines", 4, -133, -77, -51},
121 {
"Mahe 1971", 5, 41, -220, -134},
122 {
"Marco Astro", 14, -289, -124, 60},
123 {
"Massawa", 3, 639, 405, 60},
124 {
"Merchich", 5, 31, 146, 47},
125 {
"Midway Astro 1961", 14, 912, -58, 1227},
126 {
"Minna", 5, -92, -93, 122},
127 {
"NAD27 Alaska", 4, -5, 135, 172},
128 {
"NAD27 Bahamas", 4, -4, 154, 178},
129 {
"NAD27 Canada", 4, -10, 158, 187},
130 {
"NAD27 Canal Zone", 4, 0, 125, 201},
131 {
"NAD27 Caribbean", 4, -7, 152, 178},
132 {
"NAD27 Central", 4, 0, 125, 194},
133 {
"NAD27 CONUS", 4, -8, 160, 176},
134 {
"NAD27 Cuba", 4, -9, 152, 178},
135 {
"NAD27 Greenland", 4, 11, 114, 195},
136 {
"NAD27 Mexico", 4, -12, 130, 190},
137 {
"NAD27 San Salvador", 4, 1, 140, 165},
138 {
"NAD83", 11, 0, 0, 0},
139 {
"Nahrwn Masirah Ilnd", 5, -247, -148, 369},
140 {
"Nahrwn Saudi Arbia", 5, -231, -196, 482},
141 {
"Nahrwn United Arab", 5, -249, -156, 381},
142 {
"Naparima BWI", 14, -2, 374, 172},
143 {
"Observatorio 1966", 14, -425, -169, 81},
144 {
"Old Egyptian", 12, -130, 110, -13},
145 {
"Old Hawaiian", 4, 61, -285, -181},
146 {
"Oman", 5, -346, -1, 224},
147 {
"Ord Srvy Grt Britn", 0, 375, -111, 431},
148 {
"Pico De Las Nieves", 14, -307, -92, 127},
149 {
"Pitcairn Astro 1967", 14, 185, 165, 42},
150 {
"Prov So Amrican 56", 14, -288, 175, -376},
151 {
"Prov So Chilean 63", 14, 16, 196, 93},
152 {
"Puerto Rico", 4, 11, 72, -101},
153 {
"Qatar National", 14, -128, -283, 22},
154 {
"Qornoq", 14, 164, 138, -189},
155 {
"Reunion", 14, 94, -948, -1262},
156 {
"Rome 1940", 14, -225, -65, 9},
157 {
"RT 90", 3, 498, -36, 568},
158 {
"Santo (DOS)", 14, 170, 42, 84},
159 {
"Sao Braz", 14, -203, 141, 53},
160 {
"Sapper Hill 1943", 14, -355, 16, 74},
161 {
"Schwarzeck", 21, 616, 97, -251},
162 {
"South American 69", 16, -57, 1, -41},
163 {
"South Asia", 8, 7, -10, -26},
164 {
"Southeast Base", 14, -499, -249, 314},
165 {
"Southwest Base", 14, -104, 167, -38},
166 {
"Timbalai 1948", 6, -689, 691, -46},
167 {
"Tokyo", 3, -128, 481, 664},
168 {
"Tristan Astro 1968", 14, -632, 438, -609},
169 {
"Viti Levu 1916", 5, 51, 391, -36},
170 {
"Wake-Eniwetok 60", 13, 101, 52, -39},
171 {
"WGS 72", 19, 0, 0, 5},
172 {
"WGS 84", 20, 0, 0, 0},
173 {
"Zanderij", 14, -265, 120, -358},
174 {
"WGS_84", 20, 0, 0, 0},
175 {
"WGS-84", 20, 0, 0, 0},
176 {
"EUROPEAN DATUM 1950", 14, -87, -98, -121},
177 {
"European 1950 (Norway Finland)", 14, -87, -98, -121},
178 {
"ED50", 14, -87, -98, -121},
179 {
"RT90 (Sweden)", 3, 498, -36, 568},
180 {
"Monte Mario 1940", 14, -104, -49, 10},
181 {
"Ord Surv of Gr Britain 1936", 0, 375, -111, 431},
182 {
"South American 1969", 16, -57, 1, -41},
183 {
"PULKOVO 1942 (2)", 15, 25, -141, -79},
184 {
"EUROPEAN DATUM", 14, -87, -98, -121},
185 {
"BERMUDA DATUM 1957", 4, -73, 213, 296},
186 {
"COA", 14, -206, 172, -6},
187 {
"COABR", 14, -206, 172, -6},
188 {
"Roma 1940", 14, -225, -65, 9},
189 {
"ITALIENISCHES LANDESNETZ", 14, -87, -98, -121},
190 {
"HERMANSKOGEL DATUM", 3, 682, -203, 480},
191 {
"AGD66", 2, -128, -52, 153},
192 {
"ED", 14, -87, -98, -121},
193 {
"EUROPEAN 1950 (SPAIN AND PORTUGAL)", 14, -87, -98, -121},
194 {
"ED-50", 14, -87, -98, -121},
195 {
"EUROPEAN", 14, -87, -98, -121},
196 {
"POTSDAM", 14, -87, -98, -121},
197 {
"GRACIOSA SW BASE DATUM", 14, -104, 167, -38},
198 {
"WGS 1984", 20, 0, 0, 0},
199 {
"WGS 1972", 19, 0, 0, 5},
205 {
"Airy 1830", 6377563.396, 299.3249646},
206 {
"Modified Airy", 6377340.189, 299.3249646},
207 {
"Australian National", 6378160.0, 298.25},
208 {
"Bessel 1841", 6377397.155, 299.1528128},
209 {
"Clarke 1866", 6378206.4, 294.9786982},
210 {
"Clarke 1880", 6378249.145, 293.465},
211 {
"Everest (India 1830)", 6377276.345, 300.8017},
212 {
"Everest (1948)", 6377304.063, 300.8017},
213 {
"Modified Fischer 1960", 6378155.0, 298.3},
214 {
"Everest (Pakistan)", 6377309.613, 300.8017},
215 {
"Indonesian 1974", 6378160.0, 298.247},
216 {
"GRS 80", 6378137.0, 298.257222101},
217 {
"Helmert 1906", 6378200.0, 298.3},
218 {
"Hough 1960", 6378270.0, 297.0},
219 {
"International 1924", 6378388.0, 297.0},
220 {
"Krassovsky 1940", 6378245.0, 298.3},
221 {
"South American 1969", 6378160.0, 298.25},
222 {
"Everest (Malaysia 1969)", 6377295.664, 300.8017},
223 {
"Everest (Sabah Sarawak)", 6377298.556, 300.8017},
224 {
"WGS 72", 6378135.0, 298.26},
225 {
"WGS 84", 6378137.0, 298.257223563},
226 {
"Bessel 1841 (Namibia)", 6377483.865, 299.1528128},
227 {
"Everest (India 1956)", 6377301.243, 300.8017},
228 {
"Struve 1860", 6378298.3, 294.73}
232 short nDatums =
sizeof(gDatum) /
sizeof(
struct DATUM);
236 void datumParams(
short datum,
double *a,
double *es) {
237 extern struct DATUM const gDatum[];
238 extern struct ELLIPSOID const gEllipsoid[];
240 if (datum < nDatums) {
241 double f = 1.0 / gEllipsoid[gDatum[datum].ellipsoid].invf;
242 if (es) *es = 2 * f - f * f;
243 if (a) *a = gEllipsoid[gDatum[datum].ellipsoid].a;
245 double f = 1.0 / 298.257223563;
246 if (es) *es = 2 * f - f * f;
247 if (a) *a = 6378137.0;
251 static int datumNameCmp(
const char *n1,
const char *n2) {
257 else if (toupper(*n1) == toupper(*n2))
264 static int isWGS84(
int i) {
267 if (i == DATUM_INDEX_WGS84)
return i;
269 if (gDatum[i].ellipsoid != gDatum[DATUM_INDEX_WGS84].ellipsoid)
return i;
271 if (gDatum[i].dx != gDatum[DATUM_INDEX_WGS84].dx)
return i;
273 if (gDatum[i].dy != gDatum[DATUM_INDEX_WGS84].dy)
return i;
275 if (gDatum[i].dz != gDatum[DATUM_INDEX_WGS84].dz)
return i;
277 return DATUM_INDEX_WGS84;
281 int GetDatumIndex(
const char *str) {
283 while (i < (
int)nDatums) {
284 if (!datumNameCmp(str, gDatum[i].name)) {
295 void toDMS(
double a,
char *bufp,
int bufplen) {
298 int n = (int)((a - (
int)a) * 36000.0);
301 sprintf(bufp,
"%d%02d'%02d.%01d\"", (
int)(neg ? -a : a), m, s / 10, s % 10);
308 double fromDMS(
char *dms) {
311 char buf[20] = {
'\0'};
313 sscanf(dms,
"%d%[ ]%d%[ ']%lf%[ \"NSWEnswe]", &d, buf, &m, buf, &s, buf);
315 s = (double)(abs(d)) + ((
double)m + s / 60.0) / 60.0;
317 if (d >= 0 && strpbrk(buf,
"SWsw") == NULL)
return s;
326 void todmm(
int flag,
double a,
char *bufp,
int bufplen) {
330 int m = (int)((a - (
int)a) * 60000.0);
333 snprintf(bufp, bufplen,
"%d %02d.%03d'", (
int)a, m / 1000, m % 1000);
336 snprintf(bufp, bufplen,
"%02d %02d.%03d %c", (
int)a, m / 1000, (m % 1000),
338 }
else if (flag == 2) {
339 snprintf(bufp, bufplen,
"%03d %02d.%03d %c", (
int)a, m / 1000, (m % 1000),
345 void toDMM(
double a,
char *bufp,
int bufplen) {
346 todmm(0, a, bufp, bufplen);
355 void toSM(
double lat,
double lon,
double lat0,
double lon0,
double *x,
361 if ((lon * lon0 < 0.) && (fabs(lon - lon0) > 180.)) {
362 lon < 0.0 ? xlon += 360.0 : xlon -= 360.0;
365 const double z = WGS84_semimajor_axis_meters * mercator_k0;
367 *x = (xlon - lon0) * DEGREE * z;
370 const double s = sin(lat * DEGREE);
371 const double y3 = (.5 * log((1 + s) / (1 - s))) * z;
373 const double s0 = sin(lat0 * DEGREE);
374 const double y30 = (.5 * log((1 + s0) / (1 - s0))) * z;
378 double toSMcache_y30(
double lat0) {
379 const double z = WGS84_semimajor_axis_meters * mercator_k0;
380 const double s0 = sin(lat0 * DEGREE);
381 const double y30 = (.5 * log((1 + s0) / (1 - s0))) * z;
385 void toSMcache(
double lat,
double lon,
double y30,
double lon0,
double *x,
391 if ((lon * lon0 < 0.) && (fabs(lon - lon0) > 180.)) {
392 lon < 0.0 ? xlon += 360.0 : xlon -= 360.0;
395 const double z = WGS84_semimajor_axis_meters * mercator_k0;
397 *x = (xlon - lon0) * DEGREE * z;
400 const double s = sin(lat * DEGREE);
401 const double y3 = (.5 * log((1 + s) / (1 - s))) * z;
406 void fromSM(
double x,
double y,
double lat0,
double lon0,
double *lat,
408 const double z = WGS84_semimajor_axis_meters * mercator_k0;
421 const double s0 = sin(lat0 * DEGREE);
422 const double y0 = (.5 * log((1 + s0) / (1 - s0))) * z;
424 *lat = (2.0 * atan(exp((y0 + y) / z)) - PI / 2.) / DEGREE;
427 *lon = lon0 + (x / (DEGREE * z));
430 void fromSMR(
double x,
double y,
double lat0,
double lon0,
double axis_meters,
431 double *lat,
double *lon) {
432 const double s0 = sin(lat0 * DEGREE);
433 const double y0 = (.5 * log((1 + s0) / (1 - s0))) * axis_meters;
435 *lat = (2.0 * atan(exp((y0 + y) / axis_meters)) - PI / 2.) / DEGREE;
438 *lon = lon0 + (x / (DEGREE * axis_meters));
441 void toSM_ECC(
double lat,
double lon,
double lat0,
double lon0,
double *x,
443 const double f = 1.0 / WGSinvf;
444 const double e2 = 2 * f - f * f;
445 const double e = sqrt(e2);
447 const double z = WGS84_semimajor_axis_meters * mercator_k0;
449 *x = (lon - lon0) * DEGREE * z;
452 const double s = sin(lat * DEGREE);
455 const double s0 = sin(lat0 * DEGREE);
460 const double falsen = z * log(tan(PI / 4 + lat0 * DEGREE / 2) *
461 pow((1. - e * s0) / (1. + e * s0), e / 2.));
462 const double test = z * log(tan(PI / 4 + lat * DEGREE / 2) *
463 pow((1. - e * s) / (1. + e * s), e / 2.));
467 void fromSM_ECC(
double x,
double y,
double lat0,
double lon0,
double *lat,
469 const double f = 1.0 / WGSinvf;
470 const double es = 2 * f - f * f;
471 const double e = sqrt(es);
473 const double z = WGS84_semimajor_axis_meters * mercator_k0;
475 *lon = lon0 + (x / (DEGREE * z));
477 const double s0 = sin(lat0 * DEGREE);
479 const double falsen = z * log(tan(PI / 4 + lat0 * DEGREE / 2) *
480 pow((1. - e * s0) / (1. + e * s0), e / 2.));
481 const double t = exp((y + falsen) / (z));
482 const double xi = (PI / 2.) - 2.0 * atan(t);
486 double esf = (es / 2. + (5 * es * es / 24.) + (es * es * es / 12.) +
487 (13.0 * es * es * es * es / 360.)) *
489 esf += ((7. * es * es / 48.) + (29. * es * es * es / 240.) +
490 (811. * es * es * es * es / 11520.)) *
492 esf += ((7. * es * es * es / 120.) + (81 * es * es * es * es / 1120.) +
493 (4279. * es * es * es * es / 161280.)) *
496 *lat = -(xi + esf) / DEGREE;
505 void toPOLY(
double lat,
double lon,
double lat0,
double lon0,
double *x,
507 const double z = WGS84_semimajor_axis_meters * mercator_k0;
509 if (fabs((lat - lat0) * DEGREE) <= TOL) {
510 *x = (lon - lon0) * DEGREE * z;
514 const double E = (lon - lon0) * DEGREE * sin(lat * DEGREE);
515 const double cot = 1. / tan(lat * DEGREE);
517 *y = (lat * DEGREE) - (lat0 * DEGREE) + cot * (1. - cos(E));
533 void fromPOLY(
double x,
double y,
double lat0,
double lon0,
double *lat,
535 const double z = WGS84_semimajor_axis_meters * mercator_k0;
537 double yp = y - (lat0 * DEGREE * z);
538 if (fabs(yp) <= TOL) {
539 *lon = lon0 + (x / (DEGREE * z));
543 const double xp = x / z;
546 const double B = (xp * xp) + (yp * yp);
550 double tp = tan(lat3);
551 dphi = (yp * (lat3 * tp + 1.) - lat3 - .5 * (lat3 * lat3 + B) * tp) /
552 ((lat3 - yp) / tp - 1.);
554 }
while (fabs(dphi) > CONV && --i);
559 *lon = asin(xp * tan(lat3)) / sin(lat3);
563 *lat = lat3 / DEGREE;
579 void toTM(
float lat,
float lon,
float lat0,
float lon0,
double *x,
double *y) {
581 const double f = 1.0 / WGSinvf;
582 const double a = WGS84_semimajor_axis_meters;
583 const double k0 = 1.;
585 const double eccSquared = 2 * f - f * f;
586 const double eccPrimeSquared = (eccSquared) / (1 - eccSquared);
587 const double LatRad = lat * DEGREE;
588 const double LongOriginRad = lon0 * DEGREE;
589 const double LongRad = lon * DEGREE;
591 const double N = a / sqrt(1 - eccSquared * sin(LatRad) * sin(LatRad));
592 const double T = tan(LatRad) * tan(LatRad);
593 const double C = eccPrimeSquared * cos(LatRad) * cos(LatRad);
594 const double A = cos(LatRad) * (LongRad - LongOriginRad);
598 ((1 - eccSquared / 4 - 3 * eccSquared * eccSquared / 64 -
599 5 * eccSquared * eccSquared * eccSquared / 256) *
601 (3 * eccSquared / 8 + 3 * eccSquared * eccSquared / 32 +
602 45 * eccSquared * eccSquared * eccSquared / 1024) *
604 (15 * eccSquared * eccSquared / 256 +
605 45 * eccSquared * eccSquared * eccSquared / 1024) *
607 (35 * eccSquared * eccSquared * eccSquared / 3072) * sin(6 * LatRad));
610 (A + (1 - T + C) * A * A * A / 6 +
611 (5 - 18 * T + T * T + 72 * C - 58 * eccPrimeSquared) * A * A * A * A *
616 (MM + N * tan(LatRad) *
617 (A * A / 2 + (5 - T + 9 * C + 4 * C * C) * A * A * A * A / 24 +
618 (61 - 58 * T + T * T + 600 * C - 330 * eccPrimeSquared) * A *
619 A * A * A * A * A / 720)));
632 void fromTM(
double x,
double y,
double lat0,
double lon0,
double *lat,
634 const double rad2deg = 1. / DEGREE;
637 const double f = 1.0 / WGSinvf;
638 const double a = WGS84_semimajor_axis_meters;
639 const double k0 = 1.;
641 const double eccSquared = 2 * f - f * f;
643 const double eccPrimeSquared = (eccSquared) / (1 - eccSquared);
645 (1.0 - sqrt(1.0 - eccSquared)) / (1.0 + sqrt(1.0 - eccSquared));
647 const double MM = y / k0;
649 MM / (a * (1 - eccSquared / 4 - 3 * eccSquared * eccSquared / 64 -
650 5 * eccSquared * eccSquared * eccSquared / 256));
652 const double phi1Rad =
653 mu + (3 * e1 / 2 - 27 * e1 * e1 * e1 / 32) * sin(2 * mu) +
654 (21 * e1 * e1 / 16 - 55 * e1 * e1 * e1 * e1 / 32) * sin(4 * mu) +
655 (151 * e1 * e1 * e1 / 96) * sin(6 * mu);
657 const double N1 = a / sqrt(1 - eccSquared * sin(phi1Rad) * sin(phi1Rad));
658 const double T1 = tan(phi1Rad) * tan(phi1Rad);
659 const double C1 = eccPrimeSquared * cos(phi1Rad) * cos(phi1Rad);
660 const double R1 = a * (1 - eccSquared) /
661 pow(1 - eccSquared * sin(phi1Rad) * sin(phi1Rad), 1.5);
662 const double D = x / (N1 * k0);
665 (N1 * tan(phi1Rad) / R1) *
667 (5 + 3 * T1 + 10 * C1 - 4 * C1 * C1 - 9 * eccPrimeSquared) * D *
669 (61 + 90 * T1 + 298 * C1 + 45 * T1 * T1 - 252 * eccPrimeSquared -
671 D * D * D * D * D * D / 720);
672 *lat = lat0 + (*lat * rad2deg);
674 *lon = (D - (1 + 2 * T1 + C1) * D * D * D / 6 +
675 (5 - 2 * C1 + 28 * T1 - 3 * C1 * C1 + 8 * eccPrimeSquared +
677 D * D * D * D * D / 120) /
679 *lon = lon0 + *lon * rad2deg;
687 void cache_phi0(
double lat0,
double *sin_phi0,
double *cos_phi0) {
688 double phi0 = lat0 * DEGREE;
689 *sin_phi0 = sin(phi0);
690 *cos_phi0 = cos(phi0);
693 void toORTHO(
double lat,
double lon,
double sin_phi0,
double cos_phi0,
694 double lon0,
double *x,
double *y) {
695 const double z = WGS84_semimajor_axis_meters * mercator_k0;
699 if ((lon * lon0 < 0.) && (fabs(lon - lon0) > 180.))
700 lon < 0.0 ? xlon += 360.0 : xlon -= 360.0;
702 double theta = (xlon - lon0) * DEGREE;
703 double phi = lat * DEGREE;
704 double cos_phi = cos(phi);
706 double vy = sin(phi), vz = cos(theta) * cos_phi;
708 if (vy * sin_phi0 + vz * cos_phi0 < 0) {
713 double vx = sin(theta) * cos_phi;
714 double vw = vy * cos_phi0 - vz * sin_phi0;
720 void fromORTHO(
double x,
double y,
double lat0,
double lon0,
double *lat,
722 const double z = WGS84_semimajor_axis_meters * mercator_k0;
727 double phi0 = lat0 * DEGREE;
728 double d = 1 - vx * vx - vw * vw;
735 double sin_phi0 = sin(phi0), cos_phi0 = cos(phi0);
736 double vy = vw * cos_phi0 + sqrt(d) * sin_phi0;
737 double phi = asin(vy);
739 double vz = (vy * cos_phi0 - vw) / sin_phi0;
740 double theta = atan2(vx, vz);
743 *lon = theta / DEGREE + lon0;
746 double toPOLARcache_e(
double lat0) {
747 double pole = lat0 > 0 ? 90 : -90;
748 return tan((pole - lat0) * DEGREE / 2);
751 void toPOLAR(
double lat,
double lon,
double e,
double lat0,
double lon0,
752 double *x,
double *y) {
753 const double z = WGS84_semimajor_axis_meters * mercator_k0;
757 if ((lon * lon0 < 0.) && (fabs(lon - lon0) > 180.))
758 lon < 0.0 ? xlon += 360.0 : xlon -= 360.0;
760 double theta = (xlon - lon0) * DEGREE;
761 double pole = lat0 > 0 ? 90 : -90;
763 double d = tan((pole - lat) * DEGREE / 2);
765 *x = fabs(d) * sin(theta) * z;
766 *y = (e - d * cos(theta)) * z;
769 void fromPOLAR(
double x,
double y,
double lat0,
double lon0,
double *lat,
771 const double z = WGS84_semimajor_axis_meters * mercator_k0;
772 double pole = lat0 > 0 ? 90 : -90;
774 double e = tan((pole - lat0) * DEGREE / 2);
777 double yn = e - y / z;
778 double d = sqrt(xn * xn + yn * yn);
782 *lat = pole - atan(d) * 2 / DEGREE;
784 double theta = atan2(xn, yn);
785 *lon = theta / DEGREE + lon0;
788 static inline void toSTEREO1(
double &u,
double &v,
double &w,
double lat,
789 double lon,
double sin_phi0,
double cos_phi0,
793 if ((lon * lon0 < 0.) && (fabs(lon - lon0) > 180.))
794 lon < 0.0 ? xlon += 360.0 : xlon -= 360.0;
796 double theta = (xlon - lon0) * DEGREE, phi = lat * DEGREE;
797 double cos_phi = cos(phi), v0 = sin(phi), w0 = cos(theta) * cos_phi;
799 u = sin(theta) * cos_phi;
800 v = cos_phi0 * v0 - sin_phi0 * w0;
801 w = sin_phi0 * v0 + cos_phi0 * w0;
804 static inline void fromSTEREO1(
double *lat,
double *lon,
double lat0,
805 double lon0,
double u,
double v,
double w) {
806 double phi0 = lat0 * DEGREE;
807 double v0 = sin(phi0) * w + cos(phi0) * v;
808 double w0 = cos(phi0) * w - sin(phi0) * v;
809 double phi = asin(v0);
810 double theta = atan2(u, w0);
813 *lon = theta / DEGREE + lon0;
816 void toSTEREO(
double lat,
double lon,
double sin_phi0,
double cos_phi0,
817 double lon0,
double *x,
double *y) {
818 const double z = WGS84_semimajor_axis_meters * mercator_k0;
821 toSTEREO1(u, v, w, lat, lon, sin_phi0, cos_phi0, lon0);
823 double t = 2 / (w + 1);
828 void fromSTEREO(
double x,
double y,
double lat0,
double lon0,
double *lat,
830 const double z = WGS84_semimajor_axis_meters * mercator_k0;
834 double t = (x * x + y * y) / 4 + 1;
838 double w = 2 / t - 1;
840 fromSTEREO1(lat, lon, lat0, lon0, u, v, w);
843 void toGNO(
double lat,
double lon,
double sin_phi0,
double cos_phi0,
844 double lon0,
double *x,
double *y) {
845 const double z = WGS84_semimajor_axis_meters * mercator_k0;
848 toSTEREO1(u, v, w, lat, lon, sin_phi0, cos_phi0, lon0);
859 void fromGNO(
double x,
double y,
double lat0,
double lon0,
double *lat,
861 const double z = WGS84_semimajor_axis_meters * mercator_k0;
865 double w = 1 / sqrt(x * x + y * y + 1);
869 fromSTEREO1(lat, lon, lat0, lon0, u, v, w);
872 void toEQUIRECT(
double lat,
double lon,
double lat0,
double lon0,
double *x,
874 const double z = WGS84_semimajor_axis_meters * mercator_k0;
877 if ((lon * lon0 < 0.) && (fabs(lon - lon0) > 180.))
878 lon < 0.0 ? xlon += 360.0 : xlon -= 360.0;
880 *x = (xlon - lon0) * DEGREE * z;
881 *y = (lat - lat0) * DEGREE * z;
884 void fromEQUIRECT(
double x,
double y,
double lat0,
double lon0,
double *lat,
886 const double z = WGS84_semimajor_axis_meters * mercator_k0;
888 *lat = lat0 + (y / (DEGREE * z));
890 *lon = lon0 + (x / (DEGREE * z));
910 void MolodenskyTransform(
double lat,
double lon,
double *to_lat,
double *to_lon,
911 int from_datum_index,
int to_datum_index) {
915 if (from_datum_index < nDatums) {
916 const double from_lat = lat * DEGREE;
917 const double from_lon = lon * DEGREE;
918 const double from_f =
920 gEllipsoid[gDatum[from_datum_index].ellipsoid].invf;
921 const double from_esq = 2 * from_f - from_f * from_f;
922 const double from_a =
923 gEllipsoid[gDatum[from_datum_index].ellipsoid].a;
924 const double dx = gDatum[from_datum_index].dx;
925 const double dy = gDatum[from_datum_index].dy;
926 const double dz = gDatum[from_datum_index].dz;
928 1.0 / gEllipsoid[gDatum[to_datum_index].ellipsoid].invf;
930 gEllipsoid[gDatum[to_datum_index].ellipsoid].a;
931 const double da = to_a - from_a;
932 const double df = to_f - from_f;
933 const double from_h = 0;
935 const double slat = sin(from_lat);
936 const double clat = cos(from_lat);
937 const double slon = sin(from_lon);
938 const double clon = cos(from_lon);
939 const double ssqlat = slat * slat;
940 const double adb = 1.0 / (1.0 - from_f);
942 const double rn = from_a / sqrt(1.0 - from_esq * ssqlat);
944 from_a * (1. - from_esq) / pow((1.0 - from_esq * ssqlat), 1.5);
946 dlat = (((((-dx * slat * clon - dy * slat * slon) + dz * clat) +
947 (da * ((rn * from_esq * slat * clat) / from_a))) +
948 (df * (rm * adb + rn / adb) * slat * clat))) /
951 dlon = (-dx * slon + dy * clon) / ((rn + from_h) * clat);
954 *to_lon = lon + dlon / DEGREE;
955 *to_lat = lat + dlat / DEGREE;
993 #define HALFPI 1.5707963267948966
994 #define SPI 3.14159265359
995 #define TWOPI 6.2831853071795864769
996 #define ONEPI 3.14159265358979323846
997 #define MERI_TOL 1e-9
999 double adjlon(
double lon) {
1000 if (fabs(lon) <= SPI)
return (lon);
1002 lon -= TWOPI * floor(lon / TWOPI);
1016 void ll_gc_ll(
double lat,
double lon,
double brg,
double dist,
double *dlat,
1018 double th1, costh1, sinth1, sina12, cosa12, M, N, c1, c2, D, P, s1;
1021 if((brg == 90.) || (brg == 180.)){
1029 double phi1, lam1, phi2, lam2;
1034 double es, onef, f, f4;
1037 phi1 = lat * DEGREE;
1038 lam1 = lon * DEGREE;
1039 al12 = brg * DEGREE;
1040 geod_S = dist * 1852.0;
1049 geod_a = WGS84_semimajor_axis_meters;
1052 onef = sqrt(1. - es);
1056 al12 = adjlon(al12);
1057 signS = fabs(al12) > HALFPI ? 1 : 0;
1058 th1 = ellipse ? atan(onef * tan(phi1)) : phi1;
1061 if ((merid = fabs(sina12 = sin(al12)) < MERI_TOL)) {
1063 cosa12 = fabs(al12) < HALFPI ? 1. : -1.;
1067 M = costh1 * sina12;
1069 N = costh1 * cosa12;
1079 c2 = f4 * (1. - M * M);
1080 D = (1. - c2) * (1. - c2 - c1 * M);
1081 P = (1. + .5 * c1 * M) * c2 / D;
1087 s1 = (fabs(M) >= 1.) ? 0. : acos(M);
1088 s1 = sinth1 / sin(s1);
1089 s1 = (fabs(s1) >= 1.) ? 0. : acos(s1);
1095 double d, sind, u, V, X, ds, cosds, sinds, ss, de;
1100 d = geod_S / (D * geod_a);
1104 X = c2 * c2 * (sind = sin(d)) * cos(d) * (2. * V * V - 1.);
1105 ds = d + X - 2. * P * V * (1. - 2. * P * cos(u)) * sind;
1108 ds = geod_S / geod_a;
1109 if (signS) ds = -ds;
1113 if (signS) sinds = -sinds;
1114 al21 = N * cosds - sinth1 * sinds;
1116 phi2 = atan(tan(HALFPI + s1 - ds) / onef);
1134 al21 = atan(M / al21);
1135 if (al21 > 0) al21 += PI;
1136 if (al12 < 0.) al21 -= PI;
1137 al21 = adjlon(al21);
1138 phi2 = atan(-(sinth1 * cosds + N * sinds) * sin(al21) /
1139 (ellipse ? onef * M : M));
1140 de = atan2(sinds * sina12, (costh1 * cosds - sinth1 * sinds * cosa12));
1143 de += c1 * ((1. - c2) * ds + c2 * sinds * cos(ss));
1145 de -= c1 * ((1. - c2) * ds - c2 * sinds * cos(ss));
1148 lam2 = adjlon(lam1 + de);
1151 *dlat = phi2 / DEGREE;
1152 *dlon = lam2 / DEGREE;
1155 void ll_gc_ll_reverse(
double lat1,
double lon1,
double lat2,
double lon2,
1156 double *bearing,
double *dist) {
1159 if ((fabs(lon2 - lon1) < 0.1) && (fabs(lat2 - lat1) < 0.1)) {
1160 DistanceBearingMercator(lat2, lon2, lat1, lon1, bearing, dist);
1167 double phi1, lam1, phi2, lam2;
1172 double es, onef, f, f64, f2, f4;
1175 phi1 = lat1 * DEGREE;
1176 lam1 = lon1 * DEGREE;
1177 phi2 = lat2 * DEGREE;
1178 lam2 = lon2 * DEGREE;
1182 double th1, th2, thm, dthm, dlamm, dlam, sindlamm, costhm, sinthm,
1183 cosdthm, sindthm, L, E, cosd, d, X, Y, T, sind, tandlammp, u, v, D, A,
1192 geod_a = WGS84_semimajor_axis_meters;
1195 onef = sqrt(1. - es);
1199 f64 = geod_f * geod_f / 64;
1202 th1 = atan(onef * tan(phi1));
1203 th2 = atan(onef * tan(phi2));
1208 thm = .5 * (th1 + th2);
1209 dthm = .5 * (th2 - th1);
1210 dlamm = .5 * (dlam = adjlon(lam2 - lam1));
1211 if (fabs(dlam) < DTOL && fabs(dthm) < DTOL) {
1212 al12 = al21 = geod_S = 0.;
1213 if (bearing) *bearing = 0.;
1214 if (dist) *dist = 0.;
1217 sindlamm = sin(dlamm);
1220 cosdthm = cos(dthm);
1221 sindthm = sin(dthm);
1222 L = sindthm * sindthm +
1223 (cosdthm * cosdthm - sinthm * sinthm) * sindlamm * sindlamm;
1224 d = acos(cosd = 1 - L - L);
1228 Y = sinthm * cosdthm;
1229 Y *= (Y + Y) / (1. - L);
1230 T = sindthm * costhm;
1238 geod_S = geod_a * sind *
1239 (T - f4 * (T * X - Y) +
1240 f64 * (X * (A + (T - .5 * (A - E)) * X) - Y * (B + E * Y) +
1243 tan(.5 * (dlam - .25 * (Y + Y - E * (4. - X)) *
1244 (f2 * T + f64 * (32. * T - (20. * T - A) * X -
1248 geod_S = geod_a * d;
1249 tandlammp = tan(dlamm);
1251 u = atan2(sindthm, (tandlammp * costhm));
1252 v = atan2(cosdthm, (tandlammp * sinthm));
1253 al12 = adjlon(TWOPI + v - u);
1254 al21 = adjlon(TWOPI - v - u);
1255 if (al12 < 0) al12 += 2 * PI;
1256 if (bearing) *bearing = al12 / DEGREE;
1257 if (dist) *dist = geod_S / 1852.0;
1262 void PositionBearingDistanceMercator(
double lat,
double lon,
double brg,
1263 double dist,
double *dlat,
double *dlon) {
1264 ll_gc_ll(lat, lon, brg, dist, dlat, dlon);
1267 double DistLoxodrome(
double slat,
double slon,
double dlat,
double dlon) {
1269 60 * sqrt(pow(slat - dlat, 2) +
1270 pow((slon - dlon) * cos((slat + dlat) / 2 * DEGREE), 2));
1272 if (slon * dlon < 0){
1273 if (slon < 0) slon += 360.;
1274 else if (dlon < 0) dlon += 360.;
1276 60 * sqrt(pow(slat - dlat, 2) +
1277 pow((slon - dlon) * cos((slat + dlat) / 2 * DEGREE), 2));
1278 return (std::min)(dist, distrtw);
1293 double DistGreatCircle(
double slat,
double slon,
double dlat,
double dlon) {
1295 d5 = DistLoxodrome(slat, slon, dlat, dlon);
1303 double phi1, lam1, phi2, lam2;
1308 double es, onef, f, f64, f4;
1310 phi1 = slat * DEGREE;
1311 lam1 = slon * DEGREE;
1312 phi2 = dlat * DEGREE;
1313 lam2 = dlon * DEGREE;
1317 double th1, th2, thm, dthm, dlamm, dlam, sindlamm, costhm, sinthm, cosdthm,
1318 sindthm, L, E, cosd, d, X, Y, T, sind, D, A, B;
1327 geod_a = WGS84_semimajor_axis_meters;
1330 onef = sqrt(1. - es);
1334 f64 = geod_f * geod_f / 64;
1337 th1 = atan(onef * tan(phi1));
1338 th2 = atan(onef * tan(phi2));
1343 thm = .5 * (th1 + th2);
1344 dthm = .5 * (th2 - th1);
1345 dlamm = .5 * (dlam = adjlon(lam2 - lam1));
1346 if (fabs(dlam) < DTOL && fabs(dthm) < DTOL) {
1349 sindlamm = sin(dlamm);
1352 cosdthm = cos(dthm);
1353 sindthm = sin(dthm);
1354 L = sindthm * sindthm +
1355 (cosdthm * cosdthm - sinthm * sinthm) * sindlamm * sindlamm;
1356 d = acos(cosd = 1 - L - L);
1361 Y = sinthm * cosdthm;
1362 Y *= (Y + Y) / (1. - L);
1363 T = sindthm * costhm;
1371 geod_S = geod_a * sind *
1372 (T - f4 * (T * X - Y) +
1373 f64 * (X * (A + (T - .5 * (A - E)) * X) - Y * (B + E * Y) +
1380 geod_S = geod_a * d;
1389 d5 = geod_S / 1852.0;
1394 void DistanceBearingMercator(
double lat1,
double lon1,
double lat0,
double lon0,
1395 double *brg,
double *dist) {
1397 double latm = (lat0 + lat1) / 2 * DEGREE;
1398 double delta_lat = (lat1 - lat0);
1399 double delta_lon = (lon1 - lon0);
1400 double ex_lat0, ex_lat1;
1401 double bearing, distance;
1403 if (delta_lon < -180) delta_lon += 360;
1404 if (delta_lon > 180) delta_lon -= 360;
1408 else if (fabs(delta_lat) != .0)
1409 bearing = atan(delta_lon * cos(latm) / delta_lat);
1413 distance = sqrt(pow(delta_lat, 2) + pow(delta_lon * cos(latm), 2));
1418 if (delta_lat != 0.) {
1419 ex_lat0 = 10800 / PI * log(tan(PI / 4 + lat0 * DEGREE / 2));
1420 ex_lat1 = 10800 / PI * log(tan(PI / 4 + lat1 * DEGREE / 2));
1421 bearing = atan(delta_lon * 60 / (ex_lat1 - ex_lat0));
1422 distance = fabs(delta_lat / cos(bearing));
1428 bearing = fabs(bearing);
1430 if (delta_lon < 0) bearing = 2 * PI - bearing;
1434 bearing = PI - bearing;
1436 bearing = PI + bearing;
1439 if (brg) *brg = bearing * RADIAN;
1440 if (dist) *dist = distance * 60;
1466 double my_fit_function(
double tx,
double ty,
int n_par,
double *p) {
1467 double ret = p[0] + p[1] * tx;
1469 if (n_par > 2) ret += p[2] * ty;
1471 ret += p[3] * tx * tx;
1472 ret += p[4] * tx * ty;
1473 ret += p[5] * ty * ty;
1476 ret += p[6] * tx * tx * tx;
1477 ret += p[7] * tx * tx * ty;
1478 ret += p[8] * tx * ty * ty;
1479 ret += p[9] * ty * ty * ty;
1485 int Georef_Calculate_Coefficients_Onedir(
int n_points,
int n_par,
double *tx,
1486 double *ty,
double *y,
double *p,
1487 double hintp0,
double hintp1,
1502 lm_initialize_control(&control);
1504 for (i = 0; i < 12; i++) p[i] = 0.;
1513 data.user_func = my_fit_function;
1518 data.print_flag = 0;
1522 lm_minimize(n_points, n_par, p, lm_evaluate_default, lm_print_default, &data,
1530 return control.info;
1533 int Georef_Calculate_Coefficients(
struct GeoRef *cp,
int nlin_lon) {
1535 for (
int i = 0; i < 10; ++i)
1536 cp->pwx[i] = cp->pwy[i] = cp->wpx[i] = cp->wpy[i] = 0.0;
1540 switch (cp->order) {
1555 const int mp_lat = mp;
1558 const int mp_lon = nlin_lon ? 2 : mp;
1561 std::vector<double> pnull(cp->count, 1.0);
1566 int r1 = Georef_Calculate_Coefficients_Onedir(
1567 cp->count, mp_lon, cp->tx, cp->ty, cp->lon, cp->pwx,
1569 (cp->txmin * (cp->lonmax - cp->lonmin) / (cp->txmax - cp->txmin)),
1570 (cp->lonmax - cp->lonmin) / (cp->txmax - cp->txmin), 0.);
1574 double *px = nlin_lon ? &pnull[0] : cp->tx;
1576 int r2 = Georef_Calculate_Coefficients_Onedir(
1577 cp->count, mp_lat, px, cp->ty, cp->lat, cp->pwy,
1579 (cp->tymin * (cp->latmax - cp->latmin) / (cp->tymax - cp->tymin)),
1580 0., (cp->latmax - cp->latmin) / (cp->tymax - cp->tymin));
1585 int r3 = Georef_Calculate_Coefficients_Onedir(
1586 cp->count, mp_lon, cp->lon, cp->lat, cp->tx, cp->wpx,
1588 ((cp->txmax - cp->txmin) * cp->lonmin / (cp->lonmax - cp->lonmin)),
1589 (cp->txmax - cp->txmin) / (cp->lonmax - cp->lonmin), 0.0);
1593 px = nlin_lon ? &pnull[0] : cp->lon;
1595 int r4 = Georef_Calculate_Coefficients_Onedir(
1596 cp->count, mp_lat, &pnull[0] , cp->lat, cp->ty, cp->wpy,
1598 ((cp->tymax - cp->tymin) * cp->latmin / (cp->latmax - cp->latmin)),
1599 0.0, (cp->tymax - cp->tymin) / (cp->latmax - cp->latmin));
1601 if ((r1) && (r1 < 4) && (r2) && (r2 < 4) && (r3) && (r3 < 4) && (r4) &&
1608 int Georef_Calculate_Coefficients_Proj(
struct GeoRef *cp) {
1613 cp->pwx[6] = cp->pwy[6] = cp->wpx[6] = cp->wpy[6] = 0.;
1614 cp->pwx[7] = cp->pwy[7] = cp->wpx[7] = cp->wpy[7] = 0.;
1615 cp->pwx[8] = cp->pwy[8] = cp->wpx[8] = cp->wpy[8] = 0.;
1616 cp->pwx[9] = cp->pwy[9] = cp->wpx[9] = cp->wpy[9] = 0.;
1617 cp->pwx[3] = cp->pwy[3] = cp->wpx[3] = cp->wpy[3] = 0.;
1618 cp->pwx[4] = cp->pwy[4] = cp->wpx[4] = cp->wpy[4] = 0.;
1619 cp->pwx[5] = cp->pwy[5] = cp->wpx[5] = cp->wpy[5] = 0.;
1620 cp->pwx[0] = cp->pwy[0] = cp->wpx[0] = cp->wpy[0] = 0.;
1621 cp->pwx[1] = cp->pwy[1] = cp->wpx[1] = cp->wpy[1] = 0.;
1622 cp->pwx[2] = cp->pwy[2] = cp->wpx[2] = cp->wpy[2] = 0.;
1629 r1 = Georef_Calculate_Coefficients_Onedir(
1630 cp->count, mp, cp->tx, cp->ty, cp->lon, cp->pwx,
1632 (cp->txmin * (cp->lonmax - cp->lonmin) / (cp->txmax - cp->txmin)),
1633 (cp->lonmax - cp->lonmin) / (cp->txmax - cp->txmin), 0.);
1635 r2 = Georef_Calculate_Coefficients_Onedir(
1636 cp->count, mp, cp->tx, cp->ty, cp->lat, cp->pwy,
1638 (cp->tymin * (cp->latmax - cp->latmin) / (cp->tymax - cp->tymin)),
1639 0., (cp->latmax - cp->latmin) / (cp->tymax - cp->tymin));
1644 r3 = Georef_Calculate_Coefficients_Onedir(
1645 cp->count, mp, cp->lon, cp->lat, cp->tx, cp->wpx,
1647 ((cp->txmax - cp->txmin) * cp->lonmin / (cp->lonmax - cp->lonmin)),
1648 (cp->txmax - cp->txmin) / (cp->lonmax - cp->lonmin), 0.0);
1650 r4 = Georef_Calculate_Coefficients_Onedir(
1651 cp->count, mp, cp->lon, cp->lat, cp->ty, cp->wpy,
1653 ((cp->tymax - cp->tymin) * cp->latmin / (cp->latmax - cp->latmin)),
1654 0.0, (cp->tymax - cp->tymin) / (cp->latmax - cp->latmin));
1656 if ((r1) && (r1 < 4) && (r2) && (r2 < 4) && (r3) && (r3 < 4) && (r4) &&
1669 void lm_evaluate_default(
double *par,
int m_dat,
double *fvec,
void *data,
1689 for (
int i = 0; i < m_dat; i++) {
1690 fvec[i] = mydata->user_y[i] - mydata->user_func(mydata->user_tx[i],
1692 mydata->n_par, par);
1698 void lm_print_default(
int n_par,
double *par,
int m_dat,
double *fvec,
1699 void *data,
int iflag,
int iter,
int nfev)
1710 if (mydata->print_flag) {
1712 printf(
"trying step in gradient direction\n");
1713 }
else if (iflag == 1) {
1714 printf(
"determining gradient (iteration %d)\n", iter);
1715 }
else if (iflag == 0) {
1716 printf(
"starting minimization\n");
1717 }
else if (iflag == -1) {
1718 printf(
"terminated after %d evaluations\n", nfev);
1722 for (
int i = 0; i < n_par; ++i) printf(
" %12g", par[i]);
1723 printf(
" => norm: %12g\n", lm_enorm(m_dat, fvec));
1726 printf(
" fitting data as follows:\n");
1727 for (
int i = 0; i < m_dat; ++i) {
1728 const double tx = (mydata->user_tx)[i];
1729 const double ty = (mydata->user_ty)[i];
1730 const double y = (mydata->user_y)[i];
1731 const double f = mydata->user_func(tx, ty, mydata->n_par, par);
1733 " tx[%2d]=%8g ty[%2d]=%8g y=%12g fit=%12g "
1735 i, tx, i, ty, y, f, y - f);
1746 control->maxcall = 100;
1747 control->epsilon = 1.e-10;
1748 control->stepbound = 100;
1749 control->ftol = 1.e-14;
1750 control->xtol = 1.e-14;
1751 control->gtol = 1.e-14;
1754 void lm_minimize(
int m_dat,
int n_par,
double *par, lm_evaluate_ftype *evaluate,
1755 lm_print_ftype *printout,
void *data,
1757 std::vector<double> fvec(m_dat);
1758 std::vector<double> diag(n_par);
1759 std::vector<double> qtf(n_par);
1760 std::vector<double> fjac(n_par * m_dat);
1761 std::vector<double> wa1(n_par);
1762 std::vector<double> wa2(n_par);
1763 std::vector<double> wa3(n_par);
1764 std::vector<double> wa4(m_dat);
1765 std::vector<int> ipvt(n_par);
1773 lm_lmdif(m_dat, n_par, par, &fvec[0], control->ftol, control->xtol,
1774 control->gtol, control->maxcall * (n_par + 1), control->epsilon,
1775 &diag[0], 1, control->stepbound, &(control->info), &(control->nfev),
1776 &fjac[0], &ipvt[0], &qtf[0], &wa1[0], &wa2[0], &wa3[0], &wa4[0],
1777 evaluate, printout, data);
1779 (*printout)(n_par, par, m_dat, &fvec[0], data, -1, 0, control->nfev);
1780 control->fnorm = lm_enorm(m_dat, &fvec[0]);
1781 if (control->info < 0) control->info = 10;
1786 const char *lm_infmsg[] = {
1787 "improper input parameters",
1788 "the relative error in the sum of squares is at most tol",
1789 "the relative error between x and the solution is at most tol",
1790 "both errors are at most tol",
1791 "fvec is orthogonal to the columns of the jacobian to machine precision",
1792 "number of calls to fcn has reached or exceeded 200*(n+1)",
1793 "ftol is too small. no further reduction in the sum of squares is possible",
1794 "xtol too small. no further improvement in approximate solution x possible",
1795 "gtol too small. no further improvement in approximate solution x possible",
1796 "not enough memory",
1797 "break requested within function evaluation"};
1799 const char *lm_shortmsg[] = {
"invalid input",
"success (f)",
"success (p)",
1800 "success (f,p)",
"degenerate",
"call limit",
1801 "failed (f)",
"failed (p)",
"failed (o)",
1802 "no memory",
"user break"};
1815 #define LM_MACHEP 1.2e-16
1816 #define LM_DWARF 1.0e-38
1823 #define LM_SQRT_DWARF 3.834e-20
1824 #define LM_SQRT_GIANT 1.304e19
1826 void lm_qrfac(
int m,
int n,
double *a,
int pivot,
int *ipvt,
double *rdiag,
1827 double *acnorm,
double *wa);
1828 void lm_qrsolv(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
double *qtb,
1829 double *x,
double *sdiag,
double *wa);
1830 void lm_lmpar(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
double *qtb,
1831 double delta,
double *par,
double *x,
double *sdiag,
double *wa1,
1834 #define MIN(a, b) (((a) <= (b)) ? (a) : (b))
1835 #define MAX(a, b) (((a) >= (b)) ? (a) : (b))
1836 #define SQR(x) (x) * (x)
1840 void lm_lmdif(
int m,
int n,
double *x,
double *fvec,
double ftol,
double xtol,
1841 double gtol,
int maxfev,
double epsfcn,
double *diag,
int mode,
1842 double factor,
int *info,
int *nfev,
double *fjac,
int *ipvt,
1843 double *qtf,
double *wa1,
double *wa2,
double *wa3,
double *wa4,
1844 lm_evaluate_ftype *evaluate, lm_print_ftype *printout,
2026 if ((n <= 0) || (m < n) || (ftol < 0.) || (xtol < 0.) || (gtol < 0.) ||
2027 (maxfev <= 0) || (factor <= 0.)) {
2034 for (
int j = 0; j < n; j++)
2036 if (diag[j] <= 0.0) {
2049 (*evaluate)(x, m, fvec, data, info);
2050 (*printout)(n, x, m, fvec, data, 0, 0, ++(*nfev));
2051 if (*info < 0)
return;
2058 double fnorm = lm_enorm(m, fvec);
2059 const double eps = sqrt(MAX(
2065 printf(
"lmdif/ outer loop iter=%d nfev=%d fnorm=%.10e\n", iter, *nfev,
2071 for (
int j = 0; j < n; j++) {
2072 const double temp = x[j];
2073 const double step = eps * fabs(temp) == 0.0 ? eps : eps * fabs(temp);
2077 (*evaluate)(x, m, wa4, data, info);
2078 (*printout)(n, x, m, wa4, data, 1, iter, ++(*nfev));
2079 if (*info < 0)
return;
2081 for (
int i = 0; i < m; i++) fjac[j * m + i] = (wa4[i] - fvec[i]) / step;
2085 for (i = 0; i < m; i++) {
2086 for (j = 0; j < n; j++) printf(
"%.5e ", y[j * m + i]);
2093 lm_qrfac(m, n, fjac, 1, ipvt, wa1, wa2, wa3);
2102 for (
int j = 0; j < n; j++) {
2104 if (wa2[j] == 0.) diag[j] = 1.;
2111 for (
int j = 0; j < n; j++) wa3[j] = diag[j] * x[j];
2113 xnorm = lm_enorm(n, wa3);
2114 delta = factor * xnorm;
2115 if (delta == 0.) delta = factor;
2120 for (
int i = 0; i < m; i++) wa4[i] = fvec[i];
2122 for (
int j = 0; j < n; j++) {
2123 double temp3 = fjac[j * m + j];
2126 for (
int i = j; i < m; i++) sum += fjac[j * m + i] * wa4[i];
2127 double temp = -sum / temp3;
2128 for (
int i = j; i < m; i++) wa4[i] += fjac[j * m + i] * temp;
2130 fjac[j * m + j] = wa1[j];
2138 for (
int j = 0; j < n; j++) {
2139 if (wa2[ipvt[j]] == 0)
continue;
2142 for (
int i = 0; i <= j; i++) sum += fjac[j * m + i] * qtf[i] / fnorm;
2143 gnorm = MAX(gnorm, fabs(sum / wa2[ipvt[j]]));
2147 if (gnorm <= gtol) {
2155 for (
int j = 0; j < n; j++) diag[j] = MAX(diag[j], wa2[j]);
2159 const double kP0001 = 1.0e-4;
2163 printf(
"lmdif/ inner loop iter=%d nfev=%d\n", iter, *nfev);
2168 lm_lmpar(n, fjac, m, ipvt, diag, qtf, delta, &par, wa1, wa2, wa3, wa4);
2172 for (
int j = 0; j < n; j++) {
2174 wa2[j] = x[j] + wa1[j];
2175 wa3[j] = diag[j] * wa1[j];
2177 double pnorm = lm_enorm(n, wa3);
2182 delta = MIN(delta, pnorm);
2187 (*evaluate)(wa2, m, wa4, data, info);
2188 (*printout)(n, x, m, wa4, data, 2, iter, ++(*nfev));
2189 if (*info < 0)
return;
2191 double fnorm1 = lm_enorm(m, wa4);
2194 "lmdif/ pnorm %.10e fnorm1 %.10e fnorm %.10e"
2195 " delta=%.10e par=%.10e\n",
2196 pnorm, fnorm1, fnorm, delta, par);
2200 const double kP1 = 0.1;
2201 const double actred = kP1 * fnorm1 < fnorm ? 1 - SQR(fnorm1 / fnorm) : -1;
2206 for (
int j = 0; j < n; j++) {
2208 for (
int i = 0; i <= j; i++) wa3[i] += fjac[j * m + i] * wa1[ipvt[j]];
2210 const double temp1 = lm_enorm(n, wa3) / fnorm;
2211 const double temp2 = sqrt(par) * pnorm / fnorm;
2212 const double prered = SQR(temp1) + 2 * SQR(temp2);
2213 const double dirder = -(SQR(temp1) + SQR(temp2));
2217 ratio = prered != 0 ? actred / prered : 0.0;
2220 "lmdif/ actred=%.10e prered=%.10e ratio=%.10e"
2221 " sq(1)=%.10e sq(2)=%.10e dd=%.10e\n",
2222 actred, prered, prered != 0 ? ratio : 0., SQR(temp1), SQR(temp2),
2227 const double kP5 = 0.5;
2228 const double kP25 = 0.25;
2229 const double kP75 = 0.75;
2231 if (ratio <= kP25) {
2233 actred >= 0.0 ? kP5 : kP5 * dirder / (dirder + kP5 * actred);
2234 if (kP1 * fnorm1 >= fnorm || temp < kP1) temp = kP1;
2235 delta = temp * MIN(delta, pnorm / kP1);
2237 }
else if (par == 0. || ratio >= kP75) {
2238 delta = pnorm / kP5;
2243 if (ratio >= kP0001) {
2246 for (
int j = 0; j < n; j++) {
2248 wa2[j] = diag[j] * x[j];
2250 for (
int i = 0; i < m; i++) fvec[i] = wa4[i];
2251 xnorm = lm_enorm(n, wa2);
2257 printf(
"ATTN: iteration considered unsuccessful\n");
2264 if (fabs(actred) <= ftol && prered <= ftol && kP5 * ratio <= 1) *info = 1;
2265 if (delta <= xtol * xnorm) *info += 2;
2266 if (*info != 0)
return;
2270 if (*nfev >= maxfev) *info = 5;
2271 if (fabs(actred) <= LM_MACHEP && prered <= LM_MACHEP && kP5 * ratio <= 1)
2273 if (delta <= LM_MACHEP * xnorm) *info = 7;
2274 if (gnorm <= LM_MACHEP) *info = 8;
2275 if (*info != 0)
return;
2279 }
while (ratio < kP0001);
2286 void lm_lmpar(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
double *qtb,
2287 double delta,
double *par,
double *x,
double *sdiag,
double *wa1,
2374 for (
int j = 0; j < n; j++) {
2376 if (r[j * ldr + j] == 0 && nsing == n) nsing = j;
2377 if (nsing < n) wa1[j] = 0;
2380 printf(
"nsing %d ", nsing);
2382 for (
int j = nsing - 1; j >= 0; j--) {
2383 wa1[j] = wa1[j] / r[j + ldr * j];
2384 const double temp = wa1[j];
2385 for (
int i = 0; i < j; i++) wa1[i] -= r[j * ldr + i] * temp;
2388 for (
int j = 0; j < n; j++) x[ipvt[j]] = wa1[j];
2395 for (
int j = 0; j < n; j++) wa2[j] = diag[j] * x[j];
2397 double dxnorm = lm_enorm(n, wa2);
2398 double fp = dxnorm - delta;
2399 const double kP1 = 0.1;
2400 if (fp <= kP1 * delta) {
2402 printf(
"lmpar/ terminate (fp<delta/10\n");
2414 for (
int j = 0; j < n; j++) wa1[j] = diag[ipvt[j]] * wa2[ipvt[j]] / dxnorm;
2416 for (
int j = 0; j < n; j++) {
2418 for (
int i = 0; i < j; i++) sum += r[j * ldr + i] * wa1[i];
2419 wa1[j] = (wa1[j] - sum) / r[j + ldr * j];
2421 const double temp = lm_enorm(n, wa1);
2422 parl = fp / delta / temp / temp;
2427 for (
int j = 0; j < n; j++) {
2429 for (
int i = 0; i <= j; i++) sum += r[j * ldr + i] * qtb[i];
2430 wa1[j] = sum / diag[ipvt[j]];
2432 const double gnorm = lm_enorm(n, wa1);
2434 gnorm / delta == 0.0 ? LM_DWARF / MIN(delta, kP1) : gnorm / delta;
2439 *par = MAX(*par, parl);
2440 *par = MIN(*par, paru);
2441 if (*par == 0.) *par = gnorm / dxnorm;
2443 printf(
"lmpar/ parl %.4e par %.4e paru %.4e\n", parl, *par, paru);
2450 const double kP001 = 0.001;
2451 if (*par == 0.) *par = MAX(LM_DWARF, kP001 * paru);
2452 double temp = sqrt(*par);
2453 for (
int j = 0; j < n; j++) wa1[j] = temp * diag[j];
2454 lm_qrsolv(n, r, ldr, ipvt, wa1, qtb, x, sdiag, wa2);
2455 for (
int j = 0; j < n; j++) wa2[j] = diag[j] * x[j];
2456 dxnorm = lm_enorm(n, wa2);
2457 const double fp_old = fp;
2458 fp = dxnorm - delta;
2464 if (fabs(fp) <= kP1 * delta ||
2465 (parl == 0. && fp <= fp_old && fp_old < 0.) || iter == 10)
2470 for (
int j = 0; j < n; j++) wa1[j] = diag[ipvt[j]] * wa2[ipvt[j]] / dxnorm;
2472 for (
int j = 0; j < n; j++) {
2473 wa1[j] = wa1[j] / sdiag[j];
2474 for (
int i = j + 1; i < n; i++) wa1[i] -= r[j * ldr + i] * wa1[j];
2476 temp = lm_enorm(n, wa1);
2477 double parc = fp / delta / temp / temp;
2482 parl = MAX(parl, *par);
2484 paru = MIN(paru, *par);
2489 *par = MAX(parl, *par + parc);
2493 void lm_qrfac(
int m,
int n,
double *a,
int pivot,
int *ipvt,
double *rdiag,
2494 double *acnorm,
double *wa) {
2550 for (
int j = 0; j < n; j++) {
2551 acnorm[j] = lm_enorm(m, &a[j * m]);
2552 rdiag[j] = acnorm[j];
2554 if (pivot) ipvt[j] = j;
2562 const int minmn = MIN(m, n);
2563 for (
int j = 0; j < minmn; j++) {
2565 if (!pivot)
goto pivot_ok;
2569 for (k = j + 1; k < n; k++)
2570 if (rdiag[k] > rdiag[kmax]) kmax = k;
2571 if (kmax == j)
goto pivot_ok;
2573 for (
int i = 0; i < m; i++) {
2574 std::swap(a[j * m + i], a[kmax * m + i]);
2579 rdiag[kmax] = rdiag[j];
2581 std::swap(ipvt[j], ipvt[kmax]);
2591 double ajnorm = lm_enorm(m - j, &a[j * m + j]);
2597 if (a[j * m + j] < 0.) ajnorm = -ajnorm;
2598 for (
int i = j; i < m; i++) a[j * m + i] /= ajnorm;
2604 for (k = j + 1; k < n; k++) {
2607 for (
int i = j; i < m; i++) sum += a[j * m + i] * a[k * m + i];
2609 double temp = sum / a[j + m * j];
2611 for (
int i = j; i < m; i++) a[k * m + i] -= temp * a[j * m + i];
2613 if (pivot && rdiag[k] != 0.) {
2614 temp = a[m * k + j] / rdiag[k];
2615 temp = MAX(0., 1 - temp * temp);
2616 rdiag[k] *= sqrt(temp);
2617 temp = rdiag[k] / wa[k];
2618 const double kP05 = 0.05;
2619 if (kP05 * SQR(temp) <= LM_MACHEP) {
2620 rdiag[k] = lm_enorm(m - j - 1, &a[m * k + j + 1]);
2630 void lm_qrsolv(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
double *qtb,
2631 double *x,
double *sdiag,
double *wa) {
2699 for (
int j = 0; j < n; j++) {
2700 for (
int i = j; i < n; i++) r[j * ldr + i] = r[i * ldr + j];
2701 x[j] = r[j * ldr + j];
2710 for (
int j = 0; j < n; j++) {
2716 if (diag[ipvt[j]] == 0.)
goto L90;
2717 for (
int k = j; k < n; k++) sdiag[k] = 0.;
2718 sdiag[j] = diag[ipvt[j]];
2724 for (
int k = j; k < n; k++) {
2725 const double p25 = 0.25;
2726 const double p5 = 0.5;
2731 if (sdiag[k] == 0.)
continue;
2732 const int kk = k + ldr * k;
2735 if (fabs(r[kk]) < fabs(sdiag[k])) {
2736 const double cotan = r[kk] / sdiag[k];
2737 sin = p5 / sqrt(p25 + p25 * SQR(cotan));
2740 const double tan = sdiag[k] / r[kk];
2741 cos = p5 / sqrt(p25 + p25 * SQR(tan));
2748 r[kk] = cos * r[kk] + sin * sdiag[k];
2749 double temp = cos * wa[k] + sin * qtbpj;
2750 qtbpj = -sin * wa[k] + cos * qtbpj;
2755 for (
int i = k + 1; i < n; i++) {
2756 temp = cos * r[k * ldr + i] + sin * sdiag[i];
2757 sdiag[i] = -sin * r[k * ldr + i] + cos * sdiag[i];
2758 r[k * ldr + i] = temp;
2766 sdiag[j] = r[j * ldr + j];
2767 r[j * ldr + j] = x[j];
2774 for (
int j = 0; j < n; j++) {
2775 if (sdiag[j] == 0. && nsing == n) nsing = j;
2776 if (nsing < n) wa[j] = 0;
2779 for (
int j = nsing - 1; j >= 0; j--) {
2781 for (
int i = j + 1; i < nsing; i++) sum += r[j * ldr + i] * wa[i];
2782 wa[j] = (wa[j] - sum) / sdiag[j];
2787 for (
int j = 0; j < n; j++) x[ipvt[j]] = wa[j];
2790 double lm_enorm(
int n,
double *x) {
2812 double s1 = 0.0, s2 = 0.0, s3 = 0.0;
2813 double x1max = 0.0, x3max = 0.0;
2814 const double agiant = LM_SQRT_GIANT / ((double)n);
2816 for (
int i = 0; i < n; i++) {
2817 double xabs = fabs(x[i]);
2818 if (xabs > LM_SQRT_DWARF && xabs < agiant) {
2824 if (xabs > LM_SQRT_DWARF) {
2827 s1 = 1 + s1 * SQR(x1max / xabs);
2830 s1 += SQR(xabs / x1max);
2836 s3 = 1 + s3 * SQR(x3max / xabs);
2840 s3 += SQR(xabs / x3max);
2847 if (s1 != 0)
return x1max * sqrt(s1 + (s2 / x1max) / x1max);
2849 if (s2 >= x3max)
return sqrt(s2 * (1 + (x3max / s2) * (x3max * s3)));
2850 return sqrt(x3max * ((s2 / x3max) + (x3max * s3)));
2853 return x3max * sqrt(s3);
2856 double lat_gc_crosses_meridian(
double lat1,
double lon1,
double lat2,
2857 double lon2,
double lon) {
2863 double dlon = lon * DEGREE;
2864 double dlat1 = lat1 * DEGREE;
2865 double dlat2 = lat2 * DEGREE;
2866 double dlon1 = lon1 * DEGREE;
2867 double dlon2 = lon2 * DEGREE;
2869 return RADIAN * atan((sin(dlat1) * cos(dlat2) * sin(dlon - dlon2) -
2870 sin(dlat2) * cos(dlat1) * sin(dlon - dlon1)) /
2871 (cos(dlat1) * cos(dlat2) * sin(dlon1 - dlon2)));
2874 double lat_rl_crosses_meridian(
double lat1,
double lon1,
double lat2,
2875 double lon2,
double lon) {
2883 DistanceBearingMercator(lat2, lon2, lat1, lon1, &brg, NULL);
2886 toSM(lat1, lon1, 0., lon, &x1, &y1);
2887 toSM(lat1, lon, 0., lon, &x, &y1);
2892 }
else if (brg >= 180.) {
2895 }
else if (brg >= 90.) {
2902 double ydelta = fabs(x1) * tan(brg * DEGREE);
2904 double crosslat, crosslon;
2905 fromSM(x, y1 + dir * ydelta, 0., lon, &crosslat, &crosslon);