int cyclic(unsigned char **sec) { int grid_template, nx, ny, res, scan, flag_3_3, no_dx, basic_ang, sub_ang; unsigned int npnts; unsigned char *gds; double dlon, units, lon1, lon2; get_nxny(sec, &nx, &ny, &npnts, &res, &scan); if (GDS_Scan_staggered(scan)) return 0; if (nx <= 1 || ny <= 0) return 0; grid_template = code_table_3_1(sec); gds = sec[3]; flag_3_3 = flag_table_3_3(sec); no_dx = 0; if (flag_3_3 != -1) { if ((flag_3_3 & 0x20) == 0) no_dx = 1; } if (no_dx) return 0; if (grid_template == 0) { basic_ang = GDS_LatLon_basic_ang(gds); sub_ang = GDS_LatLon_sub_ang(gds); units = basic_ang == 0 ? 0.000001 : (double) basic_ang / (double) sub_ang; /* dlon has to be defined */ dlon = units * GDS_LatLon_dlon(gds); return (fabs(nx*dlon-360.0) < ERROR); } if (grid_template == 10) { if (output_order != wesn) return 0; // only works with we:sn order lon1 = GDS_Mercator_lon1(gds); lon2 = GDS_Mercator_lon2(gds); if (lon2 < lon1) lon2 += 360.0; dlon = (lon2-lon1)*nx/(nx-1.0); return (fabs(dlon-360.0) < ERROR); } if (grid_template == 40) { basic_ang = GDS_Gaussian_basic_ang(gds); sub_ang = GDS_Gaussian_sub_ang(gds); units = basic_ang == 0 ? 0.000001 : (double) basic_ang / (double) sub_ang; /* dlon has to be defined */ dlon = units * GDS_Gaussian_dlon(gds); return (fabs(nx*dlon-360.0) < ERROR); } return 0; }
int small_grib(unsigned char **sec, int mode, float *data, double *lon, double *lat, unsigned int ndata, int ix0, int ix1, int iy0, int iy1, FILE *out) { int can_subset, grid_template; int nx, ny, res, scan, new_nx, new_ny, i, j; unsigned int sec3_len, new_ndata, k, npnts; unsigned char *sec3, *new_sec[9]; double units; int basic_ang, sub_ang, cyclic_grid; float *new_data; get_nxny(sec, &nx, &ny, &npnts, &res, &scan); /* get nx, ny, and scan mode of grid */ grid_template = code_table_3_1(sec); // make a copy of the gds (sec3) sec3_len = GB2_Sec3_size(sec); sec3 = (unsigned char *) malloc(sec3_len); for (k = 0; k < sec3_len; k++) sec3[k] = sec[3][k]; // make a copy of the sec[] with new sec3 new_sec[0] = sec[0]; new_sec[1] = sec[1]; new_sec[2] = sec[2]; new_sec[3] = sec3; new_sec[4] = sec[4]; new_sec[5] = sec[5]; new_sec[6] = sec[6]; new_sec[7] = sec[7]; // new_sec[8] = sec[8]; not needed by writing routines can_subset = 1; if (lat == NULL || lon == NULL) can_subset = 0; new_nx = ix1-ix0+1; new_ny = iy1-iy0+1; if (new_nx <= 0) fatal_error("small_grib, new_nx is <= 0",""); if (new_ny <= 0) fatal_error("small_grib, new_ny is <= 0",""); new_ndata = new_nx * new_ny; cyclic_grid = 0; if (can_subset) { cyclic_grid = cyclic(sec); // lat-lon grid - no thinning if ((grid_template == 0 && sec3_len == 72) || (grid_template == 1 && sec3_len == 04)) { uint_char(new_nx,sec3+30); // nx uint_char(new_ny,sec3+34); // ny basic_ang = GDS_LatLon_basic_ang(sec3); sub_ang = GDS_LatLon_sub_ang(sec3); if (basic_ang != 0) { units = (double) basic_ang / (double) sub_ang; } else { units = 0.000001; } i = lat[ idx(ix0,iy0,nx,ny,cyclic_grid) ] / units; // lat1 int_char(i,sec3+46); i = lon[ idx(ix0,iy0,nx,ny,cyclic_grid) ] / units; // lon1 int_char(i,sec3+50); i = lat[ idx(ix1,iy1,nx,ny,cyclic_grid) ] / units; // lat2 int_char(i,sec3+55); i = lon[ idx(ix1,iy1,nx,ny,cyclic_grid) ] / units; // lon2 int_char(i,sec3+59); } else if ((grid_template == 40 && sec3_len == 72)) { // full Gaussian grid uint_char(new_nx,sec3+30); // nx uint_char(new_ny,sec3+34); // ny basic_ang = GDS_Gaussian_basic_ang(sec3); sub_ang = GDS_Gaussian_sub_ang(sec3); if (basic_ang != 0) { units = (double) basic_ang / (double) sub_ang; } else { units = 0.000001; } i = lat[ idx(ix0,iy0,nx,ny,cyclic_grid) ] / units; // lat1 int_char(i,sec3+46); i = lon[ idx(ix0,iy0,nx,ny,cyclic_grid) ] / units; // lon1 int_char(i,sec3+50); i = lat[ idx(ix1,iy1,nx,ny,cyclic_grid) ] / units; // lat2 int_char(i,sec3+55); i = lon[ idx(ix1,iy1,nx,ny,cyclic_grid) ] / units; // lon2 int_char(i,sec3+59); } // polar-stereo graphic, lambert conformal , no thinning else if ((grid_template == 20 && sec3_len == 65) || // polar stereographic (grid_template == 30 && sec3_len == 81)) { // lambert conformal uint_char(new_nx,sec3+30); // nx uint_char(new_ny,sec3+34); // ny i = (int) (lat[ idx(ix0,iy0,nx,ny,cyclic_grid) ] * 1000000.0); // lat1 int_char(i,sec3+38); i = (int) (lon[ idx(ix0,iy0,nx,ny,cyclic_grid) ] * 1000000.0); // lon1 int_char(i,sec3+42); } // mercator, no thinning else if (grid_template == 10 && sec3_len == 72) { // mercator uint_char(new_nx,sec3+30); // nx uint_char(new_ny,sec3+34); // ny units = 0.000001; i = lat[ idx(ix0,iy0,nx,ny,cyclic_grid) ] / units; // lat1 int_char(i,sec3+38); i = lon[ idx(ix0,iy0,nx,ny,cyclic_grid) ] / units; // lon1 int_char(i,sec3+42); i = lat[ idx(ix1,iy1,nx,ny,cyclic_grid) ] / units; // lat2 int_char(i,sec3+51); i = lon[ idx(ix1,iy1,nx,ny,cyclic_grid) ] / units; // lon2 int_char(i,sec3+55); } else { can_subset = 0; } } // copy data to a new array if (can_subset) { uint_char(new_ndata, sec3+6); new_data = (float *) malloc(new_ndata * sizeof(float)); #pragma omp parallel for private(i,j,k) for(j = iy0; j <= iy1; j++) { k = (j-iy0)*(ix1-ix0+1); for(i = ix0; i <= ix1; i++) { new_data[(i-ix0) + k ] = data[ idx(i,j,nx,ny,cyclic_grid) ]; } } } else { new_ndata = ndata; new_data = (float *) malloc(new_ndata * sizeof(float)); for (k = 0; k < ndata; k++) new_data[k] = data[k]; new_nx = nx; new_ny = ny; } set_order(new_sec, output_order); grib_wrt(new_sec, new_data, new_ndata, new_nx, new_ny, use_scale, dec_scale, bin_scale, wanted_bits, max_bits, grib_type, out); if (flush_mode) fflush(out); free(new_data); free(sec3); return 0; }
int gauss2ll(unsigned char **sec, double **llat, double **llon) { int nlat; /* in grib, number of latitudes must be even! */ double dx, e, w, south, north, lat1, lon1, lat2, lon2, *ylat; int isouth, inorth; double units; double *lat, *lon; int basic_ang, sub_ang; int i,j, n; unsigned int k; unsigned char *gds; int nnx, nny, nres, nscan; unsigned int nnpnts; get_nxny(sec, &nnx, &nny, &nnpnts, &nres, &nscan); gds = sec[3]; nlat = 2 * GDS_Gaussian_nlat(gds); /* figure out angle units */ basic_ang = GDS_Gaussian_basic_ang(gds); sub_ang = GDS_Gaussian_sub_ang(gds); units = basic_ang == 0 ? 0.000001 : (double) basic_ang / (double) sub_ang; lat1 = GDS_Gaussian_lat1(gds) * units; lat2 = GDS_Gaussian_lat2(gds) * units; lon1 = GDS_Gaussian_lon1(gds) * units; lon2 = GDS_Gaussian_lon2(gds) * units; if (lon1 < 0.0 || lon2 < 0.0 || lon1 > 360.0 || lon2 > 360.0) fatal_error("BAD GDS lon",""); if (lat1 < -90.0 || lat2 < -90.0 || lat1 > 90.0 || lat2 > 90.0) fatal_error("BAD GDS lat",""); /* find S latitude and dy */ if (GDS_Scan_y(nscan)) { south = lat1; north = lat2; } else { south = lat2; north = lat1; } if (south > north) fatal_error("gaussian grid: lat1 and lat2 inconsistent with scan order",""); if (nny == -1) { fprintf(stderr,"Sorry code does not handle variable ny yet\n"); return 0; } if (nny > nlat || nny < 0) fatal_error_i("gauss2ll: bad ny %d",nny); if ((*llat = (double *) malloc(nnpnts * sizeof(double))) == NULL) { fatal_error("gauss2ll memory allocation failed",""); } if ((*llon = (double *) malloc(nnpnts * sizeof(double))) == NULL) { fatal_error("gauss2ll memory allocation failed",""); } lat = *llat; lon = *llon; /* do latitudes first */ ylat = (double *) malloc(sizeof(double) * nlat); /* calculate Gaussian latitudes */ gauss2lats(nlat, ylat); /* find index of south and north */ isouth = inorth = -1; for (i = 0; i < nlat; i++) { if (fabs(south - ylat[i]) < LATERR) { isouth = i; break; } } for (i = 0; i < nlat; i++) { if (fabs(north - ylat[i]) < LATERR) { inorth = i; break; } } if (isouth < 0 || inorth < 0) fatal_error("gauss2ll: lat1/lat2 not a Gaussian latitude",""); if (inorth - isouth + 1 != nny) fatal_error("gauss2ll: lat1/lat2 not consistent with ny",""); n = 0; if (nnx >= 0) { /* regular grid */ #pragma omp parallel for private(i,j) for (j = 0; j < nny; j++) { for (i = 0; i < nnx; i++) { lat[i+j*nnx] = ylat[j+isouth]; } } } else { /* quasi regular grid */ for (j = 0; j < nny; j++) { for (i = 0; i < variable_dim[j]; i++) { lat[n++] = ylat[j+isouth]; } } } free(ylat); /* now for the longitudes */ if (GDS_Scan_x(nscan)) { e = lon1; w = lon2; } else { e = lon2; w = lon1; } if (e > w) w += 360.0; if (e < 0.0) { e += 360.0; w += 360.0; } if (e >= 360.0) { e -= 360.0; w -= 360.0; } if (nnx >= 0) { dx = (w-e) / (nnx-1); #pragma omp parallel { #pragma omp for private(i) for (i = 0; i < nnx; i++) { lon[i] = e + (dx * i) >= 360.0 ? e + (dx * i) - 360.0 : e + (dx * i); } #pragma omp for private(i,j) for (j = 1; j < nny; j++) { for (i = 0; i < nnx; i++) { lon[i+j*nnx] = lon[i]; } } } } else { n = 0; for (j = 0; j < nny; j++) { dx = (w-e) / (variable_dim[j]-1); for (i = 0; i < variable_dim[j]; i++) { lon[n++] = e + (dx * i) >= 360.0 ? e + (dx * i) - 360.0 : e + (dx * i); } } } return 0; } /* end gauss2ll() */