static int input_data(struct interp_params *params, int first_row, int last_row, struct fcell_triple *points, int fdsmooth, int fdinp, int inp_rows, int inp_cols, double zmin, double inp_ns_res, double inp_ew_res) { double x, y, sm; /* input data and smoothing */ int m1, m2; /* loop counters */ static FCELL *cellinp = NULL; /* cell buffer for input data */ static FCELL *cellsmooth = NULL; /* cell buffer for smoothing */ if (!cellinp) cellinp = Rast_allocate_f_buf(); if (!cellsmooth) cellsmooth = Rast_allocate_f_buf(); for (m1 = 0; m1 <= last_row - first_row; m1++) { Rast_get_f_row(fdinp, cellinp, inp_rows - m1 - first_row); if (fdsmooth >= 0) Rast_get_f_row(fdsmooth, cellsmooth, inp_rows - m1 - first_row); y = params->y_orig + (m1 + first_row - 1 + 0.5) * inp_ns_res; for (m2 = 0; m2 < inp_cols; m2++) { x = params->x_orig + (m2 + 0.5) * inp_ew_res; /* * z = cellinp[m2]*params->zmult; */ if (fdsmooth >= 0) sm = (double)cellsmooth[m2]; else sm = 0.01; points[m1 * inp_cols + m2].x = x - params->x_orig; points[m1 * inp_cols + m2].y = y - params->y_orig; if (!Rast_is_f_null_value(cellinp + m2)) { points[m1 * inp_cols + m2].z = cellinp[m2] * params->zmult - zmin; } else { Rast_set_f_null_value(&(points[m1 * inp_cols + m2].z), 1); } /* fprintf (stdout,"sm: %f\n",sm); */ points[m1 * inp_cols + m2].smooth = sm; } } return 1; }
/*! \brief Load raster map as floating point map Calling function must have already allocated space in buff for wind->rows * wind->cols floats. This routine simply loads the map into a 2d array by repetitve calls to get_f_raster_row. \param wind current window \param map_name raster map name \param[out] buff data buffer \param[out] nullmap null map buffer \param[out] has_null indicates if raster map contains null-data \return 1 on success \return 0 on failure */ int Gs_loadmap_as_float(struct Cell_head *wind, const char *map_name, float *buff, struct BM *nullmap, int *has_null) { FILEDESC cellfile; const char *map_set; int offset, row, col; G_debug(3, "Gs_loadmap_as_float(): name=%s", map_name); map_set = G_find_raster2(map_name, ""); if (!map_set) { G_warning(_("Raster map <%s> not found"), map_name); return 0; } *has_null = 0; cellfile = Rast_open_old(map_name, map_set); G_message(_("Loading raster map <%s>..."), G_fully_qualified_name(map_name, map_set)); for (row = 0; row < wind->rows; row++) { offset = row * wind->cols; Rast_get_f_row(cellfile, &(buff[offset]), row); G_percent(row, wind->rows, 2); for (col = 0; col < wind->cols; col++) { if (Rast_is_f_null_value(buff + offset + col)) { *has_null = 1; BM_set(nullmap, col, row, 1); } /* set nm */ } } G_percent(1, 1, 1); G_debug(4, " has_null=%d", *has_null); Rast_close(cellfile); return (1); }
int main(int argc, char *argv[]) { int m1; struct FPRange range; DCELL cellmin, cellmax; FCELL *cellrow, fcellmin; struct GModule *module; struct { struct Option *input, *elev, *slope, *aspect, *pcurv, *tcurv, *mcurv, *smooth, *maskmap, *zmult, *fi, *segmax, *npmin, *res_ew, *res_ns, *overlap, *theta, *scalex; } parm; struct { struct Flag *deriv, *cprght; } flag; G_gisinit(argv[0]); module = G_define_module(); G_add_keyword(_("raster")); G_add_keyword(_("resample")); module->description = _("Reinterpolates and optionally computes topographic analysis from " "input raster map to a new raster map (possibly with " "different resolution) using regularized spline with " "tension and smoothing."); parm.input = G_define_standard_option(G_OPT_R_INPUT); parm.res_ew = G_define_option(); parm.res_ew->key = "ew_res"; parm.res_ew->type = TYPE_DOUBLE; parm.res_ew->required = YES; parm.res_ew->description = _("Desired east-west resolution"); parm.res_ns = G_define_option(); parm.res_ns->key = "ns_res"; parm.res_ns->type = TYPE_DOUBLE; parm.res_ns->required = YES; parm.res_ns->description = _("Desired north-south resolution"); parm.elev = G_define_option(); parm.elev->key = "elev"; parm.elev->type = TYPE_STRING; parm.elev->required = NO; parm.elev->gisprompt = "new,cell,raster"; parm.elev->description = _("Output z-file (elevation) map"); parm.elev->guisection = _("Output"); parm.slope = G_define_option(); parm.slope->key = "slope"; parm.slope->type = TYPE_STRING; parm.slope->required = NO; parm.slope->gisprompt = "new,cell,raster"; parm.slope->description = _("Output slope map (or fx)"); parm.slope->guisection = _("Output"); parm.aspect = G_define_option(); parm.aspect->key = "aspect"; parm.aspect->type = TYPE_STRING; parm.aspect->required = NO; parm.aspect->gisprompt = "new,cell,raster"; parm.aspect->description = _("Output aspect map (or fy)"); parm.aspect->guisection = _("Output"); parm.pcurv = G_define_option(); parm.pcurv->key = "pcurv"; parm.pcurv->type = TYPE_STRING; parm.pcurv->required = NO; parm.pcurv->gisprompt = "new,cell,raster"; parm.pcurv->description = _("Output profile curvature map (or fxx)"); parm.pcurv->guisection = _("Output"); parm.tcurv = G_define_option(); parm.tcurv->key = "tcurv"; parm.tcurv->type = TYPE_STRING; parm.tcurv->required = NO; parm.tcurv->gisprompt = "new,cell,raster"; parm.tcurv->description = _("Output tangential curvature map (or fyy)"); parm.tcurv->guisection = _("Output"); parm.mcurv = G_define_option(); parm.mcurv->key = "mcurv"; parm.mcurv->type = TYPE_STRING; parm.mcurv->required = NO; parm.mcurv->gisprompt = "new,cell,raster"; parm.mcurv->description = _("Output mean curvature map (or fxy)"); parm.mcurv->guisection = _("Output"); parm.smooth = G_define_option(); parm.smooth->key = "smooth"; parm.smooth->type = TYPE_STRING; parm.smooth->required = NO; parm.smooth->gisprompt = "old,cell,raster"; parm.smooth->description = _("Name of raster map containing smoothing"); parm.smooth->guisection = _("Settings"); parm.maskmap = G_define_option(); parm.maskmap->key = "maskmap"; parm.maskmap->type = TYPE_STRING; parm.maskmap->required = NO; parm.maskmap->gisprompt = "old,cell,raster"; parm.maskmap->description = _("Name of raster map to be used as mask"); parm.maskmap->guisection = _("Settings"); parm.overlap = G_define_option(); parm.overlap->key = "overlap"; parm.overlap->type = TYPE_INTEGER; parm.overlap->required = NO; parm.overlap->answer = OVERLAP; parm.overlap->description = _("Rows/columns overlap for segmentation"); parm.overlap->guisection = _("Settings"); parm.zmult = G_define_option(); parm.zmult->key = "zmult"; parm.zmult->type = TYPE_DOUBLE; parm.zmult->answer = ZMULT; parm.zmult->required = NO; parm.zmult->description = _("Multiplier for z-values"); parm.zmult->guisection = _("Settings"); parm.fi = G_define_option(); parm.fi->key = "tension"; parm.fi->type = TYPE_DOUBLE; parm.fi->answer = TENSION; parm.fi->required = NO; parm.fi->description = _("Spline tension value"); parm.fi->guisection = _("Settings"); parm.theta = G_define_option(); parm.theta->key = "theta"; parm.theta->type = TYPE_DOUBLE; parm.theta->required = NO; parm.theta->description = _("Anisotropy angle (in degrees)"); parm.theta->guisection = _("Anisotropy"); parm.scalex = G_define_option(); parm.scalex->key = "scalex"; parm.scalex->type = TYPE_DOUBLE; parm.scalex->required = NO; parm.scalex->description = _("Anisotropy scaling factor"); parm.scalex->guisection = _("Anisotropy"); flag.cprght = G_define_flag(); flag.cprght->key = 't'; flag.cprght->description = _("Use dnorm independent tension"); flag.deriv = G_define_flag(); flag.deriv->key = 'd'; flag.deriv->description = _("Output partial derivatives instead of topographic parameters"); flag.deriv->guisection = _("Output"); if (G_parser(argc, argv)) exit(EXIT_FAILURE); G_get_set_window(&winhd); inp_ew_res = winhd.ew_res; inp_ns_res = winhd.ns_res; inp_cols = winhd.cols; inp_rows = winhd.rows; inp_x_orig = winhd.west; inp_y_orig = winhd.south; input = parm.input->answer; smooth = parm.smooth->answer; maskmap = parm.maskmap->answer; elev = parm.elev->answer; slope = parm.slope->answer; aspect = parm.aspect->answer; pcurv = parm.pcurv->answer; tcurv = parm.tcurv->answer; mcurv = parm.mcurv->answer; cond2 = ((pcurv != NULL) || (tcurv != NULL) || (mcurv != NULL)); cond1 = ((slope != NULL) || (aspect != NULL) || cond2); deriv = flag.deriv->answer; dtens = flag.cprght->answer; ertre = 0.1; if (!G_scan_resolution(parm.res_ew->answer, &ew_res, winhd.proj)) G_fatal_error(_("Unable to read ew_res value")); if (!G_scan_resolution(parm.res_ns->answer, &ns_res, winhd.proj)) G_fatal_error(_("Unable to read ns_res value")); if (sscanf(parm.fi->answer, "%lf", &fi) != 1) G_fatal_error(_("Invalid value for tension")); if (sscanf(parm.zmult->answer, "%lf", &zmult) != 1) G_fatal_error(_("Invalid value for zmult")); if (sscanf(parm.overlap->answer, "%d", &overlap) != 1) G_fatal_error(_("Invalid value for overlap")); if (parm.theta->answer) { if (sscanf(parm.theta->answer, "%lf", &theta) != 1) G_fatal_error(_("Invalid value for theta")); } if (parm.scalex->answer) { if (sscanf(parm.scalex->answer, "%lf", &scalex) != 1) G_fatal_error(_("Invalid value for scalex")); if (!parm.theta->answer) G_fatal_error(_("When using anisotropy both theta and scalex must be specified")); } /* * G_set_embedded_null_value_mode(1); */ outhd.ew_res = ew_res; outhd.ns_res = ns_res; outhd.east = winhd.east; outhd.west = winhd.west; outhd.north = winhd.north; outhd.south = winhd.south; outhd.proj = winhd.proj; outhd.zone = winhd.zone; G_adjust_Cell_head(&outhd, 0, 0); ew_res = outhd.ew_res; ns_res = outhd.ns_res; nsizc = outhd.cols; nsizr = outhd.rows; disk = nsizc * nsizr * sizeof(int); az = G_alloc_vector(nsizc + 1); if (cond1) { adx = G_alloc_vector(nsizc + 1); ady = G_alloc_vector(nsizc + 1); if (cond2) { adxx = G_alloc_vector(nsizc + 1); adyy = G_alloc_vector(nsizc + 1); adxy = G_alloc_vector(nsizc + 1); } } if (smooth != NULL) { fdsmooth = Rast_open_old(smooth, ""); Rast_get_cellhd(smooth, "", &smhd); if ((winhd.ew_res != smhd.ew_res) || (winhd.ns_res != smhd.ns_res)) G_fatal_error(_("Map <%s> is the wrong resolution"), smooth); if (Rast_read_fp_range(smooth, "", &range) >= 0) Rast_get_fp_range_min_max(&range, &cellmin, &cellmax); fcellmin = (float)cellmin; if (Rast_is_f_null_value(&fcellmin) || fcellmin < 0.0) G_fatal_error(_("Smoothing values can not be negative or NULL")); } Rast_get_cellhd(input, "", &inphd); if ((winhd.ew_res != inphd.ew_res) || (winhd.ns_res != inphd.ns_res)) G_fatal_error(_("Input map resolution differs from current region resolution!")); fdinp = Rast_open_old(input, ""); sdisk = 0; if (elev != NULL) sdisk += disk; if (slope != NULL) sdisk += disk; if (aspect != NULL) sdisk += disk; if (pcurv != NULL) sdisk += disk; if (tcurv != NULL) sdisk += disk; if (mcurv != NULL) sdisk += disk; G_message(_("Processing all selected output files will require")); if (sdisk > 1024) { if (sdisk > 1024 * 1024) { if (sdisk > 1024 * 1024 * 1024) { G_message(_("%.2f GB of disk space for temp files."), sdisk / (1024. * 1024. * 1024.)); } else G_message(_("%.2f MB of disk space for temp files."), sdisk / (1024. * 1024.)); } else G_message(_("%.2f KB of disk space for temp files."), sdisk / 1024.); } else G_message(_("%d bytes of disk space for temp files."), sdisk); fstar2 = fi * fi / 4.; tfsta2 = fstar2 + fstar2; deltx = winhd.east - winhd.west; delty = winhd.north - winhd.south; xmin = winhd.west; xmax = winhd.east; ymin = winhd.south; ymax = winhd.north; if (smooth != NULL) smc = -9999; else smc = 0.01; if (Rast_read_fp_range(input, "", &range) >= 0) { Rast_get_fp_range_min_max(&range, &cellmin, &cellmax); } else { cellrow = Rast_allocate_f_buf(); for (m1 = 0; m1 < inp_rows; m1++) { Rast_get_f_row(fdinp, cellrow, m1); Rast_row_update_fp_range(cellrow, m1, &range, FCELL_TYPE); } Rast_get_fp_range_min_max(&range, &cellmin, &cellmax); } fcellmin = (float)cellmin; if (Rast_is_f_null_value(&fcellmin)) G_fatal_error(_("Maximum value of a raster map is NULL.")); zmin = (double)cellmin *zmult; zmax = (double)cellmax *zmult; G_debug(1, "zmin=%f, zmax=%f", zmin, zmax); if (fd4 != NULL) fprintf(fd4, "deltx,delty %f %f \n", deltx, delty); create_temp_files(); IL_init_params_2d(¶ms, NULL, 1, 1, zmult, KMIN, KMAX, maskmap, outhd.rows, outhd.cols, az, adx, ady, adxx, adyy, adxy, fi, MAXPOINTS, SCIK1, SCIK2, SCIK3, smc, elev, slope, aspect, pcurv, tcurv, mcurv, dmin, inp_x_orig, inp_y_orig, deriv, theta, scalex, Tmp_fd_z, Tmp_fd_dx, Tmp_fd_dy, Tmp_fd_xx, Tmp_fd_yy, Tmp_fd_xy, NULL, NULL, 0, NULL); /* In the above line, the penultimate argument is supposed to be a * deviations file pointer. None is obvious, so I used NULL. */ /* The 3rd and 4th argument are int-s, elatt and smatt (from the function * definition. The value 1 seemed like a good placeholder... or not. */ IL_init_func_2d(¶ms, IL_grid_calc_2d, IL_matrix_create, IL_check_at_points_2d, IL_secpar_loop_2d, IL_crst, IL_crstg, IL_write_temp_2d); G_message(_("Temporarily changing the region to desired resolution ...")); Rast_set_window(&outhd); bitmask = IL_create_bitmask(¶ms); /* change region to initial region */ G_message(_("Changing back to the original region ...")); Rast_set_window(&winhd); ertot = 0.; cursegm = 0; G_message(_("Percent complete: ")); NPOINT = IL_resample_interp_segments_2d(¶ms, bitmask, zmin, zmax, &zminac, &zmaxac, &gmin, &gmax, &c1min, &c1max, &c2min, &c2max, &ertot, nsizc, &dnorm, overlap, inp_rows, inp_cols, fdsmooth, fdinp, ns_res, ew_res, inp_ns_res, inp_ew_res, dtens); G_message(_("dnorm in mainc after grid before out1= %f"), dnorm); if (NPOINT < 0) { clean(); G_fatal_error(_("split_and_interpolate() failed")); } if (fd4 != NULL) fprintf(fd4, "max. error found = %f \n", ertot); G_free_vector(az); if (cond1) { G_free_vector(adx); G_free_vector(ady); if (cond2) { G_free_vector(adxx); G_free_vector(adyy); G_free_vector(adxy); } } G_message(_("dnorm in mainc after grid before out2= %f"), dnorm); if (IL_resample_output_2d(¶ms, zmin, zmax, zminac, zmaxac, c1min, c1max, c2min, c2max, gmin, gmax, ertot, input, &dnorm, &outhd, &winhd, smooth, NPOINT) < 0) { clean(); G_fatal_error(_("Unable to write raster maps -- try increasing cell size")); } G_free(zero_array_cell); clean(); if (fd4) fclose(fd4); Rast_close(fdinp); if (smooth != NULL) Rast_close(fdsmooth); G_done_msg(" "); exit(EXIT_SUCCESS); }
int COGRR1(double x_or, double y_or, double z_or, int n_rows, int n_cols, int n_levs, int n_points, struct quadruple *points, struct point_3d skip_point) /*C C INTERPOLATION BY FUNCTIONAL METHOD : TPS + complete regul. c */ { int secpar_loop(); static double *w2 = NULL; static double *wz2 = NULL; static double *wz1 = NULL; double amaxa; double stepix, stepiy, stepiz, RO, xx, yy, zz, xg, yg, zg, xx2; double wm, dx, dy, dz, dxx, dyy, dxy, dxz, dyz, dzz, h, bmgd1, bmgd2, etar, zcon, r, ww, wz, r2, hcell, zzcell2, etarcell, rcell, wwcell, zzcell; double x_crs,x_crsd,x_crsdd,x_crsdr2; int n1, k1, k2, k, i1, l, l1, n4, n5, m, i; int NGST, LSIZE, ngstc, nszc, ngstr, nszr, ngstl, nszl; int POINT(); int ind, ind1; static int first_time_z = 1; off_t offset, offset1, offset2; int bmask = 1; static FCELL *cell = NULL; int cond1 = (gradient != NULL) || (aspect1 != NULL) || (aspect2 != NULL); int cond2 = (ncurv != NULL) || (gcurv != NULL) || (mcurv != NULL); #define CEULER .57721566 /* C c character*32 fncdsm c normalization c */ offset1 = nsizr * nsizc; stepix = ew_res / dnorm; stepiy = ns_res / dnorm; stepiz = tb_res / dnorm; if (!w2) { if (!(w2 = (double *)G_malloc(sizeof(double) * (KMAX2 + 1)))) { clean(); G_fatal_error(_("Not enough memory for %s"), "w2"); } } if (!wz2) { if (!(wz2 = (double *)G_malloc(sizeof(double) * (KMAX2 + 1)))) { clean(); G_fatal_error(_("Not enough memory for %s"), "wz2"); } } if (!wz1) { if (!(wz1 = (double *)G_malloc(sizeof(double) * (KMAX2 + 1)))) { clean(); G_fatal_error(_("Not enough memory for %s"), "wz1"); } } if (cell == NULL) cell = Rast_allocate_f_buf(); for (i = 1; i <= n_points; i++) { points[i - 1].x = (points[i - 1].x - x_or) / dnorm; points[i - 1].y = (points[i - 1].y - y_or) / dnorm; points[i - 1].z = (points[i - 1].z - z_or) / dnorm; } if (cv) { skip_point.x = (skip_point.x - x_or) / dnorm; skip_point.y = (skip_point.y - y_or) / dnorm; skip_point.z = (skip_point.z - z_or) / dnorm; } n1 = n_points + 1; /* C C GENERATION OF MATRIX C C FIRST COLUMN C */ A[1] = 0.; for (k = 1; k <= n_points; k++) { i1 = k + 1; A[i1] = 1.; } /* C C OTHER COLUMNS C */ RO = rsm; for (k = 1; k <= n_points; k++) { k1 = k * n1 + 1; k2 = k + 1; i1 = k1 + k; if (rsm < 0.) { /*indicates variable smoothing */ A[i1] = points[k - 1].sm; } else { A[i1] = RO; /* constant smoothing */ } for (l = k2; l <= n_points; l++) { xx = points[k - 1].x - points[l - 1].x; yy = points[k - 1].y - points[l - 1].y; zz = points[k - 1].z - points[l - 1].z; r = sqrt(xx * xx + yy * yy + zz * zz); etar = (fi * r) / 2.; if (etar == 0.) { /* printf ("ident. points in segm. \n"); printf ("x[%d]=%lf,x[%d]=%lf,y[%d]=%lf,y[%d]=%lf\n", k - 1, points[k - 1].x, l - 1, points[l - 1].x, k - 1, points[k - 1].y, l - 1, points[l - 1].y); */ } i1 = k1 + l; A[i1] = crs(etar); } } /* C C SYMMETRISATION C */ amaxa = 1.; for (k = 1; k <= n1; k++) { k1 = (k - 1) * n1; k2 = k + 1; for (l = k2; l <= n1; l++) { m = (l - 1) * n1 + k; A[m] = A[k1 + l]; amaxa = amax1(A[m], amaxa); } } /* C RIGHT SIDE C */ n4 = n1 * n1 + 1; A[n4] = 0.; for (l = 1; l <= n_points; l++) { l1 = n4 + l; A[l1] = points[l - 1].w; } n5 = n1 * (n1 + 1); for (i = 1; i <= n5; i++) A[i] = A[i] / amaxa; /* SOLVING OF SYSTEM */ if (LINEQS(n1, n1, 1, &NERROR, &DETERM)) { for (k = 1; k <= n_points; k++) { l = n4 + k; b[k] = A[l]; } b[n_points + 1] = A[n4]; POINT(n_points, points, skip_point); if (cv) return 1; if (devi != NULL && sig1 == 1) return 1; /* C C INTERPOLATION * MOST INNER LOOPS ! C */ NGST = 1; LSIZE = 0; ngstc = (int)(x_or / ew_res + 0.5) + 1; nszc = ngstc + n_cols - 1; ngstr = (int)(y_or / ns_res + 0.5) + 1; nszr = ngstr + n_rows - 1; ngstl = (int)(z_or / tb_res + 0.5) + 1; nszl = ngstl + n_levs - 1; /* fprintf(stderr," Progress percentage for each segment ..." ); */ /*fprintf(stderr,"Before loops,ngstl = %d,nszl =%d\n",ngstl,nszl); */ for (i = ngstl; i <= nszl; i++) { /*fprintf(stderr,"level=%d\n",i); */ /* G_percent(i, nszl, 2); */ offset = offset1 * (i - 1); /* levels offset */ zg = (i - ngstl) * stepiz; for (m = 1; m <= n_points; m++) { wz = zg - points[m - 1].z; wz1[m] = wz; wz2[m] = wz * wz; } for (k = ngstr; k <= nszr; k++) { yg = (k - ngstr) * stepiy; for (m = 1; m <= n_points; m++) { wm = yg - points[m - 1].y; w[m] = wm; w2[m] = wm * wm; } if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) Rast_get_f_row(fdcell, cell, n_rows_in - k); for (l = ngstc; l <= nszc; l++) { LSIZE = LSIZE + 1; if (maskmap != NULL) bmask = BM_get(bitmask, l - 1, k - 1); /*bug fix 02/03/00 jh */ xg = (l - ngstc) * stepix; ww = 0.; wwcell = 0.; dx = 0.; dy = 0.; dz = 0.; dxx = 0.; dxy = 0.; dxz = 0.; dyy = 0.; dyz = 0.; dzz = 0.; /* compute everything for area which is not masked out and where cross_input map doesn't have nulls */ if (bmask == 1 && !(cell && Rast_is_f_null_value(&cell[l - 1]))) { h = b[n1]; hcell = b[n1]; for (m = 1; m <= n_points; m++) { xx = xg - points[m - 1].x; xx2 = xx * xx; if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) { zcon = (double)(cell[l - 1] * zmult - z_or) - z_orig_in * zmult; /* bug fix 02/03/00 jh */ zcon = zcon / dnorm; zzcell = zcon - points[m - 1].z; zzcell2 = zzcell * zzcell; rcell = sqrt(xx2 + w2[m] + zzcell2); etarcell = (fi * rcell) / 2.; hcell = hcell + b[m] * crs(etarcell); } r2 = xx2 + w2[m] + wz2[m]; r = sqrt(r2); etar = (fi * r) / 2.; crs_full( etar,fi, &x_crs, cond1?&x_crsd:NULL, cond2?&x_crsdr2:NULL, cond2?&x_crsdd:NULL ); h = h + b[m] * x_crs; if(cond1) { bmgd1 = b[m] * x_crsd; dx = dx + bmgd1 * xx; dy = dy + bmgd1 * w[m]; dz = dz + bmgd1 * wz1[m]; } if(cond2) { bmgd2 = b[m] * x_crsdd; bmgd1 = b[m] * x_crsdr2; dyy = dyy + bmgd2 * w2[m] + bmgd1 * w2[m]; dzz = dzz + bmgd2 * wz2[m] + bmgd1 * wz2[m]; dxy = dxy + bmgd2 * xx * w[m] + bmgd1 * xx * w[m]; dxz = dxz + bmgd2 * xx * wz1[m] + bmgd1 * xx * wz1[m]; dyz = dyz + bmgd2 * w[m] * wz1[m] + bmgd1 * w[m] * wz1[m]; } } ww = h + wmin; if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) wwcell = hcell + wmin; az[l] = ww; if (first_time_z) { first_time_z = 0; zmaxac = zminac = ww; if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) zmaxacell = zminacell = wwcell; } zmaxac = amax1(ww, zmaxac); zminac = amin1(ww, zminac); if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) { zmaxacell = amax1(wwcell, zmaxacell); zminacell = amin1(wwcell, zminacell); } if ((ww > wmax + 0.1 * (wmax - wmin)) || (ww < wmin - 0.1 * (wmax - wmin))) { static int once = 0; if (!once) { once = 1; fprintf(stderr, "WARNING:\n"); fprintf(stderr, "Overshoot -- increase in tension suggested.\n"); fprintf(stderr, "Overshoot occurs at (%d,%d,%d) cell\n", l, k, i); fprintf(stderr, "The w-value is %lf, wmin is %lf,wmax is %lf\n", ww, wmin, wmax); } } } /* skip here if you are in masked area, ww should be 0 */ az[l] = ww; adx[l] = dx; ady[l] = dy; adz[l] = dz; /* printf("\n %f", ww); */ adxx[l] = dxx; adxy[l] = dxy; adxz[l] = dxz; adyy[l] = dyy; adyz[l] = dyz; adzz[l] = dzz; if ((gradient != NULL) || (aspect1 != NULL) || (aspect2 != NULL) || (ncurv != NULL) || (gcurv != NULL) || (mcurv != NULL)) if (!(secpar_loop(ngstc, nszc, l))) { clean(); G_fatal_error(_("Secpar_loop failed")); } if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) { zero_array_cell[l - 1] = (FCELL) (wwcell); } if (outz != NULL) { zero_array1[l - 1] = (float)(az[l] * sciz); } if (gradient != NULL) { zero_array2[l - 1] = (float)(adx[l]); } if (aspect1 != NULL) { zero_array3[l - 1] = (float)(ady[l]); } if (aspect2 != NULL) { zero_array4[l - 1] = (float)(adz[l]); } if (ncurv != NULL) { zero_array5[l - 1] = (float)(adxx[l]); } if (gcurv != NULL) { zero_array6[l - 1] = (float)(adyy[l]); } if (mcurv != NULL) { zero_array7[l - 1] = (float)(adxy[l]); } } /* columns */ ind = nsizc * (k - 1) + (ngstc - 1); ind1 = ngstc - 1; offset2 = offset + ind; /* rows*cols offset */ if ((cellinp != NULL) && (cellout != NULL) && (i == ngstl)) { G_fseek(Tmp_fd_cell, ((off_t)ind * sizeof(FCELL)), 0); if (! (fwrite (zero_array_cell + ind1, sizeof(FCELL), nszc - ngstc + 1, Tmp_fd_cell))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (outz != NULL) { G_fseek(Tmp_fd_z, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array1 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_z))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (gradient != NULL) { G_fseek(Tmp_fd_dx, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array2 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_dx))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (aspect1 != NULL) { G_fseek(Tmp_fd_dy, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array3 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_dy))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (aspect2 != NULL) { G_fseek(Tmp_fd_dz, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array4 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_dz))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (ncurv != NULL) { G_fseek(Tmp_fd_xx, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array5 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_xx))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (gcurv != NULL) { G_fseek(Tmp_fd_yy, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array6 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_yy))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } if (mcurv != NULL) { G_fseek(Tmp_fd_xy, (off_t)(offset2 * sizeof(float)), 0); if (! (fwrite (zero_array7 + ind1, sizeof(float), nszc - ngstc + 1, Tmp_fd_xy))) { clean(); G_fatal_error (_("Not enough disk space--cannot write files")); } } } } } /* falls here if LINEQS() returns 0 */ /* total++; */ /*fprintf(stderr,"wminac=%lf,wmaxac=%lf\n",zminac,zmaxac); */ return 1; }
int main(int argc, char *argv[]) { char *terrainmap, *seedmap, *lakemap; int rows, cols, in_terran_fd, out_fd, lake_fd, row, col, pases, pass; int lastcount, curcount, start_col = 0, start_row = 0; double east, north, area = 0, volume = 0; FCELL **in_terran, **out_water, water_level, max_depth = 0, min_depth = 0; FCELL water_window[3][3]; struct Option *tmap_opt, *smap_opt, *wlvl_opt, *lake_opt, *sdxy_opt; struct Flag *negative_flag, *overwrite_flag; struct GModule *module; struct Colors colr; struct Cell_head window; struct History history; G_gisinit(argv[0]); module = G_define_module(); G_add_keyword(_("raster")); G_add_keyword(_("hydrology")); G_add_keyword(_("hazard")); G_add_keyword(_("flood")); module->description = _("Fills lake at given point to given level."); tmap_opt = G_define_standard_option(G_OPT_R_ELEV); wlvl_opt = G_define_option(); wlvl_opt->key = "water_level"; wlvl_opt->description = _("Water level"); wlvl_opt->type = TYPE_DOUBLE; wlvl_opt->required = YES; lake_opt = G_define_standard_option(G_OPT_R_OUTPUT); lake_opt->key = "lake"; lake_opt->required = NO; lake_opt->guisection = _("Output"); sdxy_opt = G_define_standard_option(G_OPT_M_COORDS); sdxy_opt->label = _("Seed point coordinates"); sdxy_opt->description = _("Either this coordinates pair or a seed" " map have to be specified"); sdxy_opt->required = NO; sdxy_opt->multiple = NO; sdxy_opt->guisection = _("Seed"); smap_opt = G_define_standard_option(G_OPT_R_MAP); smap_opt->key = "seed"; smap_opt->label = _("Input raster map with given starting point(s) (at least 1 cell > 0)"); smap_opt->description = _("Either this parameter or a coordinates pair have to be specified"); smap_opt->required = NO; smap_opt->guisection = _("Seed"); negative_flag = G_define_flag(); negative_flag->key = 'n'; negative_flag->description = _("Use negative depth values for lake raster map"); overwrite_flag = G_define_flag(); overwrite_flag->key = 'o'; overwrite_flag->description = _("Overwrite seed map with result (lake) map"); overwrite_flag->guisection = _("Output"); if (G_parser(argc, argv)) /* Returns 0 if successful, non-zero otherwise */ exit(EXIT_FAILURE); if (smap_opt->answer && sdxy_opt->answer) G_fatal_error(_("Both seed map and coordinates cannot be specified")); if (!smap_opt->answer && !sdxy_opt->answer) G_fatal_error(_("Seed map or seed coordinates must be set!")); if (sdxy_opt->answer && !lake_opt->answer) G_fatal_error(_("Seed coordinates and output map lake= must be set!")); if (lake_opt->answer && overwrite_flag->answer) G_fatal_error(_("Both lake and overwrite cannot be specified")); if (!lake_opt->answer && !overwrite_flag->answer) G_fatal_error(_("Output lake map or overwrite flag must be set!")); terrainmap = tmap_opt->answer; seedmap = smap_opt->answer; sscanf(wlvl_opt->answer, "%f", &water_level); lakemap = lake_opt->answer; /* If lakemap is set, write to it, else is set overwrite flag and we should write to seedmap. */ if (lakemap) lake_fd = Rast_open_new(lakemap, 1); rows = Rast_window_rows(); cols = Rast_window_cols(); /* If we use x,y as seed... */ if (sdxy_opt->answer) { G_get_window(&window); east = window.east; north = window.north; G_scan_easting(sdxy_opt->answers[0], &east, G_projection()); G_scan_northing(sdxy_opt->answers[1], &north, G_projection()); start_col = (int)Rast_easting_to_col(east, &window); start_row = (int)Rast_northing_to_row(north, &window); if (start_row < 0 || start_row > rows || start_col < 0 || start_col > cols) G_fatal_error(_("Seed point outside the current region")); } /* Open terrain map */ in_terran_fd = Rast_open_old(terrainmap, ""); /* Open seed map */ if (smap_opt->answer) out_fd = Rast_open_old(seedmap, ""); /* Pointers to rows. Row = ptr to 'col' size array. */ in_terran = (FCELL **) G_malloc(rows * sizeof(FCELL *)); out_water = (FCELL **) G_malloc(rows * sizeof(FCELL *)); if (in_terran == NULL || out_water == NULL) G_fatal_error(_("G_malloc: out of memory")); G_debug(1, "Loading maps..."); /* foo_rows[row] == array with data (2d array). */ for (row = 0; row < rows; row++) { in_terran[row] = (FCELL *) G_malloc(cols * sizeof(FCELL)); out_water[row] = (FCELL *) G_calloc(cols, sizeof(FCELL)); /* In newly created space load data from file. */ Rast_get_f_row(in_terran_fd, in_terran[row], row); if (smap_opt->answer) Rast_get_f_row(out_fd, out_water[row], row); G_percent(row + 1, rows, 5); } /* Set seed point */ if (sdxy_opt->answer) /* Check is water level higher than seed point */ if (in_terran[start_row][start_col] >= water_level) G_fatal_error(_("Given water level at seed point is below earth surface. " "Increase water level or move seed point.")); out_water[start_row][start_col] = 1; /* Close seed map for reading. */ if (smap_opt->answer) Rast_close(out_fd); /* Open output map for writing. */ if (lakemap) out_fd = lake_fd; else out_fd = Rast_open_new(seedmap, 1); /* More pases are renudant. Real pases count is controlled by altered cell count. */ pases = (int)(rows * cols) / 2; G_debug(1, "Starting lake filling at level of %8.4f in %d passes. Percent done:", water_level, pases); lastcount = 0; for (pass = 0; pass < pases; pass++) { G_debug(3, "Pass: %d", pass); curcount = 0; /* Move from left upper corner to right lower corner. */ for (row = 0; row < rows; row++) { for (col = 0; col < cols; col++) { /* Loading water data into window. */ load_window_values(out_water, water_window, rows, cols, row, col); /* Cheking presence of water. */ if (is_near_water(water_window) == 1) { if (in_terran[row][col] < water_level) { out_water[row][col] = water_level - in_terran[row][col]; curcount++; } else { out_water[row][col] = 0; /* Cell is higher than water level -> NULL. */ } } } } if (curcount == lastcount) break; /* We done. */ lastcount = curcount; curcount = 0; /* Move backwards - from lower right corner to upper left corner. */ for (row = rows - 1; row >= 0; row--) { for (col = cols - 1; col >= 0; col--) { load_window_values(out_water, water_window, rows, cols, row, col); if (is_near_water(water_window) == 1) { if (in_terran[row][col] < water_level) { out_water[row][col] = water_level - in_terran[row][col]; curcount++; } else { out_water[row][col] = 0; } } } } G_percent(pass + 1, pases, 10); if (curcount == lastcount) break; /* We done. */ lastcount = curcount; } /*pases */ G_percent(pases, pases, 10); /* Show 100%. */ save_map(out_water, out_fd, rows, cols, negative_flag->answer, &min_depth, &max_depth, &area, &volume); G_message(_("Lake depth from %f to %f (specified water level is taken as zero)"), min_depth, max_depth); G_message(_("Lake area %f square meters"), area); G_message(_("Lake volume %f cubic meters"), volume); G_important_message(_("Volume is correct only if lake depth (terrain raster map) is in meters")); /* Close all files. Lake map gets written only now. */ Rast_close(in_terran_fd); Rast_close(out_fd); /* Add blue color gradient from light bank to dark depth */ Rast_init_colors(&colr); if (negative_flag->answer == 1) { Rast_add_f_color_rule(&max_depth, 0, 240, 255, &min_depth, 0, 50, 170, &colr); } else { Rast_add_f_color_rule(&min_depth, 0, 240, 255, &max_depth, 0, 50, 170, &colr); } Rast_write_colors(lakemap, G_mapset(), &colr); Rast_short_history(lakemap, "raster", &history); Rast_command_history(&history); Rast_write_history(lakemap, &history); return EXIT_SUCCESS; }
int extract_points(int z_flag) { struct line_pnts *points = Vect_new_line_struct(); CELL *cellbuf; FCELL *fcellbuf; DCELL *dcellbuf; int row, col; double x, y; int count; switch (data_type) { case CELL_TYPE: cellbuf = Rast_allocate_c_buf(); break; case FCELL_TYPE: fcellbuf = Rast_allocate_f_buf(); break; case DCELL_TYPE: dcellbuf = Rast_allocate_d_buf(); break; } G_message(_("Extracting points...")); count = 1; for (row = 0; row < cell_head.rows; row++) { G_percent(row, n_rows, 2); y = Rast_row_to_northing((double)(row + .5), &cell_head); switch (data_type) { case CELL_TYPE: Rast_get_c_row(input_fd, cellbuf, row); break; case FCELL_TYPE: Rast_get_f_row(input_fd, fcellbuf, row); break; case DCELL_TYPE: Rast_get_d_row(input_fd, dcellbuf, row); break; } for (col = 0; col < cell_head.cols; col++) { int cat, val; double dval; x = Rast_col_to_easting((double)(col + .5), &cell_head); switch (data_type) { case CELL_TYPE: if (Rast_is_c_null_value(cellbuf + col)) continue; val = cellbuf[col]; dval = val; break; case FCELL_TYPE: if (Rast_is_f_null_value(fcellbuf + col)) continue; dval = fcellbuf[col]; break; case DCELL_TYPE: if (Rast_is_d_null_value(dcellbuf + col)) continue; dval = dcellbuf[col]; break; } /* value_flag is used only for CELL type */ cat = (value_flag) ? val : count; Vect_reset_line(points); Vect_reset_cats(Cats); Vect_cat_set(Cats, 1, cat); Vect_append_point(points, x, y, dval); Vect_write_line(&Map, GV_POINT, points, Cats); if ((driver != NULL) && !value_flag) { insert_value(cat, val, dval); } count++; } } G_percent(row, n_rows, 2); switch (data_type) { case CELL_TYPE: G_free(cellbuf); break; case FCELL_TYPE: G_free(fcellbuf); break; case DCELL_TYPE: G_free(dcellbuf); break; } Vect_destroy_line_struct(points); return (1); }