int gctpc_get_latlon(unsigned char **sec, double **lon, double **lat) { int gdt; unsigned char *gds; double r_maj; /* major axis */ double r_min; /* minor axis */ double lat1; /* first standard parallel */ double lat2; /* second standard parallel */ double c_lon; /* center longitude */ double c_lat; /* center latitude */ double false_east; /* x offset in meters */ double false_north; double dx, dy; double x0, y0; long int (*inv_fn)(); double *llat, *llon, rlon, rlat; int i, nnx, nny, nres, nscan; unsigned int nnpnts; long long_i; gdt = code_table_3_1(sec); gds = sec[3]; /* only process certain grids */ if (gdt != 10 && gdt != 20 && gdt != 30 && gdt != 31) return 1; get_nxny(sec, &nnx, &nny, &nnpnts, &nres, &nscan); /* potentially staggered */ if (nnx < 1 || nny < 1) return 1; llat = *lat; llon = *lon; if (llat != NULL) { free(llat); free(llon); *lat = *lon = llat = llon = NULL; } inv_fn = NULL; dx = dy = 0.0; if (gdt == 10) { // mercator /* get earth axis */ axes_earth(sec, &r_maj, &r_min); dy = GDS_Mercator_dy(gds); dx = GDS_Mercator_dx(gds); /* central point */ c_lon = GDS_Mercator_ori_angle(gds) * (M_PI/180.0); c_lat = GDS_Mercator_latD(gds) * (M_PI/180.0); /* find the eastling and northing of of the 1st grid point */ false_east = false_north = 0.0; long_i = merforint(r_maj,r_min,c_lon,c_lat,false_east,false_north); rlon = GDS_Mercator_lon1(gds) * (M_PI/180.0); rlat = GDS_Mercator_lat1(gds) * (M_PI/180.0); long_i = merfor(rlon, rlat, &x0, &y0); /* initialize for 1st grid point */ x0 = -x0; y0 = -y0; long_i = merinvint(r_maj,r_min,c_lon,c_lat,x0,y0); inv_fn = &merinv; } else if (gdt == 20) { // polar stereographic /* get earth axis */ axes_earth(sec, &r_maj, &r_min); dy = GDS_Polar_dy(gds); dx = GDS_Polar_dx(gds); /* central point */ c_lon = GDS_Polar_lov(gds) * (M_PI/180.0); c_lat = GDS_Polar_lad(gds) * (M_PI/180.0); /* find the eastling and northing of of the 1st grid point */ false_east = false_north = 0.0; long_i = psforint(r_maj,r_min,c_lon,c_lat,false_east,false_north); rlon = GDS_Polar_lon1(gds) * (M_PI/180.0); rlat = GDS_Polar_lat1(gds) * (M_PI/180.0); long_i = psfor(rlon, rlat, &x0, &y0); /* initialize for 1st grid point */ x0 = -x0; y0 = -y0; long_i = psinvint(r_maj,r_min,c_lon,c_lat,x0,y0); inv_fn = &psinv; } else if (gdt == 30) { // lambert conformal conic /* get earth axis */ axes_earth(sec, &r_maj, &r_min); dy = GDS_Lambert_dy(gds); dx = GDS_Lambert_dx(gds); //printf(">>> gctpc dx %lf, dy %lf\n", dx, dy); /* latitudes of tangent/intersection */ lat1 = GDS_Lambert_Latin1(gds) * (M_PI/180.0); lat2 = GDS_Lambert_Latin2(gds) * (M_PI/180.0); /* central point */ c_lon = GDS_Lambert_Lov(gds) * (M_PI/180.0); c_lat = GDS_Lambert_LatD(gds) * (M_PI/180.0); /* find the eastling and northing of of the 1st grid point */ false_east = false_north = 0.0; long_i = lamccforint(r_maj,r_min,lat1,lat2,c_lon,c_lat,false_east,false_north); rlon = GDS_Lambert_Lo1(gds) * (M_PI/180.0); rlat = GDS_Lambert_La1(gds) * (M_PI/180.0); long_i = lamccfor(rlon, rlat, &x0, &y0); /* initialize for 1st grid point */ x0 = -x0; y0 = -y0; long_i = lamccinvint(r_maj,r_min,lat1,lat2,c_lon,c_lat,x0,y0); inv_fn = &lamccinv; } else if (gdt == 31) { // albers equal area /* get earth axis */ axes_earth(sec, &r_maj, &r_min); dy = GDS_Albers_dy(gds); dx = GDS_Albers_dx(gds); /* latitudes of tangent/intersection */ lat1 = GDS_Albers_Latin1(gds) * (M_PI/180.0); lat2 = GDS_Albers_Latin2(gds) * (M_PI/180.0); /* central point */ c_lon = GDS_Albers_Lov(gds) * (M_PI/180.0); c_lat = GDS_Albers_LatD(gds) * (M_PI/180.0); /* find the eastling and northing of of the 1st grid point */ false_east = false_north = 0.0; long_i = alberforint(r_maj,r_min,lat1,lat2,c_lon,c_lat,false_east,false_north); rlon = GDS_Albers_Lo1(gds) * (M_PI/180.0); rlat = GDS_Albers_La1(gds) * (M_PI/180.0); long_i = alberfor(rlon, rlat, &x0, &y0); /* initialize for 1st grid point */ x0 = -x0; y0 = -y0; long_i = alberinvint(r_maj,r_min,lat1,lat2,c_lon,c_lat,x0,y0); inv_fn = &alberinv; } if (inv_fn == NULL) return 1; if ((*lat = llat = (double *) malloc(nnpnts * sizeof(double))) == NULL) { fatal_error("gctpc_get_latlon memory allocation failed",""); } if ((*lon = llon = (double *) malloc(nnpnts * sizeof(double))) == NULL) { fatal_error("gctpc_get_latlon memory allocation failed",""); } /* put x[] and y[] values in lon and lat */ if (stagger(sec, nnpnts, llon, llat)) fatal_error("gctpc: stagger problem",""); printf(">> stagger gctpc x00 %lf y00 %lf\n",llon[0], llat[0]); #pragma omp parallel for schedule(static) private(i) for (i = 0; i < nnpnts; i++) { inv_fn(llon[i]*dx, llat[i]*dy, llon+i, llat+i); llat[i] *= (180.0 / M_PI); llon[i] *= (180.0 / M_PI); if (llon[i] < 0.0) llon[i] += 360.0; } return 0; }
int proj4_initialize(unsigned char **sec, struct proj4_struct *projection) { unsigned int gdt; unsigned char *gds; int has_np; double r_maj, r_min, c_lon, c_lat, x_0, y_0, lat1, lon1; double latsp1, latsp2; char proj4_def[1000]; gdt = code_table_3_1(sec); gds = sec[3]; axes_earth(sec, &r_maj, &r_min); projection->radius_minor = r_min; projection->radius_minor = r_min; if (gdt == 0) { projection->proj_is_nop = 1; projection->x_0 = projection->y_0 = projection->lat_0 = projection->lon_0 = 0.0; return 0; } projection->proj_is_nop = 0; sprintf(proj4_def,"+proj=latlong +a=%lf +b=%lf",r_maj, r_min); if ( (projection->pj_latlon = pj_init_plus(proj4_def)) == NULL) fatal_error("proj4_initialize: pj_init_plus %s failed", proj4_def); if (gdt == 10 && (GDS_Mercator_ori_angle(gds) == 0.0) ) { // mercator no rotation /* central point */ c_lon = 0.0; c_lat = GDS_Mercator_latD(gds); sprintf(proj4_def,"+proj=merc +lat_ts=%lf +lat_0=0 +lon_0=0 +x_0=0 +y_0=0 +a=%lf +b=%lf", c_lat, r_maj, r_min); if ((projection->pj_grid = pj_init_plus(proj4_def)) == NULL) fatal_error("proj4_initialize: pj_init_plus %s failed", proj4_def); /* longitude, latitude of first grid point */ projection->lat_0 = lat1 = GDS_Mercator_lat1(gds); projection->lon_0 = lon1 = GDS_Mercator_lon1(gds); x_0 = lon1 * DEG_TO_RAD; y_0 = lat1 * DEG_TO_RAD; if ( pj_transform(projection->pj_latlon, projection->pj_grid, 1, 1, &x_0, &y_0, NULL) != 0 ) fatal_error("proj4_initialize: Proj4 transform to lat-lon",""); projection->x_0 = x_0; projection->y_0 = y_0; } else if (gdt == 20) { // polar stereographic /* central point */ c_lon = GDS_Polar_lov(gds); c_lat = GDS_Polar_lad(gds); /* strange but np/sp flag is used by proj4 but not gctpc */ has_np = ((flag_table_3_5(sec) & 128) == 0); sprintf(proj4_def,"+proj=stere +lat_ts=%lf +lat_0=%s +lon_0=%lf +k_0=1 +x_0=0 +y_0=0 +a=%lf +b=%lf", c_lat, has_np ? "90" : "-90", c_lon, r_maj,r_min); if ((projection->pj_grid = pj_init_plus(proj4_def)) == NULL) fatal_error("Proj4: pj_init_plus %s failed", proj4_def); /* longitude, latitude of first grid point */ projection->lon_0 = lon1 = GDS_Polar_lon1(gds); projection->lat_0 = lat1 = GDS_Polar_lat1(gds); x_0 = lon1 * DEG_TO_RAD; y_0 = lat1 * DEG_TO_RAD; if ( pj_transform(projection->pj_latlon, projection->pj_grid, 1, 1, &x_0, &y_0, NULL) != 0 ) fatal_error("proj4_init: Proj4 transform to lat-lon",""); projection->x_0 = x_0; projection->y_0 = y_0; } else if (gdt == 30) { // lambert conformal conic /* latitudes of tangent/intersection */ latsp1 = GDS_Lambert_Latin1(gds); latsp2 = GDS_Lambert_Latin2(gds); /* central point */ c_lon = GDS_Lambert_Lov(gds); c_lat = GDS_Lambert_LatD(gds); sprintf(proj4_def,"+proj=lcc +lon_0=%lf +lat_0=%lf +lat_1=%lf +lat_2=%lf +a=%lf +b=%lf",c_lon, c_lat,latsp1,latsp2,r_maj,r_min); if ((projection->pj_grid = pj_init_plus(proj4_def)) == NULL) fatal_error("Proj4: pj_init_plus %s failed", proj4_def); /* longitude, latitude of first grid point */ lon1 = GDS_Lambert_Lo1(gds); lat1 = GDS_Lambert_La1(gds); x_0 = lon1 * DEG_TO_RAD; y_0 = lat1 * DEG_TO_RAD; /* longitude, latitude of first grid point */ projection->lon_0 = lon1 = GDS_Lambert_Lo1(gds); projection->lat_0 = lat1 = GDS_Lambert_La1(gds); x_0 = lon1 * DEG_TO_RAD; y_0 = lat1 * DEG_TO_RAD; if ( pj_transform(projection->pj_latlon, projection->pj_grid, 1, 1, &x_0, &y_0, NULL) != 0 ) fatal_error("proj4_init: Proj4 transform to lat-lon",""); projection->x_0 = x_0; projection->y_0 = y_0; } else { return 1; } return 0; }
int lambert2ll(unsigned char **sec, double **llat, double **llon) { double n; double *lat, *lon; double dx, dy, lat1r, lon1r, lon2d, lon2r, latin1r, latin2r; double lond, latd, d_lon; double f, rho, rhoref, theta, startx, starty; int j, nnx, nny, nres, nscan; double x, y, tmp; unsigned char *gds; double latDr; double earth_radius; unsigned int nnpnts; get_nxny(sec, &nnx, &nny, &nnpnts, &nres, &nscan); if (nnx <= 0 || nny <= 0) { fprintf(stderr,"Sorry code does not handle variable nx/ny yet\n"); return 0; } earth_radius = radius_earth(sec); gds = sec[3]; dy = GDS_Lambert_dy(gds); dx = GDS_Lambert_dx(gds); lat1r = GDS_Lambert_La1(gds) * (M_PI / 180.0); lon1r = GDS_Lambert_Lo1(gds) * (M_PI / 180.0); lon2d = GDS_Lambert_Lov(gds); lon2r = lon2d * (M_PI / 180.0); latin1r = GDS_Lambert_Latin1(gds) * (M_PI/180.0); latin2r = GDS_Lambert_Latin2(gds) * (M_PI/180.0); // fix for theta start value crossing 0 longitude // if ((lon1r - lon2r) > 0) lon2r = lon2r + 2*M_PI; // // Latitude of "false origin" where scales are defined. // It is used to estimate "reference_R", rhoref. // Often latDr == latin1r == latin2r and non-modified code is true and works fine. // But could be different if intersection latitudes latin1r and latin2r are different. // Usually latDr must be latin1r <= latDr <= latin2r, other could be strange. // latDr = GDS_Lambert_LatD(gds) * (M_PI/180.0); if (lon1r < 0) fatal_error("bad GDS, lon1r < 0.0",""); if ( fabs(latin1r - latin2r) < 1E-09 ) { n = sin(latin1r); } else { n = log(cos(latin1r)/cos(latin2r)) / log(tan(M_PI_4 + latin2r/2.0) / tan(M_PI_4 + latin1r/2.0)); } f = (cos(latin1r) * pow(tan(M_PI_4 + latin1r/2.0), n)) / n; rho = earth_radius * f * pow(tan(M_PI_4 + lat1r/2.0),-n); // old rhoref = earth_radius * f * pow(tan(M_PI_4 + latin1r/2.0),-n); rhoref = earth_radius * f * pow(tan(M_PI_4 + latDr/2.0),-n); // 2/2009 .. new code d_lon = lon1r - lon2r; if (d_lon > M_PI) d_lon -= 2*M_PI; if (d_lon < -M_PI) d_lon += 2*M_PI; theta = n * d_lon; // 2/2009 theta = n * (lon1r - lon2r); startx = rho * sin(theta); starty = rhoref - rho * cos(theta); if ((*llat = (double *) malloc(nnpnts * sizeof(double))) == NULL) { fatal_error("lambert2ll memory allocation failed",""); } if ((*llon = (double *) malloc(nnpnts * sizeof(double))) == NULL) { fatal_error("lambert2ll memory allocation failed",""); } lat = *llat; lon = *llon; /* put x[] and y[] values in lon[] and lat[] */ if (stagger(sec, nnpnts, lon, lat)) fatal_error("geo: stagger problem",""); dx = fabs(dx); dy = fabs(dy); #pragma omp parallel for private(j,x,y,tmp,theta,rho,lond,latd) for (j = 0; j < nnpnts; j++) { y = starty + lat[j]*dy; x = startx + lon[j]*dx; tmp = rhoref - y; theta = atan(x / tmp); rho = sqrt(x * x + tmp*tmp); rho = n > 0 ? rho : -rho; lond = lon2d + todegrees(theta/n); latd = todegrees(2.0 * atan(pow(earth_radius * f/rho,1.0/n)) - M_PI_2); lond = lond >= 360.0 ? lond - 360.0 : lond; lond = lond < 0.0 ? lond + 360.0 : lond; lon[j] = lond; lat[j] = latd; } return 0; } /* end lambert2ll() */