1 // sg_geodesy.cxx -- routines to convert between geodetic and geocentric
4 // Copied and adapted directly from LaRCsim/ls_geodesy.c
6 // See below for the complete original LaRCsim comments.
10 #include <simgear/compiler.h>
12 #ifdef SG_HAVE_STD_INCLUDES
22 #include <simgear/constants.h>
23 #include <simgear/debug/logstream.hxx>
25 #include "point3d.hxx"
26 #include "sg_geodesy.hxx"
27 #include "localconsts.hxx"
30 #ifndef SG_HAVE_NATIVE_SGI_COMPILERS
35 #define DOMAIN_ERR_DEBUG 1
38 // sgGeocToGeod(lat_geoc, radius, *lat_geod, *alt, *sea_level_r)
40 // lat_geoc Geocentric latitude, radians, + = North
41 // radius C.G. radius to earth center (meters)
44 // lat_geod Geodetic latitude, radians, + = North
45 // alt C.G. altitude above mean sea level (meters)
46 // sea_level_r radius from earth center to sea level at
47 // local vertical (surface normal) of C.G. (meters)
50 void sgGeocToGeod( double lat_geoc, double radius, double
51 *lat_geod, double *alt, double *sea_level_r )
53 #ifdef DOMAIN_ERR_DEBUG
54 errno = 0; // start with error zero'd
57 double t_lat, x_alpha, mu_alpha, delt_mu, r_alpha, l_point, rho_alpha;
58 double sin_mu_a, denom,delt_lambda, lambda_sl, sin_lambda_sl;
60 if( ( (SGD_PI_2 - lat_geoc) < SG_ONE_SECOND ) // near North pole
61 || ( (SGD_PI_2 + lat_geoc) < SG_ONE_SECOND ) ) // near South pole
64 *sea_level_r = SG_EQUATORIAL_RADIUS_M*E;
65 *alt = radius - *sea_level_r;
67 // cout << " lat_geoc = " << lat_geoc << endl;
68 t_lat = tan(lat_geoc);
69 // cout << " tan(t_lat) = " << t_lat << endl;
70 x_alpha = E*SG_EQUATORIAL_RADIUS_M/sqrt(t_lat*t_lat + E*E);
71 #ifdef DOMAIN_ERR_DEBUG
73 perror("fgGeocToGeod()");
74 SG_LOG( SG_GENERAL, SG_ALERT, "sqrt(" << t_lat*t_lat + E*E << ")" );
77 // cout << " x_alpha = " << x_alpha << endl;
78 double tmp = sqrt(SG_EQ_RAD_SQUARE_M - x_alpha * x_alpha);
79 if ( tmp < 0.0 ) { tmp = 0.0; }
80 #ifdef DOMAIN_ERR_DEBUG
82 perror("fgGeocToGeod()");
83 SG_LOG( SG_GENERAL, SG_ALERT, "sqrt(" << SG_EQ_RAD_SQUARE_M - x_alpha * x_alpha
87 mu_alpha = atan2(tmp,E*x_alpha);
88 if (lat_geoc < 0) mu_alpha = - mu_alpha;
89 sin_mu_a = sin(mu_alpha);
90 delt_lambda = mu_alpha - lat_geoc;
91 r_alpha = x_alpha/cos(lat_geoc);
92 l_point = radius - r_alpha;
93 *alt = l_point*cos(delt_lambda);
95 denom = sqrt(1-EPS*EPS*sin_mu_a*sin_mu_a);
96 #ifdef DOMAIN_ERR_DEBUG
98 perror("fgGeocToGeod()");
99 SG_LOG( SG_GENERAL, SG_ALERT, "sqrt(" <<
100 1-EPS*EPS*sin_mu_a*sin_mu_a << ")" );
103 rho_alpha = SG_EQUATORIAL_RADIUS_M*(1-EPS)/
105 delt_mu = atan2(l_point*sin(delt_lambda),rho_alpha + *alt);
106 *lat_geod = mu_alpha - delt_mu;
107 lambda_sl = atan( E*E * tan(*lat_geod) ); // SL geoc. latitude
108 sin_lambda_sl = sin( lambda_sl );
110 sqrt(SG_EQ_RAD_SQUARE_M / (1 + ((1/(E*E))-1)*sin_lambda_sl*sin_lambda_sl));
111 #ifdef DOMAIN_ERR_DEBUG
113 perror("fgGeocToGeod()");
114 SG_LOG( SG_GENERAL, SG_ALERT, "sqrt(" <<
115 SG_EQ_RAD_SQUARE_M / (1 + ((1/(E*E))-1)*sin_lambda_sl*sin_lambda_sl)
125 // sgGeodToGeoc( lat_geod, alt, *sl_radius, *lat_geoc )
127 // lat_geod Geodetic latitude, radians, + = North
128 // alt C.G. altitude above mean sea level (meters)
131 // sl_radius SEA LEVEL radius to earth center (meters)
132 // (add Altitude to get true distance from earth center.
133 // lat_geoc Geocentric latitude, radians, + = North
137 void sgGeodToGeoc( double lat_geod, double alt, double *sl_radius,
140 double lambda_sl, sin_lambda_sl, cos_lambda_sl, sin_mu, cos_mu, px, py;
142 #ifdef DOMAIN_ERR_DEBUG
146 lambda_sl = atan( E*E * tan(lat_geod) ); // sea level geocentric latitude
147 sin_lambda_sl = sin( lambda_sl );
148 cos_lambda_sl = cos( lambda_sl );
149 sin_mu = sin(lat_geod); // Geodetic (map makers') latitude
150 cos_mu = cos(lat_geod);
152 sqrt(SG_EQ_RAD_SQUARE_M / (1 + ((1/(E*E))-1)*sin_lambda_sl*sin_lambda_sl));
153 #ifdef DOMAIN_ERR_DEBUG
155 perror("fgGeodToGeoc()");
156 SG_LOG( SG_GENERAL, SG_ALERT, "sqrt(" <<
157 SG_EQ_RAD_SQUARE_M / (1 + ((1/(E*E))-1)*sin_lambda_sl*sin_lambda_sl)
161 py = *sl_radius*sin_lambda_sl + alt*sin_mu;
162 px = *sl_radius*cos_lambda_sl + alt*cos_mu;
163 *lat_geoc = atan2( py, px );
167 // Direct and inverse distance functions
169 // Proceedings of the 7th International Symposium on Geodetic
170 // Computations, 1985
172 // "The Nested Coefficient Method for Accurate Solutions of Direct and
173 // Inverse Geodetic Problems With Any Length"
178 // modified for FlightGear to use WGS84 only -- Norman Vine
180 #define GEOD_INV_PI SGD_PI
185 // for WGS_84 a = 6378137.000, rf = 298.257223563;
187 static double M0( double e2 ) {
189 return GEOD_INV_PI*(1.0 - e2*( 1.0/4.0 + e2*( 3.0/64.0 +
190 e2*(5.0/256.0) )))/2.0;
194 // given, alt, lat1, lon1, az1 and distance (s), calculate lat2, lon2
195 // and az2. Lat, lon, and azimuth are in degrees. distance in meters
196 int geo_direct_wgs_84 ( double alt, double lat1, double lon1, double az1,
197 double s, double *lat2, double *lon2, double *az2 )
199 double a = 6378137.000, rf = 298.257223563;
200 double RADDEG = (GEOD_INV_PI)/180.0, testv = 1.0E-10;
201 double f = ( rf > 0.0 ? 1.0/rf : 0.0 );
202 double b = a*(1.0-f);
203 double e2 = f*(2.0-f);
204 double phi1 = lat1*RADDEG, lam1 = lon1*RADDEG;
205 double sinphi1 = sin(phi1), cosphi1 = cos(phi1);
206 double azm1 = az1*RADDEG;
207 double sinaz1 = sin(azm1), cosaz1 = cos(azm1);
210 if( fabs(s) < 0.01 ) { // distance < centimeter => congruency
214 if( *az2 > 360.0 ) *az2 -= 360.0;
216 } else if( cosphi1 ) { // non-polar origin
217 // u1 is reduced latitude
218 double tanu1 = sqrt(1.0-e2)*sinphi1/cosphi1;
219 double sig1 = atan2(tanu1,cosaz1);
220 double cosu1 = 1.0/sqrt( 1.0 + tanu1*tanu1 ), sinu1 = tanu1*cosu1;
221 double sinaz = cosu1*sinaz1, cos2saz = 1.0-sinaz*sinaz;
222 double us = cos2saz*e2/(1.0-e2);
225 double ta = 1.0+us*(4096.0+us*(-768.0+us*(320.0-175.0*us)))/16384.0,
226 tb = us*(256.0+us*(-128.0+us*(74.0-47.0*us)))/1024.0,
229 // FIRST ESTIMATE OF SIGMA (SIG)
230 double first = s/(b*ta); // !!
232 double c2sigm, sinsig,cossig, temp,denom,rnumer, dlams, dlam;
234 c2sigm = cos(2.0*sig1+sig);
235 sinsig = sin(sig); cossig = cos(sig);
238 tb*sinsig*(c2sigm+tb*(cossig*(-1.0+2.0*c2sigm*c2sigm) -
239 tb*c2sigm*(-3.0+4.0*sinsig*sinsig)
240 *(-3.0+4.0*c2sigm*c2sigm)/6.0)
242 } while( fabs(sig-temp) > testv);
244 // LATITUDE OF POINT 2
245 // DENOMINATOR IN 2 PARTS (TEMP ALSO USED LATER)
246 temp = sinu1*sinsig-cosu1*cossig*cosaz1;
247 denom = (1.0-f)*sqrt(sinaz*sinaz+temp*temp);
250 rnumer = sinu1*cossig+cosu1*sinsig*cosaz1;
251 *lat2 = atan2(rnumer,denom)/RADDEG;
253 // DIFFERENCE IN LONGITUDE ON AUXILARY SPHERE (DLAMS )
254 rnumer = sinsig*sinaz1;
255 denom = cosu1*cossig-sinu1*sinsig*cosaz1;
256 dlams = atan2(rnumer,denom);
259 tc = f*cos2saz*(4.0+f*(4.0-3.0*cos2saz))/16.0;
261 // DIFFERENCE IN LONGITUDE
262 dlam = dlams-(1.0-tc)*f*sinaz*(sig+tc*sinsig*
266 *lon2 = (lam1+dlam)/RADDEG;
267 if (*lon2 > 180.0 ) *lon2 -= 360.0;
268 if (*lon2 < -180.0 ) *lon2 += 360.0;
270 // AZIMUTH - FROM NORTH
271 *az2 = atan2(-sinaz,temp)/RADDEG;
272 if ( fabs(*az2) < testv ) *az2 = 0.0;
273 if( *az2 < 0.0) *az2 += 360.0;
275 } else { // phi1 == 90 degrees, polar origin
276 double dM = a*M0(e2) - s;
277 double paz = ( phi1 < 0.0 ? 180.0 : 0.0 );
278 return geo_direct_wgs_84( alt, 0.0, lon1, paz, dM,lat2,lon2,az2 );
283 // given alt, lat1, lon1, lat2, lon2, calculate starting and ending
284 // az1, az2 and distance (s). Lat, lon, and azimuth are in degrees.
285 // distance in meters
286 int geo_inverse_wgs_84( double alt, double lat1, double lon1, double lat2,
287 double lon2, double *az1, double *az2, double *s )
289 double a = 6378137.000, rf = 298.257223563;
291 double RADDEG = (GEOD_INV_PI)/180.0, testv = 1.0E-10;
292 double f = ( rf > 0.0 ? 1.0/rf : 0.0 );
293 double b = a*(1.0-f);
294 // double e2 = f*(2.0-f); // unused in this routine
295 double phi1 = lat1*RADDEG, lam1 = lon1*RADDEG;
296 double sinphi1 = sin(phi1), cosphi1 = cos(phi1);
297 double phi2 = lat2*RADDEG, lam2 = lon2*RADDEG;
298 double sinphi2 = sin(phi2), cosphi2 = cos(phi2);
300 if( (fabs(lat1-lat2) < testv &&
301 ( fabs(lon1-lon2) < testv) || fabs(lat1-90.0) < testv ) )
303 // TWO STATIONS ARE IDENTICAL : SET DISTANCE & AZIMUTHS TO ZERO */
304 *az1 = 0.0; *az2 = 0.0; *s = 0.0;
306 } else if( fabs(cosphi1) < testv ) {
307 // initial point is polar
308 int k = geo_inverse_wgs_84( alt, lat2,lon2,lat1,lon1, az1,az2,s );
309 k = k; // avoid compiler error since return result is unused
310 b = *az1; *az1 = *az2; *az2 = b;
312 } else if( fabs(cosphi2) < testv ) {
313 // terminal point is polar
314 int k = geo_inverse_wgs_84( alt, lat1,lon1,lat1,lon1+180.0,
316 k = k; // avoid compiler error since return result is unused
319 if( *az2 > 360.0 ) *az2 -= 360.0;
321 } else if( (fabs( fabs(lon1-lon2) - 180 ) < testv) &&
322 (fabs(lat1+lat2) < testv) )
324 // Geodesic passes through the pole (antipodal)
326 geo_inverse_wgs_84( alt, lat1,lon1, lat1,lon2, az1,az2, &s1 );
327 geo_inverse_wgs_84( alt, lat2,lon2, lat1,lon2, az1,az2, &s2 );
332 // antipodal and polar points don't get here
333 double dlam = lam2 - lam1, dlams = dlam;
334 double sdlams,cdlams, sig,sinsig,cossig, sinaz,
336 double tc,temp, us,rnumer,denom, ta,tb;
337 double cosu1,sinu1, sinu2,cosu2;
340 temp = (1.0-f)*sinphi1/cosphi1;
341 cosu1 = 1.0/sqrt(1.0+temp*temp);
343 temp = (1.0-f)*sinphi2/cosphi2;
344 cosu2 = 1.0/sqrt(1.0+temp*temp);
348 sdlams = sin(dlams), cdlams = cos(dlams);
349 sinsig = sqrt(cosu2*cosu2*sdlams*sdlams+
350 (cosu1*sinu2-sinu1*cosu2*cdlams)*
351 (cosu1*sinu2-sinu1*cosu2*cdlams));
352 cossig = sinu1*sinu2+cosu1*cosu2*cdlams;
354 sig = atan2(sinsig,cossig);
355 sinaz = cosu1*cosu2*sdlams/sinsig;
356 cos2saz = 1.0-sinaz*sinaz;
357 c2sigm = (sinu1 == 0.0 || sinu2 == 0.0 ? cossig :
358 cossig-2.0*sinu1*sinu2/cos2saz);
359 tc = f*cos2saz*(4.0+f*(4.0-3.0*cos2saz))/16.0;
361 dlams = dlam+(1.0-tc)*f*sinaz*
363 (c2sigm+tc*cossig*(-1.0+2.0*c2sigm*c2sigm)));
364 if (fabs(dlams) > GEOD_INV_PI && iter++ > 50) {
367 } while ( fabs(temp-dlams) > testv);
369 us = cos2saz*(a*a-b*b)/(b*b); // !!
370 // BACK AZIMUTH FROM NORTH
371 rnumer = -(cosu1*sdlams);
372 denom = sinu1*cosu2-cosu1*sinu2*cdlams;
373 *az2 = atan2(rnumer,denom)/RADDEG;
374 if( fabs(*az2) < testv ) *az2 = 0.0;
375 if(*az2 < 0.0) *az2 += 360.0;
377 // FORWARD AZIMUTH FROM NORTH
378 rnumer = cosu2*sdlams;
379 denom = cosu1*sinu2-sinu1*cosu2*cdlams;
380 *az1 = atan2(rnumer,denom)/RADDEG;
381 if( fabs(*az1) < testv ) *az1 = 0.0;
382 if(*az1 < 0.0) *az1 += 360.0;
385 ta = 1.0+us*(4096.0+us*(-768.0+us*(320.0-175.0*us)))/
387 tb = us*(256.0+us*(-128.0+us*(74.0-47.0*us)))/1024.0;
390 *s = b*ta*(sig-tb*sinsig*
391 (c2sigm+tb*(cossig*(-1.0+2.0*c2sigm*c2sigm)-tb*
392 c2sigm*(-3.0+4.0*sinsig*sinsig)*
393 (-3.0+4.0*c2sigm*c2sigm)/6.0)/
400 /***************************************************************************
404 ----------------------------------------------------------------------------
406 FUNCTION: Converts geocentric coordinates to geodetic positions
408 ----------------------------------------------------------------------------
410 MODULE STATUS: developmental
412 ----------------------------------------------------------------------------
414 GENEALOGY: Written as part of LaRCSim project by E. B. Jackson
416 ----------------------------------------------------------------------------
418 DESIGNED BY: E. B. Jackson
420 CODED BY: E. B. Jackson
422 MAINTAINED BY: E. B. Jackson
424 ----------------------------------------------------------------------------
426 MODIFICATION HISTORY:
430 930208 Modified to avoid singularity near polar region. EBJ
431 930602 Moved backwards calcs here from ls_step. EBJ
432 931214 Changed erroneous Latitude and Altitude variables to
433 *lat_geod and *alt in routine ls_geoc_to_geod. EBJ
434 940111 Changed header files from old ls_eom.h style to ls_types,
435 and ls_constants. Also replaced old DATA type with new
441 * Revision 1.5 1994/01/11 18:47:05 bjax
442 * Changed include files to use types and constants, not ls_eom.h
443 * Also changed DATA type to SCALAR type.
445 * Revision 1.4 1993/12/14 21:06:47 bjax
446 * Removed global variable references Altitude and Latitude. EBJ
448 * Revision 1.3 1993/06/02 15:03:40 bjax
449 * Made new subroutine for calculating geodetic to geocentric; changed name
450 * of forward conversion routine from ls_geodesy to ls_geoc_to_geod.
453 ----------------------------------------------------------------------------
457 [ 1] Stevens, Brian L.; and Lewis, Frank L.: "Aircraft
458 Control and Simulation", Wiley and Sons, 1992.
462 ----------------------------------------------------------------------------
466 ----------------------------------------------------------------------------
470 ----------------------------------------------------------------------------
473 lat_geoc Geocentric latitude, radians, + = North
474 radius C.G. radius to earth center, ft
476 ----------------------------------------------------------------------------
479 lat_geod Geodetic latitude, radians, + = North
480 alt C.G. altitude above mean sea level, ft
481 sea_level_r radius from earth center to sea level at
482 local vertical (surface normal) of C.G.
484 --------------------------------------------------------------------------*/