/*----------------------------------------------------------------------------------------*/ double P_Mean_Calc(struct Cell_head *Elaboration, struct Point *obs, int npoints) { int i, mean_count = 0; double mean = 0.0; struct bound_box mean_box; Vect_region_box(Elaboration, &mean_box); mean_box.W -= CONTOUR; mean_box.E += CONTOUR; mean_box.N += CONTOUR; mean_box.S -= CONTOUR; for (i = 0; i < npoints; i++) { /* */ if (Vect_point_in_box (obs[i].coordX, obs[i].coordY, obs[i].coordZ, &mean_box)) { mean_count++; mean += obs[i].coordZ; } } if (mean_count == 0) mean = .0; else mean /= (double)mean_count; return mean; }
double P_estimate_splinestep(struct Map_info *Map, double *dens, double *dist) { int type, npoints = 0; double xmin = 0, xmax = 0, ymin = 0, ymax = 0; double x, y, z; struct line_pnts *points; struct line_cats *categories; struct bound_box region_box; struct Cell_head orig; G_get_set_window(&orig); Vect_region_box(&orig, ®ion_box); points = Vect_new_line_struct(); categories = Vect_new_cats_struct(); Vect_rewind(Map); while ((type = Vect_read_next_line(Map, points, categories)) > 0) { if (!(type & GV_POINT)) continue; x = points->x[0]; y = points->y[0]; if (points->z != NULL) z = points->z[0]; else z = 0.0; /* only use points in current region */ if (Vect_point_in_box(x, y, z, ®ion_box)) { npoints++; if (npoints > 1) { if (xmin > x) xmin = x; else if (xmax < x) xmax = x; if (ymin > y) ymin = y; else if (ymax < y) ymax = y; } else { xmin = xmax = x; ymin = ymax = y; } } } if (npoints > 0) { /* estimated average distance between points in map units */ *dist = sqrt(((xmax - xmin) * (ymax - ymin)) / npoints); /* estimated point density as number of points per square map unit */ *dens = npoints / ((xmax - xmin) * (ymax - ymin)); return 0; } else { return -1; } }
int loadSiteCoordinates(struct Map_info *Map, struct Point **points, int region, struct Cell_head *window, int field, struct cat_list *cat_list) { int i, pointIdx; struct line_pnts *sites; struct line_cats *cats; struct bound_box box; int type; sites = Vect_new_line_struct(); cats = Vect_new_cats_struct(); *points = NULL; pointIdx = 0; /* copy window to box */ Vect_region_box(window, &box); while ((type = Vect_read_next_line(Map, sites, cats)) > -1) { if (type != GV_POINT && !(type & GV_LINES)) continue; if (field > 0 && !Vect_cats_in_constraint(cats, field, cat_list)) continue; for (i = 0; i < sites->n_points; i++) { G_debug(4, "Point: %f|%f|%f", sites->x[i], sites->y[i], sites->z[i]); if (region && !Vect_point_in_box(sites->x[i], sites->y[i], sites->z[i], &box)) continue; G_debug(4, "Point in the box"); if ((pointIdx % ALLOC_CHUNK) == 0) *points = (struct Point *) G_realloc(*points, (pointIdx + ALLOC_CHUNK) * sizeof(struct Point)); (*points)[pointIdx].x = sites->x[i]; (*points)[pointIdx].y = sites->y[i]; (*points)[pointIdx].z = sites->z[i]; pointIdx++; } } if (pointIdx > 0) *points = (struct Point *)G_realloc(*points, (pointIdx + 1) * sizeof(struct Point)); return pointIdx; }
int main(int argc, char *argv[]) { int i, j, nlines, type, field, cat; int fd; /* struct Categories RCats; *//* TODO */ struct Cell_head window; RASTER_MAP_TYPE out_type; CELL *cell; DCELL *dcell; double drow, dcol; char buf[2000]; struct Option *vect_opt, *rast_opt, *field_opt, *col_opt, *where_opt; int Cache_size; struct order *cache; int cur_row; struct GModule *module; struct Map_info Map; struct line_pnts *Points; struct line_cats *Cats; int point; int point_cnt; /* number of points in cache */ int outside_cnt; /* points outside region */ int nocat_cnt; /* points inside region but without category */ int dupl_cnt; /* duplicate categories */ struct bound_box box; int *catexst, *cex; struct field_info *Fi; dbString stmt; dbDriver *driver; int select, norec_cnt, update_cnt, upderr_cnt, col_type; G_gisinit(argv[0]); module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("raster")); G_add_keyword(_("position")); G_add_keyword(_("querying")); G_add_keyword(_("attribute table")); module->description = _("Uploads raster values at positions of vector points to the table."); vect_opt = G_define_standard_option(G_OPT_V_INPUT); vect_opt->key = "vector"; vect_opt->description = _("Name of input vector points map for which to edit attribute table"); rast_opt = G_define_standard_option(G_OPT_R_INPUT); rast_opt->key = "raster"; rast_opt->description = _("Name of existing raster map to be queried"); field_opt = G_define_standard_option(G_OPT_V_FIELD); col_opt = G_define_option(); col_opt->key = "column"; col_opt->type = TYPE_STRING; col_opt->required = YES; col_opt->description = _("Column name (will be updated by raster values)"); where_opt = G_define_standard_option(G_OPT_DB_WHERE); if (G_parser(argc, argv)) exit(EXIT_FAILURE); field = atoi(field_opt->answer); db_init_string(&stmt); Points = Vect_new_line_struct(); Cats = Vect_new_cats_struct(); G_get_window(&window); Vect_region_box(&window, &box); /* T and B set to +/- PORT_DOUBLE_MAX */ /* Open vector */ Vect_set_open_level(2); Vect_open_old(&Map, vect_opt->answer, ""); Fi = Vect_get_field(&Map, field); if (Fi == NULL) G_fatal_error(_("Database connection not defined for layer %d"), field); /* Open driver */ driver = db_start_driver_open_database(Fi->driver, Fi->database); if (driver == NULL) { G_fatal_error(_("Unable to open database <%s> by driver <%s>"), Fi->database, Fi->driver); } /* Open raster */ fd = Rast_open_old(rast_opt->answer, ""); out_type = Rast_get_map_type(fd); /* TODO: Later possibly category labels */ /* if ( Rast_read_cats (name, "", &RCats) < 0 ) G_fatal_error ( "Cannot read category file"); */ /* Check column type */ col_type = db_column_Ctype(driver, Fi->table, col_opt->answer); if (col_type == -1) G_fatal_error(_("Column <%s> not found"), col_opt->answer); if (col_type != DB_C_TYPE_INT && col_type != DB_C_TYPE_DOUBLE) G_fatal_error(_("Column type not supported")); if (out_type == CELL_TYPE && col_type == DB_C_TYPE_DOUBLE) G_warning(_("Raster type is integer and column type is float")); if (out_type != CELL_TYPE && col_type == DB_C_TYPE_INT) G_warning(_("Raster type is float and column type is integer, some data lost!!")); /* Read vector points to cache */ Cache_size = Vect_get_num_primitives(&Map, GV_POINT); /* Note: Some space may be wasted (outside region or no category) */ cache = (struct order *)G_calloc(Cache_size, sizeof(struct order)); point_cnt = outside_cnt = nocat_cnt = 0; nlines = Vect_get_num_lines(&Map); G_debug(1, "Reading %d vector features fom map", nlines); for (i = 1; i <= nlines; i++) { type = Vect_read_line(&Map, Points, Cats, i); G_debug(4, "line = %d type = %d", i, type); /* check type */ if (!(type & GV_POINT)) continue; /* Points only */ /* check region */ if (!Vect_point_in_box(Points->x[0], Points->y[0], 0.0, &box)) { outside_cnt++; continue; } Vect_cat_get(Cats, field, &cat); if (cat < 0) { /* no category of given field */ nocat_cnt++; continue; } G_debug(4, " cat = %d", cat); /* Add point to cache */ drow = Rast_northing_to_row(Points->y[0], &window); dcol = Rast_easting_to_col(Points->x[0], &window); /* a special case. * if north falls at southern edge, or east falls on eastern edge, * the point will appear outside the window. * So, for these edges, bring the point inside the window */ if (drow == window.rows) drow--; if (dcol == window.cols) dcol--; cache[point_cnt].row = (int)drow; cache[point_cnt].col = (int)dcol; cache[point_cnt].cat = cat; cache[point_cnt].count = 1; point_cnt++; } Vect_set_db_updated(&Map); Vect_hist_command(&Map); Vect_close(&Map); G_debug(1, "Read %d vector points", point_cnt); /* Cache may contain duplicate categories, sort by cat, find and remove duplicates * and recalc count and decrease point_cnt */ qsort(cache, point_cnt, sizeof(struct order), by_cat); G_debug(1, "Points are sorted, starting duplicate removal loop"); for (i = 0, j = 1; j < point_cnt; j++) if (cache[i].cat != cache[j].cat) cache[++i] = cache[j]; else cache[i].count++; point_cnt = i + 1; G_debug(1, "%d vector points left after removal of duplicates", point_cnt); /* Report number of points not used */ if (outside_cnt) G_warning(_("%d points outside current region were skipped"), outside_cnt); if (nocat_cnt) G_warning(_("%d points without category were skipped"), nocat_cnt); /* Sort cache by current region row */ qsort(cache, point_cnt, sizeof(struct order), by_row); /* Allocate space for raster row */ if (out_type == CELL_TYPE) cell = Rast_allocate_c_buf(); else dcell = Rast_allocate_d_buf(); /* Extract raster values from file and store in cache */ G_debug(1, "Extracting raster values"); cur_row = -1; for (point = 0; point < point_cnt; point++) { if (cache[point].count > 1) continue; /* duplicate cats */ if (cur_row != cache[point].row) { if (out_type == CELL_TYPE) Rast_get_c_row(fd, cell, cache[point].row); else Rast_get_d_row(fd, dcell, cache[point].row); } cur_row = cache[point].row; if (out_type == CELL_TYPE) { cache[point].value = cell[cache[point].col]; } else { cache[point].dvalue = dcell[cache[point].col]; } } /* point loop */ /* Update table from cache */ G_debug(1, "Updating db table"); /* select existing categories to array (array is sorted) */ select = db_select_int(driver, Fi->table, Fi->key, NULL, &catexst); db_begin_transaction(driver); norec_cnt = update_cnt = upderr_cnt = dupl_cnt = 0; for (point = 0; point < point_cnt; point++) { if (cache[point].count > 1) { G_warning(_("More points (%d) of category %d, value set to 'NULL'"), cache[point].count, cache[point].cat); dupl_cnt++; } /* category exist in DB ? */ cex = (int *)bsearch((void *)&(cache[point].cat), catexst, select, sizeof(int), srch_cat); if (cex == NULL) { /* cat does not exist in DB */ norec_cnt++; G_warning(_("No record for category %d in table <%s>"), cache[point].cat, Fi->table); continue; } sprintf(buf, "update %s set %s = ", Fi->table, col_opt->answer); db_set_string(&stmt, buf); if (out_type == CELL_TYPE) { if (cache[point].count > 1 || Rast_is_c_null_value(&cache[point].value)) { sprintf(buf, "NULL"); } else { sprintf(buf, "%d ", cache[point].value); } } else { /* FCELL or DCELL */ if (cache[point].count > 1 || Rast_is_d_null_value(&cache[point].dvalue)) { sprintf(buf, "NULL"); } else { sprintf(buf, "%.10f", cache[point].dvalue); } } db_append_string(&stmt, buf); sprintf(buf, " where %s = %d", Fi->key, cache[point].cat); db_append_string(&stmt, buf); /* user provides where condition: */ if (where_opt->answer) { sprintf(buf, " AND %s", where_opt->answer); db_append_string(&stmt, buf); } G_debug(3, db_get_string(&stmt)); /* Update table */ if (db_execute_immediate(driver, &stmt) == DB_OK) { update_cnt++; } else { upderr_cnt++; } } G_debug(1, "Committing DB transaction"); db_commit_transaction(driver); G_free(catexst); db_close_database_shutdown_driver(driver); db_free_string(&stmt); /* Report */ G_message(_("%d categories loaded from table"), select); G_message(_("%d categories loaded from vector"), point_cnt); G_message(_("%d categories from vector missing in table"), norec_cnt); G_message(_("%d duplicate categories in vector"), dupl_cnt); if (!where_opt->answer) G_message(_("%d records updated"), update_cnt); G_message(_("%d update errors"), upderr_cnt); exit(EXIT_SUCCESS); }
/*----------------------------------------------------------------------------------------------------------*/ int main(int argc, char *argv[]) { /* Declarations */ int dim_vect, nparameters, BW, npoints; int nsply, nsplx, nsplx_adj, nsply_adj; int nsubregion_col, nsubregion_row; int subregion = 0, nsubregions = 0; const char *dvr, *db, *mapset; char table_name[GNAME_MAX]; char xname[GNAME_MAX], xmapset[GMAPSET_MAX]; double lambda, mean, stepN, stepE, HighThresh, LowThresh; double N_extension, E_extension, edgeE, edgeN; int i, nterrain, count_terrain; int last_row, last_column, flag_auxiliar = FALSE; int *lineVect; double *TN, *Q, *parVect; /* Interpolating and least-square vectors */ double **N, **obsVect, **obsVect_all; /* Interpolation and least-square matrix */ struct Map_info In, Out, Terrain; struct Option *in_opt, *out_opt, *out_terrain_opt, *stepE_opt, *stepN_opt, *lambda_f_opt, *Thresh_A_opt, *Thresh_B_opt; struct Flag *spline_step_flag; struct GModule *module; struct Cell_head elaboration_reg, original_reg; struct Reg_dimens dims; struct bound_box general_box, overlap_box; struct Point *observ; struct lidar_cat *lcat; dbDriver *driver; /*----------------------------------------------------------------------------------------------------------*/ /* Options' declaration */ module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("LIDAR")); module->description = _("Corrects the v.lidar.growing output. It is the last of the three algorithms for LIDAR filtering."); spline_step_flag = G_define_flag(); spline_step_flag->key = 'e'; spline_step_flag->label = _("Estimate point density and distance"); spline_step_flag->description = _("Estimate point density and distance for the input vector points within the current region extends and quit"); in_opt = G_define_standard_option(G_OPT_V_INPUT); in_opt->description = _("Input observation vector map name (v.lidar.growing output)"); out_opt = G_define_standard_option(G_OPT_V_OUTPUT); out_opt->description = _("Output classified vector map name"); out_terrain_opt = G_define_option(); out_terrain_opt->key = "terrain"; out_terrain_opt->type = TYPE_STRING; out_terrain_opt->key_desc = "name"; out_terrain_opt->required = YES; out_terrain_opt->gisprompt = "new,vector,vector"; out_terrain_opt->description = _("Only 'terrain' points output vector map"); stepE_opt = G_define_option(); stepE_opt->key = "ew_step"; stepE_opt->type = TYPE_DOUBLE; stepE_opt->required = NO; stepE_opt->answer = "25"; stepE_opt->description = _("Length of each spline step in the east-west direction"); stepE_opt->guisection = _("Settings"); stepN_opt = G_define_option(); stepN_opt->key = "ns_step"; stepN_opt->type = TYPE_DOUBLE; stepN_opt->required = NO; stepN_opt->answer = "25"; stepN_opt->description = _("Length of each spline step in the north-south direction"); stepN_opt->guisection = _("Settings"); lambda_f_opt = G_define_option(); lambda_f_opt->key = "lambda_c"; lambda_f_opt->type = TYPE_DOUBLE; lambda_f_opt->required = NO; lambda_f_opt->description = _("Regularization weight in reclassification evaluation"); lambda_f_opt->answer = "1"; Thresh_A_opt = G_define_option(); Thresh_A_opt->key = "tch"; Thresh_A_opt->type = TYPE_DOUBLE; Thresh_A_opt->required = NO; Thresh_A_opt->description = _("High threshold for object to terrain reclassification"); Thresh_A_opt->answer = "2"; Thresh_B_opt = G_define_option(); Thresh_B_opt->key = "tcl"; Thresh_B_opt->type = TYPE_DOUBLE; Thresh_B_opt->required = NO; Thresh_B_opt->description = _("Low threshold for terrain to object reclassification"); Thresh_B_opt->answer = "1"; /* Parsing */ G_gisinit(argv[0]); if (G_parser(argc, argv)) exit(EXIT_FAILURE); stepN = atof(stepN_opt->answer); stepE = atof(stepE_opt->answer); lambda = atof(lambda_f_opt->answer); HighThresh = atof(Thresh_A_opt->answer); LowThresh = atof(Thresh_B_opt->answer); if (!(db = G_getenv_nofatal2("DB_DATABASE", G_VAR_MAPSET))) G_fatal_error(_("Unable to read name of database")); if (!(dvr = G_getenv_nofatal2("DB_DRIVER", G_VAR_MAPSET))) G_fatal_error(_("Unable to read name of driver")); /* Setting auxiliar table's name */ if (G_name_is_fully_qualified(out_opt->answer, xname, xmapset)) { sprintf(table_name, "%s_aux", xname); } else sprintf(table_name, "%s_aux", out_opt->answer); /* Something went wrong in a previous v.lidar.correction execution */ if (db_table_exists(dvr, db, table_name)) { /* Start driver and open db */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); db_set_error_handler_driver(driver); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Old auxiliar table could not be dropped")); db_close_database_shutdown_driver(driver); } /* Checking vector names */ Vect_check_input_output_name(in_opt->answer, out_opt->answer, G_FATAL_EXIT); /* Open input vector */ if ((mapset = G_find_vector2(in_opt->answer, "")) == NULL) G_fatal_error(_("Vector map <%s> not found"), in_opt->answer); Vect_set_open_level(1); /* without topology */ if (1 > Vect_open_old(&In, in_opt->answer, mapset)) G_fatal_error(_("Unable to open vector map <%s>"), in_opt->answer); /* Input vector must be 3D */ if (!Vect_is_3d(&In)) G_fatal_error(_("Input vector map <%s> is not 3D!"), in_opt->answer); /* Estimate point density and mean distance for current region */ if (spline_step_flag->answer) { double dens, dist; if (P_estimate_splinestep(&In, &dens, &dist) == 0) { G_message("Estimated point density: %.4g", dens); G_message("Estimated mean distance between points: %.4g", dist); } else G_warning(_("No points in current region!")); Vect_close(&In); exit(EXIT_SUCCESS); } /* Open output vector */ if (0 > Vect_open_new(&Out, out_opt->answer, WITH_Z)) { Vect_close(&In); G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); } if (0 > Vect_open_new(&Terrain, out_terrain_opt->answer, WITH_Z)) { Vect_close(&In); Vect_close(&Out); G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); } /* Copy vector Head File */ Vect_copy_head_data(&In, &Out); Vect_hist_copy(&In, &Out); Vect_hist_command(&Out); Vect_copy_head_data(&In, &Terrain); Vect_hist_copy(&In, &Terrain); Vect_hist_command(&Terrain); /* Start driver and open db */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); db_set_error_handler_driver(driver); /* Create auxiliar table */ if ((flag_auxiliar = P_Create_Aux2_Table(driver, table_name)) == FALSE) { Vect_close(&In); Vect_close(&Out); Vect_close(&Terrain); exit(EXIT_FAILURE); } db_create_index2(driver, table_name, "ID"); /* sqlite likes that ??? */ db_close_database_shutdown_driver(driver); driver = db_start_driver_open_database(dvr, db); /* Setting regions and boxes */ G_get_set_window(&original_reg); G_get_set_window(&elaboration_reg); Vect_region_box(&elaboration_reg, &overlap_box); Vect_region_box(&elaboration_reg, &general_box); /*------------------------------------------------------------------ | Subdividing and working with tiles: | Each original region will be divided into several subregions. | Each one will be overlaped by its neighbouring subregions. | The overlapping is calculated as a fixed OVERLAP_SIZE times | the largest spline step plus 2 * edge ----------------------------------------------------------------*/ /* Fixing parameters of the elaboration region */ P_zero_dim(&dims); nsplx_adj = NSPLX_MAX; nsply_adj = NSPLY_MAX; if (stepN > stepE) dims.overlap = OVERLAP_SIZE * stepN; else dims.overlap = OVERLAP_SIZE * stepE; P_get_edge(P_BILINEAR, &dims, stepE, stepN); P_set_dim(&dims, stepE, stepN, &nsplx_adj, &nsply_adj); G_verbose_message(n_("adjusted EW spline %d", "adjusted EW splines %d", nsplx_adj), nsplx_adj); G_verbose_message(n_("adjusted NS spline %d", "adjusted NS splines %d", nsply_adj), nsply_adj); /* calculate number of subregions */ edgeE = dims.ew_size - dims.overlap - 2 * dims.edge_v; edgeN = dims.sn_size - dims.overlap - 2 * dims.edge_h; N_extension = original_reg.north - original_reg.south; E_extension = original_reg.east - original_reg.west; nsubregion_col = ceil(E_extension / edgeE) + 0.5; nsubregion_row = ceil(N_extension / edgeN) + 0.5; if (nsubregion_col < 0) nsubregion_col = 0; if (nsubregion_row < 0) nsubregion_row = 0; nsubregions = nsubregion_row * nsubregion_col; elaboration_reg.south = original_reg.north; last_row = FALSE; while (last_row == FALSE) { /* For each row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_ROW); if (elaboration_reg.north > original_reg.north) { /* First row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_ROW); } if (elaboration_reg.south <= original_reg.south) { /* Last row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_ROW); last_row = TRUE; } nsply = ceil((elaboration_reg.north - elaboration_reg.south) / stepN) + 0.5; /* if (nsply > NSPLY_MAX) { nsply = NSPLY_MAX; } */ G_debug(1, _("nsply = %d"), nsply); elaboration_reg.east = original_reg.west; last_column = FALSE; while (last_column == FALSE) { /* For each column */ subregion++; if (nsubregions > 1) G_message(_("subregion %d of %d"), subregion, nsubregions); P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_COLUMN); if (elaboration_reg.west < original_reg.west) { /* First column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_COLUMN); } if (elaboration_reg.east >= original_reg.east) { /* Last column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_COLUMN); last_column = TRUE; } nsplx = ceil((elaboration_reg.east - elaboration_reg.west) / stepE) + 0.5; /* if (nsplx > NSPLX_MAX) { nsplx = NSPLX_MAX; } */ G_debug(1, _("nsplx = %d"), nsplx); dim_vect = nsplx * nsply; G_debug(1, _("read vector region map")); observ = P_Read_Vector_Correction(&In, &elaboration_reg, &npoints, &nterrain, dim_vect, &lcat); G_debug(5, _("npoints = %d, nterrain = %d"), npoints, nterrain); if (npoints > 0) { /* If there is any point falling into elaboration_reg. */ count_terrain = 0; nparameters = nsplx * nsply; /* Mean calculation */ G_debug(3, _("Mean calculation")); mean = P_Mean_Calc(&elaboration_reg, observ, npoints); /*Least Squares system */ BW = P_get_BandWidth(P_BILINEAR, nsply); /* Bilinear interpolation */ N = G_alloc_matrix(nparameters, BW); /* Normal matrix */ TN = G_alloc_vector(nparameters); /* vector */ parVect = G_alloc_vector(nparameters); /* Bilinear parameters vector */ obsVect = G_alloc_matrix(nterrain + 1, 3); /* Observation vector with terrain points */ obsVect_all = G_alloc_matrix(npoints + 1, 3); /* Observation vector with all points */ Q = G_alloc_vector(nterrain + 1); /* "a priori" var-cov matrix */ lineVect = G_alloc_ivector(npoints + 1); /* Setting obsVect vector & Q matrix */ G_debug(3, _("Only TERRAIN points")); for (i = 0; i < npoints; i++) { if (observ[i].cat == TERRAIN_SINGLE) { obsVect[count_terrain][0] = observ[i].coordX; obsVect[count_terrain][1] = observ[i].coordY; obsVect[count_terrain][2] = observ[i].coordZ - mean; Q[count_terrain] = 1; /* Q=I */ count_terrain++; } lineVect[i] = observ[i].lineID; obsVect_all[i][0] = observ[i].coordX; obsVect_all[i][1] = observ[i].coordY; obsVect_all[i][2] = observ[i].coordZ - mean; } G_free(observ); G_verbose_message(_("Bilinear interpolation")); normalDefBilin(N, TN, Q, obsVect, stepE, stepN, nsplx, nsply, elaboration_reg.west, elaboration_reg.south, nterrain, nparameters, BW); nCorrectGrad(N, lambda, nsplx, nsply, stepE, stepN); G_math_solver_cholesky_sband(N, parVect, TN, nparameters, BW); G_free_matrix(N); G_free_vector(TN); G_free_vector(Q); G_free_matrix(obsVect); G_verbose_message( _("Correction and creation of terrain vector")); P_Sparse_Correction(&In, &Out, &Terrain, &elaboration_reg, general_box, overlap_box, obsVect_all, lcat, parVect, lineVect, stepN, stepE, dims.overlap, HighThresh, LowThresh, nsplx, nsply, npoints, driver, mean, table_name); G_free_vector(parVect); G_free_matrix(obsVect_all); G_free_ivector(lineVect); } else { G_free(observ); G_warning(_("No data within this subregion. " "Consider changing the spline step.")); } G_free(lcat); } /*! END WHILE; last_column = TRUE */ } /*! END WHILE; last_row = TRUE */ /* Dropping auxiliar table */ if (npoints > 0) { G_debug(1, _("Dropping <%s>"), table_name); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Auxiliar table could not be dropped")); } db_close_database_shutdown_driver(driver); Vect_close(&In); Vect_close(&Out); Vect_close(&Terrain); G_done_msg(" "); exit(EXIT_SUCCESS); } /*! END MAIN */
/*--------------------------------------------------------------------*/ int main(int argc, char *argv[]) { /* Variable declarations */ int nsply, nsplx, nrows, ncols, nsplx_adj, nsply_adj; int nsubregion_col, nsubregion_row, subregion_row, subregion_col; int subregion = 0, nsubregions = 0; int last_row, last_column, grid, bilin, ext, flag_auxiliar, cross; /* booleans */ double stepN, stepE, lambda, mean; double N_extension, E_extension, edgeE, edgeN; const char *mapset, *drv, *db, *vector, *map; char table_name[GNAME_MAX], title[64]; char xname[GNAME_MAX], xmapset[GMAPSET_MAX]; int dim_vect, nparameters, BW; int *lineVect; /* Vector restoring primitive's ID */ double *TN, *Q, *parVect; /* Interpolating and least-square vectors */ double **N, **obsVect; /* Interpolation and least-square matrix */ SEGMENT out_seg, mask_seg; const char *out_file, *mask_file; int out_fd, mask_fd; double seg_size; int seg_mb, segments_in_memory; int have_mask; /* Structs declarations */ int raster; struct Map_info In, In_ext, Out; struct History history; struct GModule *module; struct Option *in_opt, *in_ext_opt, *out_opt, *out_map_opt, *stepE_opt, *stepN_opt, *lambda_f_opt, *type_opt, *dfield_opt, *col_opt, *mask_opt, *memory_opt, *solver, *error, *iter; struct Flag *cross_corr_flag, *spline_step_flag; struct Reg_dimens dims; struct Cell_head elaboration_reg, original_reg; struct bound_box general_box, overlap_box, original_box; struct Point *observ; struct line_cats *Cats; dbCatValArray cvarr; int with_z; int nrec, ctype = 0; struct field_info *Fi; dbDriver *driver, *driver_cats; /*----------------------------------------------------------------*/ /* Options declarations */ module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("surface")); G_add_keyword(_("interpolation")); G_add_keyword(_("LIDAR")); module->description = _("Performs bicubic or bilinear spline interpolation with Tykhonov regularization."); cross_corr_flag = G_define_flag(); cross_corr_flag->key = 'c'; cross_corr_flag->description = _("Find the best Tykhonov regularizing parameter using a \"leave-one-out\" cross validation method"); spline_step_flag = G_define_flag(); spline_step_flag->key = 'e'; spline_step_flag->label = _("Estimate point density and distance"); spline_step_flag->description = _("Estimate point density and distance for the input vector points within the current region extends and quit"); in_opt = G_define_standard_option(G_OPT_V_INPUT); in_opt->label = _("Name of input vector point map"); dfield_opt = G_define_standard_option(G_OPT_V_FIELD); dfield_opt->guisection = _("Settings"); col_opt = G_define_standard_option(G_OPT_DB_COLUMN); col_opt->required = NO; col_opt->label = _("Name of the attribute column with values to be used for approximation"); col_opt->description = _("If not given and input is 3D vector map then z-coordinates are used."); col_opt->guisection = _("Settings"); in_ext_opt = G_define_standard_option(G_OPT_V_INPUT); in_ext_opt->key = "sparse_input"; in_ext_opt->required = NO; in_ext_opt->label = _("Name of input vector map with sparse points"); out_opt = G_define_standard_option(G_OPT_V_OUTPUT); out_opt->required = NO; out_opt->guisection = _("Outputs"); out_map_opt = G_define_standard_option(G_OPT_R_OUTPUT); out_map_opt->key = "raster_output"; out_map_opt->required = NO; out_map_opt->guisection = _("Outputs"); mask_opt = G_define_standard_option(G_OPT_R_INPUT); mask_opt->key = "mask"; mask_opt->label = _("Raster map to use for masking (applies to raster output only)"); mask_opt->description = _("Only cells that are not NULL and not zero are interpolated"); mask_opt->required = NO; stepE_opt = G_define_option(); stepE_opt->key = "ew_step"; stepE_opt->type = TYPE_DOUBLE; stepE_opt->required = NO; stepE_opt->answer = "4"; stepE_opt->description = _("Length of each spline step in the east-west direction"); stepE_opt->guisection = _("Settings"); stepN_opt = G_define_option(); stepN_opt->key = "ns_step"; stepN_opt->type = TYPE_DOUBLE; stepN_opt->required = NO; stepN_opt->answer = "4"; stepN_opt->description = _("Length of each spline step in the north-south direction"); stepN_opt->guisection = _("Settings"); type_opt = G_define_option(); type_opt->key = "method"; type_opt->description = _("Spline interpolation algorithm"); type_opt->type = TYPE_STRING; type_opt->options = "bilinear,bicubic"; type_opt->answer = "bilinear"; type_opt->guisection = _("Settings"); G_asprintf((char **) &(type_opt->descriptions), "bilinear;%s;bicubic;%s", _("Bilinear interpolation"), _("Bicubic interpolation")); lambda_f_opt = G_define_option(); lambda_f_opt->key = "lambda_i"; lambda_f_opt->type = TYPE_DOUBLE; lambda_f_opt->required = NO; lambda_f_opt->description = _("Tykhonov regularization parameter (affects smoothing)"); lambda_f_opt->answer = "0.01"; lambda_f_opt->guisection = _("Settings"); solver = N_define_standard_option(N_OPT_SOLVER_SYMM); solver->options = "cholesky,cg"; solver->answer = "cholesky"; iter = N_define_standard_option(N_OPT_MAX_ITERATIONS); error = N_define_standard_option(N_OPT_ITERATION_ERROR); memory_opt = G_define_option(); memory_opt->key = "memory"; memory_opt->type = TYPE_INTEGER; memory_opt->required = NO; memory_opt->answer = "300"; memory_opt->label = _("Maximum memory to be used (in MB)"); memory_opt->description = _("Cache size for raster rows"); /*----------------------------------------------------------------*/ /* Parsing */ G_gisinit(argv[0]); if (G_parser(argc, argv)) exit(EXIT_FAILURE); vector = out_opt->answer; map = out_map_opt->answer; if (vector && map) G_fatal_error(_("Choose either vector or raster output, not both")); if (!vector && !map && !cross_corr_flag->answer) G_fatal_error(_("No raster or vector or cross-validation output")); if (!strcmp(type_opt->answer, "linear")) bilin = P_BILINEAR; else bilin = P_BICUBIC; stepN = atof(stepN_opt->answer); stepE = atof(stepE_opt->answer); lambda = atof(lambda_f_opt->answer); flag_auxiliar = FALSE; drv = db_get_default_driver_name(); if (!drv) { if (db_set_default_connection() != DB_OK) G_fatal_error(_("Unable to set default DB connection")); drv = db_get_default_driver_name(); } db = db_get_default_database_name(); if (!db) G_fatal_error(_("No default DB defined")); /* Set auxiliary table's name */ if (vector) { if (G_name_is_fully_qualified(out_opt->answer, xname, xmapset)) { sprintf(table_name, "%s_aux", xname); } else sprintf(table_name, "%s_aux", out_opt->answer); } /* Something went wrong in a previous v.surf.bspline execution */ if (db_table_exists(drv, db, table_name)) { /* Start driver and open db */ driver = db_start_driver_open_database(drv, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), drv); db_set_error_handler_driver(driver); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Old auxiliary table could not be dropped")); db_close_database_shutdown_driver(driver); } /* Open input vector */ if ((mapset = G_find_vector2(in_opt->answer, "")) == NULL) G_fatal_error(_("Vector map <%s> not found"), in_opt->answer); Vect_set_open_level(1); /* WITHOUT TOPOLOGY */ if (1 > Vect_open_old(&In, in_opt->answer, mapset)) G_fatal_error(_("Unable to open vector map <%s> at the topological level"), in_opt->answer); bspline_field = 0; /* assume 3D input */ bspline_column = col_opt->answer; with_z = !bspline_column && Vect_is_3d(&In); if (Vect_is_3d(&In)) { if (!with_z) G_verbose_message(_("Input is 3D: using attribute values instead of z-coordinates for approximation")); else G_verbose_message(_("Input is 3D: using z-coordinates for approximation")); } else { /* 2D */ if (!bspline_column) G_fatal_error(_("Input vector map is 2D. Parameter <%s> required."), col_opt->key); } if (!with_z) { bspline_field = Vect_get_field_number(&In, dfield_opt->answer); } /* Estimate point density and mean distance for current region */ if (spline_step_flag->answer) { double dens, dist; if (P_estimate_splinestep(&In, &dens, &dist) == 0) { fprintf(stdout, _("Estimated point density: %.4g"), dens); fprintf(stdout, _("Estimated mean distance between points: %.4g"), dist); } else { fprintf(stdout, _("No points in current region")); } Vect_close(&In); exit(EXIT_SUCCESS); } /*----------------------------------------------------------------*/ /* Cross-correlation begins */ if (cross_corr_flag->answer) { G_debug(1, "CrossCorrelation()"); cross = cross_correlation(&In, stepE, stepN); if (cross != TRUE) G_fatal_error(_("Cross validation didn't finish correctly")); else { G_debug(1, "Cross validation finished correctly"); Vect_close(&In); G_done_msg(_("Cross validation finished for ew_step = %f and ns_step = %f"), stepE, stepN); exit(EXIT_SUCCESS); } } /* Open input ext vector */ ext = FALSE; if (in_ext_opt->answer) { ext = TRUE; G_message(_("Vector map <%s> of sparse points will be interpolated"), in_ext_opt->answer); if ((mapset = G_find_vector2(in_ext_opt->answer, "")) == NULL) G_fatal_error(_("Vector map <%s> not found"), in_ext_opt->answer); Vect_set_open_level(1); /* WITHOUT TOPOLOGY */ if (1 > Vect_open_old(&In_ext, in_ext_opt->answer, mapset)) G_fatal_error(_("Unable to open vector map <%s> at the topological level"), in_opt->answer); } /* Open output map */ /* vector output */ if (vector && !map) { if (strcmp(drv, "dbf") == 0) G_fatal_error(_("Sorry, the <%s> driver is not compatible with " "the vector output of this module. " "Try with raster output or another driver."), drv); Vect_check_input_output_name(in_opt->answer, out_opt->answer, G_FATAL_EXIT); grid = FALSE; if (0 > Vect_open_new(&Out, out_opt->answer, WITH_Z)) G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); /* Copy vector Head File */ if (ext == FALSE) { Vect_copy_head_data(&In, &Out); Vect_hist_copy(&In, &Out); } else { Vect_copy_head_data(&In_ext, &Out); Vect_hist_copy(&In_ext, &Out); } Vect_hist_command(&Out); G_verbose_message(_("Points in input vector map <%s> will be interpolated"), vector); } /* read z values from attribute table */ if (bspline_field > 0) { G_message(_("Reading values from attribute table...")); db_CatValArray_init(&cvarr); Fi = Vect_get_field(&In, bspline_field); if (Fi == NULL) G_fatal_error(_("Cannot read layer info")); driver_cats = db_start_driver_open_database(Fi->driver, Fi->database); /*G_debug (0, _("driver=%s db=%s"), Fi->driver, Fi->database); */ if (driver_cats == NULL) G_fatal_error(_("Unable to open database <%s> by driver <%s>"), Fi->database, Fi->driver); db_set_error_handler_driver(driver_cats); nrec = db_select_CatValArray(driver_cats, Fi->table, Fi->key, col_opt->answer, NULL, &cvarr); G_debug(3, "nrec = %d", nrec); ctype = cvarr.ctype; if (ctype != DB_C_TYPE_INT && ctype != DB_C_TYPE_DOUBLE) G_fatal_error(_("Column type not supported")); if (nrec < 0) G_fatal_error(_("Unable to select data from table")); G_verbose_message(_("%d records selected from table"), nrec); db_close_database_shutdown_driver(driver_cats); } /*----------------------------------------------------------------*/ /* Interpolation begins */ G_debug(1, "Interpolation()"); /* Open driver and database */ driver = db_start_driver_open_database(drv, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. " "Run db.connect."), drv); db_set_error_handler_driver(driver); /* Create auxiliary table */ if (vector) { if ((flag_auxiliar = P_Create_Aux4_Table(driver, table_name)) == FALSE) { P_Drop_Aux_Table(driver, table_name); G_fatal_error(_("Interpolation: Creating table: " "It was impossible to create table <%s>."), table_name); } /* db_create_index2(driver, table_name, "ID"); */ /* sqlite likes that ??? */ db_close_database_shutdown_driver(driver); driver = db_start_driver_open_database(drv, db); } /* raster output */ raster = -1; Rast_set_fp_type(DCELL_TYPE); if (!vector && map) { grid = TRUE; raster = Rast_open_fp_new(out_map_opt->answer); G_verbose_message(_("Cells for raster map <%s> will be interpolated"), map); } /* Setting regions and boxes */ G_debug(1, "Interpolation: Setting regions and boxes"); G_get_window(&original_reg); G_get_window(&elaboration_reg); Vect_region_box(&original_reg, &original_box); Vect_region_box(&elaboration_reg, &overlap_box); Vect_region_box(&elaboration_reg, &general_box); nrows = Rast_window_rows(); ncols = Rast_window_cols(); /* Alloc raster matrix */ have_mask = 0; out_file = mask_file = NULL; out_fd = mask_fd = -1; if (grid == TRUE) { int row; DCELL *drastbuf; seg_mb = atoi(memory_opt->answer); if (seg_mb < 3) G_fatal_error(_("Memory in MB must be >= 3")); if (mask_opt->answer) seg_size = sizeof(double) + sizeof(char); else seg_size = sizeof(double); seg_size = (seg_size * SEGSIZE * SEGSIZE) / (1 << 20); segments_in_memory = seg_mb / seg_size + 0.5; G_debug(1, "%d %dx%d segments held in memory", segments_in_memory, SEGSIZE, SEGSIZE); out_file = G_tempfile(); out_fd = creat(out_file, 0666); if (Segment_format(out_fd, nrows, ncols, SEGSIZE, SEGSIZE, sizeof(double)) != 1) G_fatal_error(_("Can not create temporary file")); close(out_fd); out_fd = open(out_file, 2); if (Segment_init(&out_seg, out_fd, segments_in_memory) != 1) G_fatal_error(_("Can not initialize temporary file")); /* initialize output */ G_message(_("Initializing output...")); drastbuf = Rast_allocate_buf(DCELL_TYPE); Rast_set_d_null_value(drastbuf, ncols); for (row = 0; row < nrows; row++) { G_percent(row, nrows, 2); Segment_put_row(&out_seg, drastbuf, row); } G_percent(row, nrows, 2); if (mask_opt->answer) { int row, col, maskfd; DCELL dval, *drastbuf; char mask_val; G_message(_("Load masking map")); mask_file = G_tempfile(); mask_fd = creat(mask_file, 0666); if (Segment_format(mask_fd, nrows, ncols, SEGSIZE, SEGSIZE, sizeof(char)) != 1) G_fatal_error(_("Can not create temporary file")); close(mask_fd); mask_fd = open(mask_file, 2); if (Segment_init(&mask_seg, mask_fd, segments_in_memory) != 1) G_fatal_error(_("Can not initialize temporary file")); maskfd = Rast_open_old(mask_opt->answer, ""); drastbuf = Rast_allocate_buf(DCELL_TYPE); for (row = 0; row < nrows; row++) { G_percent(row, nrows, 2); Rast_get_d_row(maskfd, drastbuf, row); for (col = 0; col < ncols; col++) { dval = drastbuf[col]; if (Rast_is_d_null_value(&dval) || dval == 0) mask_val = 0; else mask_val = 1; Segment_put(&mask_seg, &mask_val, row, col); } } G_percent(row, nrows, 2); G_free(drastbuf); Rast_close(maskfd); have_mask = 1; } } /*------------------------------------------------------------------ | Subdividing and working with tiles: | Each original region will be divided into several subregions. | Each one will be overlaped by its neighbouring subregions. | The overlapping is calculated as a fixed OVERLAP_SIZE times | the largest spline step plus 2 * edge ----------------------------------------------------------------*/ /* Fixing parameters of the elaboration region */ P_zero_dim(&dims); /* Set dim struct to zero */ nsplx_adj = NSPLX_MAX; nsply_adj = NSPLY_MAX; if (stepN > stepE) dims.overlap = OVERLAP_SIZE * stepN; else dims.overlap = OVERLAP_SIZE * stepE; P_get_edge(bilin, &dims, stepE, stepN); P_set_dim(&dims, stepE, stepN, &nsplx_adj, &nsply_adj); G_verbose_message(_("Adjusted EW splines %d"), nsplx_adj); G_verbose_message(_("Adjusted NS splines %d"), nsply_adj); /* calculate number of subregions */ edgeE = dims.ew_size - dims.overlap - 2 * dims.edge_v; edgeN = dims.sn_size - dims.overlap - 2 * dims.edge_h; N_extension = original_reg.north - original_reg.south; E_extension = original_reg.east - original_reg.west; nsubregion_col = ceil(E_extension / edgeE) + 0.5; nsubregion_row = ceil(N_extension / edgeN) + 0.5; if (nsubregion_col < 0) nsubregion_col = 0; if (nsubregion_row < 0) nsubregion_row = 0; nsubregions = nsubregion_row * nsubregion_col; /* Creating line and categories structs */ Cats = Vect_new_cats_struct(); Vect_cat_set(Cats, 1, 0); subregion_row = 0; elaboration_reg.south = original_reg.north; last_row = FALSE; while (last_row == FALSE) { /* For each subregion row */ subregion_row++; P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_ROW); if (elaboration_reg.north > original_reg.north) { /* First row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_ROW); } if (elaboration_reg.south <= original_reg.south) { /* Last row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_ROW); last_row = TRUE; } nsply = ceil((elaboration_reg.north - elaboration_reg.south) / stepN) + 0.5; G_debug(1, "Interpolation: nsply = %d", nsply); /* if (nsply > NSPLY_MAX) nsply = NSPLY_MAX; */ elaboration_reg.east = original_reg.west; last_column = FALSE; subregion_col = 0; /* TODO: process each subregion using its own thread (via OpenMP or pthreads) */ /* I'm not sure about pthreads, but you can tell OpenMP to start all at the same time and it will keep num_workers supplied with the next job as free cpus become available */ while (last_column == FALSE) { /* For each subregion column */ int npoints = 0; /* needed for sparse points interpolation */ int npoints_ext, *lineVect_ext = NULL; double **obsVect_ext; /*, mean_ext = .0; */ struct Point *observ_ext; subregion_col++; subregion++; if (nsubregions > 1) G_message(_("Processing subregion %d of %d..."), subregion, nsubregions); P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_COLUMN); if (elaboration_reg.west < original_reg.west) { /* First column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_COLUMN); } if (elaboration_reg.east >= original_reg.east) { /* Last column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_COLUMN); last_column = TRUE; } nsplx = ceil((elaboration_reg.east - elaboration_reg.west) / stepE) + 0.5; G_debug(1, "Interpolation: nsplx = %d", nsplx); /* if (nsplx > NSPLX_MAX) nsplx = NSPLX_MAX; */ G_debug(1, "Interpolation: (%d,%d): subregion bounds", subregion_row, subregion_col); G_debug(1, "Interpolation: \t\tNORTH:%.2f\t", elaboration_reg.north); G_debug(1, "Interpolation: WEST:%.2f\t\tEAST:%.2f", elaboration_reg.west, elaboration_reg.east); G_debug(1, "Interpolation: \t\tSOUTH:%.2f", elaboration_reg.south); #ifdef DEBUG_SUBREGIONS fprintf(stdout, "B 5\n"); fprintf(stdout, " %.11g %.11g\n", elaboration_reg.east, elaboration_reg.north); fprintf(stdout, " %.11g %.11g\n", elaboration_reg.west, elaboration_reg.north); fprintf(stdout, " %.11g %.11g\n", elaboration_reg.west, elaboration_reg.south); fprintf(stdout, " %.11g %.11g\n", elaboration_reg.east, elaboration_reg.south); fprintf(stdout, " %.11g %.11g\n", elaboration_reg.east, elaboration_reg.north); fprintf(stdout, "C 1 1\n"); fprintf(stdout, " %.11g %.11g\n", (elaboration_reg.west + elaboration_reg.east) / 2, (elaboration_reg.south + elaboration_reg.north) / 2); fprintf(stdout, " 1 %d\n", subregion); #endif /* reading points in interpolation region */ dim_vect = nsplx * nsply; observ_ext = NULL; if (grid == FALSE && ext == TRUE) { observ_ext = P_Read_Vector_Region_Map(&In_ext, &elaboration_reg, &npoints_ext, dim_vect, 1); } else npoints_ext = 1; if (grid == TRUE && have_mask) { /* any unmasked cells in general region ? */ mean = 0; observ_ext = P_Read_Raster_Region_masked(&mask_seg, &original_reg, original_box, general_box, &npoints_ext, dim_vect, mean); } observ = NULL; if (npoints_ext > 0) { observ = P_Read_Vector_Region_Map(&In, &elaboration_reg, &npoints, dim_vect, bspline_field); } else npoints = 1; G_debug(1, "Interpolation: (%d,%d): Number of points in <elaboration_box> is %d", subregion_row, subregion_col, npoints); if (npoints > 0) G_verbose_message(_("%d points found in this subregion"), npoints); /* only interpolate if there are any points in current subregion */ if (npoints > 0 && npoints_ext > 0) { int i; nparameters = nsplx * nsply; BW = P_get_BandWidth(bilin, nsply); /* Least Squares system */ N = G_alloc_matrix(nparameters, BW); /* Normal matrix */ TN = G_alloc_vector(nparameters); /* vector */ parVect = G_alloc_vector(nparameters); /* Parameters vector */ obsVect = G_alloc_matrix(npoints, 3); /* Observation vector */ Q = G_alloc_vector(npoints); /* "a priori" var-cov matrix */ lineVect = G_alloc_ivector(npoints); /* */ for (i = 0; i < npoints; i++) { /* Setting obsVect vector & Q matrix */ double dval; Q[i] = 1; /* Q=I */ lineVect[i] = observ[i].lineID; obsVect[i][0] = observ[i].coordX; obsVect[i][1] = observ[i].coordY; /* read z coordinates from attribute table */ if (bspline_field > 0) { int cat, ival, ret; cat = observ[i].cat; if (cat < 0) continue; if (ctype == DB_C_TYPE_INT) { ret = db_CatValArray_get_value_int(&cvarr, cat, &ival); obsVect[i][2] = ival; observ[i].coordZ = ival; } else { /* DB_C_TYPE_DOUBLE */ ret = db_CatValArray_get_value_double(&cvarr, cat, &dval); obsVect[i][2] = dval; observ[i].coordZ = dval; } if (ret != DB_OK) { G_warning(_("Interpolation: (%d,%d): No record for point (cat = %d)"), subregion_row, subregion_col, cat); continue; } } /* use z coordinates of 3D vector */ else { obsVect[i][2] = observ[i].coordZ; } } /* Mean calculation for every point */ mean = P_Mean_Calc(&elaboration_reg, observ, npoints); G_debug(1, "Interpolation: (%d,%d): mean=%lf", subregion_row, subregion_col, mean); G_free(observ); for (i = 0; i < npoints; i++) obsVect[i][2] -= mean; /* Bilinear interpolation */ if (bilin) { G_debug(1, "Interpolation: (%d,%d): Bilinear interpolation...", subregion_row, subregion_col); normalDefBilin(N, TN, Q, obsVect, stepE, stepN, nsplx, nsply, elaboration_reg.west, elaboration_reg.south, npoints, nparameters, BW); nCorrectGrad(N, lambda, nsplx, nsply, stepE, stepN); } /* Bicubic interpolation */ else { G_debug(1, "Interpolation: (%d,%d): Bicubic interpolation...", subregion_row, subregion_col); normalDefBicubic(N, TN, Q, obsVect, stepE, stepN, nsplx, nsply, elaboration_reg.west, elaboration_reg.south, npoints, nparameters, BW); nCorrectGrad(N, lambda, nsplx, nsply, stepE, stepN); } if(G_strncasecmp(solver->answer, "cg", 2) == 0) G_math_solver_cg_sband(N, parVect, TN, nparameters, BW, atoi(iter->answer), atof(error->answer)); else G_math_solver_cholesky_sband(N, parVect, TN, nparameters, BW); G_free_matrix(N); G_free_vector(TN); G_free_vector(Q); if (grid == TRUE) { /* GRID INTERPOLATION ==> INTERPOLATION INTO A RASTER */ G_debug(1, "Interpolation: (%d,%d): Regular_Points...", subregion_row, subregion_col); if (!have_mask) { P_Regular_Points(&elaboration_reg, &original_reg, general_box, overlap_box, &out_seg, parVect, stepN, stepE, dims.overlap, mean, nsplx, nsply, nrows, ncols, bilin); } else { P_Sparse_Raster_Points(&out_seg, &elaboration_reg, &original_reg, general_box, overlap_box, observ_ext, parVect, stepE, stepN, dims.overlap, nsplx, nsply, npoints_ext, bilin, mean); } } else { /* OBSERVATION POINTS INTERPOLATION */ if (ext == FALSE) { G_debug(1, "Interpolation: (%d,%d): Sparse_Points...", subregion_row, subregion_col); P_Sparse_Points(&Out, &elaboration_reg, general_box, overlap_box, obsVect, parVect, lineVect, stepE, stepN, dims.overlap, nsplx, nsply, npoints, bilin, Cats, driver, mean, table_name); } else { /* FLAG_EXT == TRUE */ /* done that earlier */ /* int npoints_ext, *lineVect_ext = NULL; double **obsVect_ext; struct Point *observ_ext; observ_ext = P_Read_Vector_Region_Map(&In_ext, &elaboration_reg, &npoints_ext, dim_vect, 1); */ obsVect_ext = G_alloc_matrix(npoints_ext, 3); /* Observation vector_ext */ lineVect_ext = G_alloc_ivector(npoints_ext); for (i = 0; i < npoints_ext; i++) { /* Setting obsVect_ext vector & Q matrix */ obsVect_ext[i][0] = observ_ext[i].coordX; obsVect_ext[i][1] = observ_ext[i].coordY; obsVect_ext[i][2] = observ_ext[i].coordZ - mean; lineVect_ext[i] = observ_ext[i].lineID; } G_free(observ_ext); G_debug(1, "Interpolation: (%d,%d): Sparse_Points...", subregion_row, subregion_col); P_Sparse_Points(&Out, &elaboration_reg, general_box, overlap_box, obsVect_ext, parVect, lineVect_ext, stepE, stepN, dims.overlap, nsplx, nsply, npoints_ext, bilin, Cats, driver, mean, table_name); G_free_matrix(obsVect_ext); G_free_ivector(lineVect_ext); } /* END FLAG_EXT == TRUE */ } /* END GRID == FALSE */ G_free_vector(parVect); G_free_matrix(obsVect); G_free_ivector(lineVect); } else { if (observ) G_free(observ); if (observ_ext) G_free(observ_ext); if (npoints == 0) G_warning(_("No data within this subregion. " "Consider increasing spline step values.")); } } /*! END WHILE; last_column = TRUE */ } /*! END WHILE; last_row = TRUE */ G_verbose_message(_("Writing output...")); /* Writing the output raster map */ if (grid == TRUE) { int row, col; DCELL *drastbuf, dval; if (have_mask) { Segment_release(&mask_seg); /* release memory */ close(mask_fd); unlink(mask_file); } drastbuf = Rast_allocate_buf(DCELL_TYPE); for (row = 0; row < nrows; row++) { G_percent(row, nrows, 2); for (col = 0; col < ncols; col++) { Segment_get(&out_seg, &dval, row, col); drastbuf[col] = dval; } Rast_put_d_row(raster, drastbuf); } Rast_close(raster); Segment_release(&out_seg); /* release memory */ close(out_fd); unlink(out_file); /* set map title */ sprintf(title, "%s interpolation with Tykhonov regularization", type_opt->answer); Rast_put_cell_title(out_map_opt->answer, title); /* write map history */ Rast_short_history(out_map_opt->answer, "raster", &history); Rast_command_history(&history); Rast_write_history(out_map_opt->answer, &history); } /* Writing to the output vector map the points from the overlapping zones */ else if (flag_auxiliar == TRUE) { if (ext == FALSE) P_Aux_to_Vector(&In, &Out, driver, table_name); else P_Aux_to_Vector(&In_ext, &Out, driver, table_name); /* Drop auxiliary table */ G_debug(1, "%s: Dropping <%s>", argv[0], table_name); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Auxiliary table could not be dropped")); } db_close_database_shutdown_driver(driver); Vect_close(&In); if (ext != FALSE) Vect_close(&In_ext); if (vector) Vect_close(&Out); G_done_msg(" "); exit(EXIT_SUCCESS); } /*END MAIN */
struct Point *P_Read_Raster_Region_Map(SEGMENT *in_seg, struct Cell_head *Elaboration, struct Cell_head *Original, int *num_points, int dim_vect) { int col, row, startcol, endcol, startrow, endrow, nrows, ncols; int pippo, npoints; double x, y, z; struct Point *obs; struct bound_box elaboration_box; pippo = dim_vect; obs = (struct Point *)G_calloc(pippo, sizeof(struct Point)); /* Reading points inside elaboration zone */ Vect_region_box(Elaboration, &elaboration_box); npoints = 0; nrows = Original->rows; ncols = Original->cols; if (Original->north > Elaboration->north) startrow = (Original->north - Elaboration->north) / Original->ns_res - 1; else startrow = 0; if (Original->north > Elaboration->south) { endrow = (Original->north - Elaboration->south) / Original->ns_res + 1; if (endrow > nrows) endrow = nrows; } else endrow = nrows; if (Elaboration->west > Original->west) startcol = (Elaboration->west - Original->west) / Original->ew_res - 1; else startcol = 0; if (Elaboration->east > Original->west) { endcol = (Elaboration->east - Original->west) / Original->ew_res + 1; if (endcol > ncols) endcol = ncols; } else endcol = ncols; for (row = startrow; row < endrow; row++) { for (col = startcol; col < endcol; col++) { Segment_get(in_seg, &z, row, col); if (!Rast_is_d_null_value(&z)) { x = Rast_col_to_easting((double)(col) + 0.5, Original); y = Rast_row_to_northing((double)(row) + 0.5, Original); if (Vect_point_in_box(x, y, 0, &elaboration_box)) { npoints++; if (npoints >= pippo) { pippo += dim_vect; obs = (struct Point *)G_realloc((void *)obs, (signed int)pippo * sizeof(struct Point)); } /* Storing observation vector */ obs[npoints - 1].coordX = x; obs[npoints - 1].coordY = y; obs[npoints - 1].coordZ = z; } } } } *num_points = npoints; return obs; }
struct Point *P_Read_Vector_Region_Map(struct Map_info *Map, struct Cell_head *Elaboration, int *num_points, int dim_vect, int layer) { int line_num, pippo, npoints, cat, type; double x, y, z; struct Point *obs; struct line_pnts *points; struct line_cats *categories; struct bound_box elaboration_box; pippo = dim_vect; obs = (struct Point *)G_calloc(pippo, sizeof(struct Point)); points = Vect_new_line_struct(); categories = Vect_new_cats_struct(); /* Reading points inside elaboration zone */ Vect_region_box(Elaboration, &elaboration_box); npoints = 0; line_num = 0; Vect_rewind(Map); while ((type = Vect_read_next_line(Map, points, categories)) > 0) { if (!(type & GV_POINT)) continue; line_num++; x = points->x[0]; y = points->y[0]; if (points->z != NULL) z = points->z[0]; else z = 0.0; /* Reading and storing points only if in elaboration_reg */ if (Vect_point_in_box(x, y, z, &elaboration_box)) { npoints++; if (npoints >= pippo) { pippo += dim_vect; obs = (struct Point *)G_realloc((void *)obs, (signed int)pippo * sizeof(struct Point)); } /* Storing observation vector */ obs[npoints - 1].coordX = x; obs[npoints - 1].coordY = y; obs[npoints - 1].coordZ = z; obs[npoints - 1].lineID = line_num; /* Storing also the line's number */ Vect_cat_get(categories, layer, &cat); obs[npoints - 1].cat = cat; } } Vect_destroy_line_struct(points); Vect_destroy_cats_struct(categories); *num_points = npoints; return obs; }
/*--------------------------------------------------------------------*/ int main(int argc, char *argv[]) { /* Variables declarations */ int nsplx_adj, nsply_adj; int nsubregion_col, nsubregion_row; int subregion = 0, nsubregions = 0; double N_extension, E_extension, edgeE, edgeN; int dim_vect, nparameters, BW, npoints; double mean, lambda; const char *dvr, *db, *mapset; char table_name[GNAME_MAX]; char xname[GNAME_MAX], xmapset[GMAPSET_MAX]; int last_row, last_column, flag_auxiliar = FALSE; int filter_mode; int *lineVect; double *TN, *Q, *parVect; /* Interpolating and least-square vectors */ double **N, **obsVect; /* Interpolation and least-square matrix */ /* Structs declarations */ struct Map_info In, Out, Outlier, Qgis; struct Option *in_opt, *out_opt, *outlier_opt, *qgis_opt, *stepE_opt, *stepN_opt, *lambda_f_opt, *Thres_O_opt, *filter_opt; struct Flag *spline_step_flag; struct GModule *module; struct Reg_dimens dims; struct Cell_head elaboration_reg, original_reg; struct bound_box general_box, overlap_box; struct Point *observ; dbDriver *driver; /*----------------------------------------------------------------*/ /* Options declaration */ module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("statistics")); G_add_keyword(_("extract")); G_add_keyword(_("select")); G_add_keyword(_("filter")); module->description = _("Removes outliers from vector point data."); spline_step_flag = G_define_flag(); spline_step_flag->key = 'e'; spline_step_flag->label = _("Estimate point density and distance"); spline_step_flag->description = _("Estimate point density and distance for the input vector points within the current region extends and quit"); in_opt = G_define_standard_option(G_OPT_V_INPUT); out_opt = G_define_standard_option(G_OPT_V_OUTPUT); outlier_opt = G_define_option(); outlier_opt->key = "outlier"; outlier_opt->type = TYPE_STRING; outlier_opt->key_desc = "name"; outlier_opt->required = YES; outlier_opt->gisprompt = "new,vector,vector"; outlier_opt->description = _("Name of output outlier vector map"); qgis_opt = G_define_option(); qgis_opt->key = "qgis"; qgis_opt->type = TYPE_STRING; qgis_opt->key_desc = "name"; qgis_opt->required = NO; qgis_opt->gisprompt = "new,vector,vector"; qgis_opt->description = _("Name of vector map for visualization in QGIS"); stepE_opt = G_define_option(); stepE_opt->key = "ew_step"; stepE_opt->type = TYPE_DOUBLE; stepE_opt->required = NO; stepE_opt->answer = "10"; stepE_opt->description = _("Length of each spline step in the east-west direction"); stepE_opt->guisection = _("Settings"); stepN_opt = G_define_option(); stepN_opt->key = "ns_step"; stepN_opt->type = TYPE_DOUBLE; stepN_opt->required = NO; stepN_opt->answer = "10"; stepN_opt->description = _("Length of each spline step in the north-south direction"); stepN_opt->guisection = _("Settings"); lambda_f_opt = G_define_option(); lambda_f_opt->key = "lambda"; lambda_f_opt->type = TYPE_DOUBLE; lambda_f_opt->required = NO; lambda_f_opt->description = _("Tykhonov regularization weight"); lambda_f_opt->answer = "0.1"; lambda_f_opt->guisection = _("Settings"); Thres_O_opt = G_define_option(); Thres_O_opt->key = "threshold"; Thres_O_opt->type = TYPE_DOUBLE; Thres_O_opt->required = NO; Thres_O_opt->description = _("Threshold for the outliers"); Thres_O_opt->answer = "50"; filter_opt = G_define_option(); filter_opt->key = "filter"; filter_opt->type = TYPE_STRING; filter_opt->required = NO; filter_opt->description = _("Filtering option"); filter_opt->options = "both,positive,negative"; filter_opt->answer = "both"; /* Parsing */ G_gisinit(argv[0]); if (G_parser(argc, argv)) exit(EXIT_FAILURE); if (!(db = G_getenv_nofatal2("DB_DATABASE", G_VAR_MAPSET))) G_fatal_error(_("Unable to read name of database")); if (!(dvr = G_getenv_nofatal2("DB_DRIVER", G_VAR_MAPSET))) G_fatal_error(_("Unable to read name of driver")); stepN = atof(stepN_opt->answer); stepE = atof(stepE_opt->answer); lambda = atof(lambda_f_opt->answer); Thres_Outlier = atof(Thres_O_opt->answer); filter_mode = 0; if (strcmp(filter_opt->answer, "positive") == 0) filter_mode = 1; else if (strcmp(filter_opt->answer, "negative") == 0) filter_mode = -1; P_set_outlier_fn(filter_mode); flag_auxiliar = FALSE; /* Checking vector names */ Vect_check_input_output_name(in_opt->answer, out_opt->answer, G_FATAL_EXIT); if ((mapset = G_find_vector2(in_opt->answer, "")) == NULL) { G_fatal_error(_("Vector map <%s> not found"), in_opt->answer); } /* Setting auxiliar table's name */ if (G_name_is_fully_qualified(out_opt->answer, xname, xmapset)) { sprintf(table_name, "%s_aux", xname); } else sprintf(table_name, "%s_aux", out_opt->answer); /* Something went wrong in a previous v.outlier execution */ if (db_table_exists(dvr, db, table_name)) { /* Start driver and open db */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); db_set_error_handler_driver(driver); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Old auxiliar table could not be dropped")); db_close_database_shutdown_driver(driver); } /* Open input vector */ Vect_set_open_level(1); /* WITHOUT TOPOLOGY */ if (1 > Vect_open_old(&In, in_opt->answer, mapset)) G_fatal_error(_("Unable to open vector map <%s> at the topological level"), in_opt->answer); /* Input vector must be 3D */ if (!Vect_is_3d(&In)) G_fatal_error(_("Input vector map <%s> is not 3D!"), in_opt->answer); /* Estimate point density and mean distance for current region */ if (spline_step_flag->answer) { double dens, dist; if (P_estimate_splinestep(&In, &dens, &dist) == 0) { G_message("Estimated point density: %.4g", dens); G_message("Estimated mean distance between points: %.4g", dist); } else G_warning(_("No points in current region!")); Vect_close(&In); exit(EXIT_SUCCESS); } /* Open output vector */ if (qgis_opt->answer) if (0 > Vect_open_new(&Qgis, qgis_opt->answer, WITHOUT_Z)) G_fatal_error(_("Unable to create vector map <%s>"), qgis_opt->answer); if (0 > Vect_open_new(&Out, out_opt->answer, WITH_Z)) { Vect_close(&Qgis); G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); } if (0 > Vect_open_new(&Outlier, outlier_opt->answer, WITH_Z)) { Vect_close(&Out); Vect_close(&Qgis); G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); } /* Copy vector Head File */ Vect_copy_head_data(&In, &Out); Vect_hist_copy(&In, &Out); Vect_hist_command(&Out); Vect_copy_head_data(&In, &Outlier); Vect_hist_copy(&In, &Outlier); Vect_hist_command(&Outlier); if (qgis_opt->answer) { Vect_copy_head_data(&In, &Qgis); Vect_hist_copy(&In, &Qgis); Vect_hist_command(&Qgis); } /* Open driver and database */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); db_set_error_handler_driver(driver); /* Create auxiliar table */ if ((flag_auxiliar = P_Create_Aux2_Table(driver, table_name)) == FALSE) G_fatal_error(_("It was impossible to create <%s> table."), table_name); db_create_index2(driver, table_name, "ID"); /* sqlite likes that ??? */ db_close_database_shutdown_driver(driver); driver = db_start_driver_open_database(dvr, db); /* Setting regions and boxes */ G_get_set_window(&original_reg); G_get_set_window(&elaboration_reg); Vect_region_box(&elaboration_reg, &overlap_box); Vect_region_box(&elaboration_reg, &general_box); /*------------------------------------------------------------------ | Subdividing and working with tiles: | Each original region will be divided into several subregions. | Each one will be overlaped by its neighbouring subregions. | The overlapping is calculated as a fixed OVERLAP_SIZE times | the largest spline step plus 2 * edge ----------------------------------------------------------------*/ /* Fixing parameters of the elaboration region */ P_zero_dim(&dims); /* Set dim struct to zero */ nsplx_adj = NSPLX_MAX; nsply_adj = NSPLY_MAX; if (stepN > stepE) dims.overlap = OVERLAP_SIZE * stepN; else dims.overlap = OVERLAP_SIZE * stepE; P_get_edge(P_BILINEAR, &dims, stepE, stepN); P_set_dim(&dims, stepE, stepN, &nsplx_adj, &nsply_adj); G_verbose_message(_("Adjusted EW splines %d"), nsplx_adj); G_verbose_message(_("Adjusted NS splines %d"), nsply_adj); /* calculate number of subregions */ edgeE = dims.ew_size - dims.overlap - 2 * dims.edge_v; edgeN = dims.sn_size - dims.overlap - 2 * dims.edge_h; N_extension = original_reg.north - original_reg.south; E_extension = original_reg.east - original_reg.west; nsubregion_col = ceil(E_extension / edgeE) + 0.5; nsubregion_row = ceil(N_extension / edgeN) + 0.5; if (nsubregion_col < 0) nsubregion_col = 0; if (nsubregion_row < 0) nsubregion_row = 0; nsubregions = nsubregion_row * nsubregion_col; elaboration_reg.south = original_reg.north; last_row = FALSE; while (last_row == FALSE) { /* For each row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_ROW); if (elaboration_reg.north > original_reg.north) { /* First row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_ROW); } if (elaboration_reg.south <= original_reg.south) { /* Last row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_ROW); last_row = TRUE; } nsply = ceil((elaboration_reg.north - elaboration_reg.south) / stepN) + 0.5; /* if (nsply > NSPLY_MAX) nsply = NSPLY_MAX; */ G_debug(1, "nsply = %d", nsply); elaboration_reg.east = original_reg.west; last_column = FALSE; while (last_column == FALSE) { /* For each column */ subregion++; if (nsubregions > 1) G_message(_("Processing subregion %d of %d..."), subregion, nsubregions); else /* v.outlier -e will report mean point distance: */ G_warning(_("No subregions found! Check values for 'ew_step' and 'ns_step' parameters")); P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_COLUMN); if (elaboration_reg.west < original_reg.west) { /* First column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_COLUMN); } if (elaboration_reg.east >= original_reg.east) { /* Last column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_COLUMN); last_column = TRUE; } nsplx = ceil((elaboration_reg.east - elaboration_reg.west) / stepE) + 0.5; /* if (nsplx > NSPLX_MAX) nsplx = NSPLX_MAX; */ G_debug(1, "nsplx = %d", nsplx); /*Setting the active region */ dim_vect = nsplx * nsply; observ = P_Read_Vector_Region_Map(&In, &elaboration_reg, &npoints, dim_vect, 1); if (npoints > 0) { /* If there is any point falling into elaboration_reg... */ int i; nparameters = nsplx * nsply; /* Mean calculation */ mean = P_Mean_Calc(&elaboration_reg, observ, npoints); /* Least Squares system */ G_debug(1, "Allocation memory for bilinear interpolation"); BW = P_get_BandWidth(P_BILINEAR, nsply); /* Bilinear interpolation */ N = G_alloc_matrix(nparameters, BW); /* Normal matrix */ TN = G_alloc_vector(nparameters); /* vector */ parVect = G_alloc_vector(nparameters); /* Bicubic parameters vector */ obsVect = G_alloc_matrix(npoints, 3); /* Observation vector */ Q = G_alloc_vector(npoints); /* "a priori" var-cov matrix */ lineVect = G_alloc_ivector(npoints); /* Setting obsVect vector & Q matrix */ for (i = 0; i < npoints; i++) { obsVect[i][0] = observ[i].coordX; obsVect[i][1] = observ[i].coordY; obsVect[i][2] = observ[i].coordZ - mean; lineVect[i] = observ[i].lineID; Q[i] = 1; /* Q=I */ } G_free(observ); G_verbose_message(_("Bilinear interpolation")); normalDefBilin(N, TN, Q, obsVect, stepE, stepN, nsplx, nsply, elaboration_reg.west, elaboration_reg.south, npoints, nparameters, BW); nCorrectGrad(N, lambda, nsplx, nsply, stepE, stepN); G_math_solver_cholesky_sband(N, parVect, TN, nparameters, BW); G_free_matrix(N); G_free_vector(TN); G_free_vector(Q); G_verbose_message(_("Outlier detection")); if (qgis_opt->answer) P_Outlier(&Out, &Outlier, &Qgis, elaboration_reg, general_box, overlap_box, obsVect, parVect, mean, dims.overlap, lineVect, npoints, driver, table_name); else P_Outlier(&Out, &Outlier, NULL, elaboration_reg, general_box, overlap_box, obsVect, parVect, mean, dims.overlap, lineVect, npoints, driver, table_name); G_free_vector(parVect); G_free_matrix(obsVect); G_free_ivector(lineVect); } /*! END IF; npoints > 0 */ else { G_free(observ); G_warning(_("No data within this subregion. " "Consider increasing spline step values.")); } } /*! END WHILE; last_column = TRUE */ } /*! END WHILE; last_row = TRUE */ /* Drop auxiliar table */ if (npoints > 0) { G_debug(1, "%s: Dropping <%s>", argv[0], table_name); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Auxiliary table could not be dropped")); } db_close_database_shutdown_driver(driver); Vect_close(&In); Vect_close(&Out); Vect_close(&Outlier); if (qgis_opt->answer) { Vect_build(&Qgis); Vect_close(&Qgis); } G_done_msg(" "); exit(EXIT_SUCCESS); } /*END MAIN */
int main(int argc, char **argv) { int i; int **cats, *ncats, nfields, *fields; struct Flag *line_flag; /* struct Flag *all_flag; */ struct Option *in_opt, *out_opt; struct Flag *table_flag; struct GModule *module; struct line_pnts *Points; struct line_cats *Cats; int node, nnodes; COOR *coor; int ncoor, acoor; int line, nlines, type, ctype, area, nareas; int err_boundaries, err_centr_out, err_centr_dupl, err_nocentr; G_gisinit(argv[0]); module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("geometry")); G_add_keyword(_("triangulation")); module->description = _("Creates a Voronoi diagram from an input vector " "map containing points or centroids."); in_opt = G_define_standard_option(G_OPT_V_INPUT); out_opt = G_define_standard_option(G_OPT_V_OUTPUT); /* all_flag = G_define_flag (); all_flag->key = 'a'; all_flag->description = _("Use all points (do not limit to current region)"); */ line_flag = G_define_flag(); line_flag->key = 'l'; line_flag->description = _("Output tessellation as a graph (lines), not areas"); table_flag = G_define_flag(); table_flag->key = 't'; table_flag->description = _("Do not create attribute table"); if (G_parser(argc, argv)) exit(EXIT_FAILURE); if (line_flag->answer) Type = GV_LINE; else Type = GV_BOUNDARY; All = 0; Points = Vect_new_line_struct(); Cats = Vect_new_cats_struct(); /* open files */ Vect_set_open_level(2); Vect_open_old(&In, in_opt->answer, ""); if (Vect_open_new(&Out, out_opt->answer, 0) < 0) G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); Vect_hist_copy(&In, &Out); Vect_hist_command(&Out); /* initialize working region */ G_get_window(&Window); Vect_region_box(&Window, &Box); Box.T = 0.5; Box.B = -0.5; freeinit(&sfl, sizeof(struct Site)); G_message(_("Reading sites...")); readsites(); siteidx = 0; geominit(); triangulate = 0; plot = 0; debug = 0; G_message(_("Voronoi triangulation...")); voronoi(triangulate, nextone); /* Close free ends by current region */ Vect_build_partial(&Out, GV_BUILD_BASE); ncoor = 0; acoor = 100; coor = (COOR *) G_malloc(sizeof(COOR) * acoor); nnodes = Vect_get_num_nodes(&Out); for (node = 1; node <= nnodes; node++) { double x, y; if (Vect_get_node_n_lines(&Out, node) < 2) { /* add coordinates */ Vect_get_node_coor(&Out, node, &x, &y, NULL); if (ncoor == acoor - 5) { /* always space for 5 region corners */ acoor += 100; coor = (COOR *) G_realloc(coor, sizeof(COOR) * acoor); } coor[ncoor].x = x; coor[ncoor].y = y; ncoor++; } } /* Add region corners */ coor[ncoor].x = Box.W; coor[ncoor].y = Box.S; ncoor++; coor[ncoor].x = Box.E; coor[ncoor].y = Box.S; ncoor++; coor[ncoor].x = Box.E; coor[ncoor].y = Box.N; ncoor++; coor[ncoor].x = Box.W; coor[ncoor].y = Box.N; ncoor++; /* Sort */ qsort(coor, ncoor, sizeof(COOR), (void *)cmp); /* add last (first corner) */ coor[ncoor].x = Box.W; coor[ncoor].y = Box.S; ncoor++; for (i = 1; i < ncoor; i++) { if (coor[i].x == coor[i - 1].x && coor[i].y == coor[i - 1].y) continue; /* duplicate */ Vect_reset_line(Points); Vect_append_point(Points, coor[i].x, coor[i].y, 0.0); Vect_append_point(Points, coor[i - 1].x, coor[i - 1].y, 0.0); Vect_write_line(&Out, Type, Points, Cats); } G_free(coor); /* Copy input points as centroids */ nfields = Vect_cidx_get_num_fields(&In); cats = (int **)G_malloc(nfields * sizeof(int *)); ncats = (int *)G_malloc(nfields * sizeof(int)); fields = (int *)G_malloc(nfields * sizeof(int)); for (i = 0; i < nfields; i++) { ncats[i] = 0; cats[i] = (int *)G_malloc(Vect_cidx_get_num_cats_by_index(&In, i) * sizeof(int)); fields[i] = Vect_cidx_get_field_number(&In, i); } if (line_flag->answer) ctype = GV_POINT; else ctype = GV_CENTROID; nlines = Vect_get_num_lines(&In); G_message(_("Writing sites to output...")); for (line = 1; line <= nlines; line++) { G_percent(line, nlines, 2); type = Vect_read_line(&In, Points, Cats, line); if (!(type & GV_POINTS)) continue; if (!Vect_point_in_box(Points->x[0], Points->y[0], 0.0, &Box)) continue; Vect_write_line(&Out, ctype, Points, Cats); for (i = 0; i < Cats->n_cats; i++) { int f, j; f = -1; for (j = 0; j < nfields; j++) { /* find field */ if (fields[j] == Cats->field[i]) { f = j; break; } } if (f > -1) { cats[f][ncats[f]] = Cats->cat[i]; ncats[f]++; } } } /* Copy tables */ if (!(table_flag->answer)) { int ttype, ntabs = 0; struct field_info *IFi, *OFi; /* Number of output tabs */ for (i = 0; i < Vect_get_num_dblinks(&In); i++) { int f, j; IFi = Vect_get_dblink(&In, i); f = -1; for (j = 0; j < nfields; j++) { /* find field */ if (fields[j] == IFi->number) { f = j; break; } } if (f > -1) { if (ncats[f] > 0) ntabs++; } } if (ntabs > 1) ttype = GV_MTABLE; else ttype = GV_1TABLE; for (i = 0; i < nfields; i++) { int ret; if (fields[i] == 0) continue; G_message(_("Layer %d"), fields[i]); /* Make a list of categories */ IFi = Vect_get_field(&In, fields[i]); if (!IFi) { /* no table */ G_message(_("No table")); continue; } OFi = Vect_default_field_info(&Out, IFi->number, IFi->name, ttype); ret = db_copy_table_by_ints(IFi->driver, IFi->database, IFi->table, OFi->driver, Vect_subst_var(OFi->database, &Out), OFi->table, IFi->key, cats[i], ncats[i]); if (ret == DB_FAILED) { G_warning(_("Cannot copy table")); } else { Vect_map_add_dblink(&Out, OFi->number, OFi->name, OFi->table, IFi->key, OFi->database, OFi->driver); } } } Vect_close(&In); /* cleaning part 1: count errors */ Vect_build_partial(&Out, GV_BUILD_CENTROIDS); err_boundaries = err_centr_out = err_centr_dupl = err_nocentr = 0; nlines = Vect_get_num_lines(&Out); for (line = 1; line <= nlines; line++) { if (!Vect_line_alive(&Out, line)) continue; type = Vect_get_line_type(&Out, line); if (type == GV_BOUNDARY) { int left, right; Vect_get_line_areas(&Out, line, &left, &right); if (left == 0 || right == 0) { G_debug(3, "line = %d left = %d right = %d", line, left, right); err_boundaries++; } } if (type == GV_CENTROID) { area = Vect_get_centroid_area(&Out, line); if (area == 0) err_centr_out++; else if (area < 0) err_centr_dupl++; } } err_nocentr = 0; nareas = Vect_get_num_areas(&Out); for (area = 1; area <= nareas; area++) { if (!Vect_area_alive(&Out, area)) continue; line = Vect_get_area_centroid(&Out, area); if (line == 0) err_nocentr++; } /* cleaning part 2: snap */ if (err_nocentr || err_centr_dupl || err_centr_out) { int nmod; G_important_message(_("Output needs topological cleaning")); Vect_snap_lines(&Out, GV_BOUNDARY, 1e-7, NULL); do { Vect_break_lines(&Out, GV_BOUNDARY, NULL); Vect_remove_duplicates(&Out, GV_BOUNDARY, NULL); nmod = Vect_clean_small_angles_at_nodes(&Out, GV_BOUNDARY, NULL); } while (nmod > 0); err_boundaries = 0; nlines = Vect_get_num_lines(&Out); for (line = 1; line <= nlines; line++) { if (!Vect_line_alive(&Out, line)) continue; type = Vect_get_line_type(&Out, line); if (type == GV_BOUNDARY) { int left, right; Vect_get_line_areas(&Out, line, &left, &right); if (left == 0 || right == 0) { G_debug(3, "line = %d left = %d right = %d", line, left, right); err_boundaries++; } } } } /* cleaning part 3: remove remaining incorrect boundaries */ if (err_boundaries) { G_important_message(_("Removing incorrect boundaries from output")); nlines = Vect_get_num_lines(&Out); for (line = 1; line <= nlines; line++) { if (!Vect_line_alive(&Out, line)) continue; type = Vect_get_line_type(&Out, line); if (type == GV_BOUNDARY) { int left, right; Vect_get_line_areas(&Out, line, &left, &right); /* &&, not ||, no typo */ if (left == 0 && right == 0) { G_debug(3, "line = %d left = %d right = %d", line, left, right); Vect_delete_line(&Out, line); } } } } /* build clean topology */ Vect_build_partial(&Out, GV_BUILD_NONE); Vect_build(&Out); Vect_close(&Out); G_done_msg(" "); exit(EXIT_SUCCESS); }
int main(int argc, char *argv[]) { /* Variables' declarations */ int nsplx_adj, nsply_adj; int nsubregion_col, nsubregion_row, subregion = 0, nsubregions = 0; double N_extension, E_extension, edgeE, edgeN; int dim_vect, nparameters, BW, npoints; double lambda_B, lambda_F, grad_H, grad_L, alpha, mean; const char *dvr, *db, *mapset; char table_interpolation[GNAME_MAX], table_name[GNAME_MAX]; char xname[GNAME_MAX], xmapset[GMAPSET_MAX]; int last_row, last_column, flag_auxiliar = FALSE; int *lineVect; double *TN, *Q, *parVect_bilin, *parVect_bicub; /* Interpolating and least-square vectors */ double **N, **obsVect; /* Interpolation and least-square matrix */ /* Structs' declarations */ struct Map_info In, Out; struct Option *in_opt, *out_opt, *stepE_opt, *stepN_opt, *lambdaF_opt, *lambdaB_opt, *gradH_opt, *gradL_opt, *alfa_opt; struct Flag *spline_step_flag; struct GModule *module; struct Cell_head elaboration_reg, original_reg; struct Reg_dimens dims; struct bound_box general_box, overlap_box; struct Point *observ; dbDriver *driver; /*------------------------------------------------------------------------------------------*/ /* Options' declaration */ module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("LIDAR")); G_add_keyword(_("edges")); module->description = _("Detects the object's edges from a LIDAR data set."); spline_step_flag = G_define_flag(); spline_step_flag->key = 'e'; spline_step_flag->label = _("Estimate point density and distance"); spline_step_flag->description = _("Estimate point density and distance for the input vector points within the current region extends and quit"); in_opt = G_define_standard_option(G_OPT_V_INPUT); out_opt = G_define_standard_option(G_OPT_V_OUTPUT); stepE_opt = G_define_option(); stepE_opt->key = "see"; stepE_opt->type = TYPE_DOUBLE; stepE_opt->required = NO; stepE_opt->answer = "4"; stepE_opt->description = _("Interpolation spline step value in east direction"); stepE_opt->guisection = _("Settings"); stepN_opt = G_define_option(); stepN_opt->key = "sen"; stepN_opt->type = TYPE_DOUBLE; stepN_opt->required = NO; stepN_opt->answer = "4"; stepN_opt->description = _("Interpolation spline step value in north direction"); stepN_opt->guisection = _("Settings"); lambdaB_opt = G_define_option(); lambdaB_opt->key = "lambda_g"; lambdaB_opt->type = TYPE_DOUBLE; lambdaB_opt->required = NO; lambdaB_opt->description = _("Regularization weight in gradient evaluation"); lambdaB_opt->answer = "0.01"; lambdaB_opt->guisection = _("Settings"); gradH_opt = G_define_option(); gradH_opt->key = "tgh"; gradH_opt->type = TYPE_DOUBLE; gradH_opt->required = NO; gradH_opt->description = _("High gradient threshold for edge classification"); gradH_opt->answer = "6"; gradH_opt->guisection = _("Settings"); gradL_opt = G_define_option(); gradL_opt->key = "tgl"; gradL_opt->type = TYPE_DOUBLE; gradL_opt->required = NO; gradL_opt->description = _("Low gradient threshold for edge classification"); gradL_opt->answer = "3"; gradL_opt->guisection = _("Settings"); alfa_opt = G_define_option(); alfa_opt->key = "theta_g"; alfa_opt->type = TYPE_DOUBLE; alfa_opt->required = NO; alfa_opt->description = _("Angle range for same direction detection"); alfa_opt->answer = "0.26"; alfa_opt->guisection = _("Settings"); lambdaF_opt = G_define_option(); lambdaF_opt->key = "lambda_r"; lambdaF_opt->type = TYPE_DOUBLE; lambdaF_opt->required = NO; lambdaF_opt->description = _("Regularization weight in residual evaluation"); lambdaF_opt->answer = "2"; lambdaF_opt->guisection = _("Settings"); /* Parsing */ G_gisinit(argv[0]); if (G_parser(argc, argv)) exit(EXIT_FAILURE); line_out_counter = 1; stepN = atof(stepN_opt->answer); stepE = atof(stepE_opt->answer); lambda_F = atof(lambdaF_opt->answer); lambda_B = atof(lambdaB_opt->answer); grad_H = atof(gradH_opt->answer); grad_L = atof(gradL_opt->answer); alpha = atof(alfa_opt->answer); grad_L = grad_L * grad_L; grad_H = grad_H * grad_H; if (!(db = G__getenv2("DB_DATABASE", G_VAR_MAPSET))) G_fatal_error(_("Unable to read name of database")); if (!(dvr = G__getenv2("DB_DRIVER", G_VAR_MAPSET))) G_fatal_error(_("Unable to read name of driver")); /* Setting auxiliar table's name */ if (G_name_is_fully_qualified(out_opt->answer, xname, xmapset)) { sprintf(table_name, "%s_aux", xname); sprintf(table_interpolation, "%s_edge_Interpolation", xname); } else { sprintf(table_name, "%s_aux", out_opt->answer); sprintf(table_interpolation, "%s_edge_Interpolation", out_opt->answer); } /* Something went wrong in a previous v.lidar.edgedetection execution */ if (db_table_exists(dvr, db, table_name)) { /* Start driver and open db */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_fatal_error(_("Old auxiliar table could not be dropped")); db_close_database_shutdown_driver(driver); } /* Something went wrong in a previous v.lidar.edgedetection execution */ if (db_table_exists(dvr, db, table_interpolation)) { /* Start driver and open db */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); if (P_Drop_Aux_Table(driver, table_interpolation) != DB_OK) G_fatal_error(_("Old auxiliar table could not be dropped")); db_close_database_shutdown_driver(driver); } /* Checking vector names */ Vect_check_input_output_name(in_opt->answer, out_opt->answer, G_FATAL_EXIT); if ((mapset = G_find_vector2(in_opt->answer, "")) == NULL) { G_fatal_error(_("Vector map <%s> not found"), in_opt->answer); } Vect_set_open_level(1); /* Open input vector */ if (1 > Vect_open_old(&In, in_opt->answer, mapset)) G_fatal_error(_("Unable to open vector map <%s>"), in_opt->answer); /* Input vector must be 3D */ if (!Vect_is_3d(&In)) G_fatal_error(_("Input vector map <%s> is not 3D!"), in_opt->answer); /* Estimate point density and mean distance for current region */ if (spline_step_flag->answer) { double dens, dist; if (P_estimate_splinestep(&In, &dens, &dist) == 0) { G_message("Estimated point density: %.4g", dens); G_message("Estimated mean distance between points: %.4g", dist); } else G_warning(_("No points in current region!")); Vect_close(&In); exit(EXIT_SUCCESS); } /* Open output vector */ if (0 > Vect_open_new(&Out, out_opt->answer, WITH_Z)) G_fatal_error(_("Unable to create vector map <%s>"), out_opt->answer); /* Copy vector Head File */ Vect_copy_head_data(&In, &Out); Vect_hist_copy(&In, &Out); Vect_hist_command(&Out); /* Start driver and open db */ driver = db_start_driver_open_database(dvr, db); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), dvr); db_set_error_handler_driver(driver); /* Create auxiliar and interpolation table */ if ((flag_auxiliar = P_Create_Aux4_Table(driver, table_name)) == FALSE) G_fatal_error(_("It was impossible to create <%s>."), table_name); if (P_Create_Aux2_Table(driver, table_interpolation) == FALSE) G_fatal_error(_("It was impossible to create <%s> interpolation table in database."), out_opt->answer); db_create_index2(driver, table_name, "ID"); db_create_index2(driver, table_interpolation, "ID"); /* sqlite likes that ??? */ db_close_database_shutdown_driver(driver); driver = db_start_driver_open_database(dvr, db); /* Setting regions and boxes */ G_get_set_window(&original_reg); G_get_set_window(&elaboration_reg); Vect_region_box(&elaboration_reg, &overlap_box); Vect_region_box(&elaboration_reg, &general_box); /*------------------------------------------------------------------ | Subdividing and working with tiles: | Each original region will be divided into several subregions. | Each one will be overlaped by its neighbouring subregions. | The overlapping is calculated as a fixed OVERLAP_SIZE times | the largest spline step plus 2 * edge ----------------------------------------------------------------*/ /* Fixing parameters of the elaboration region */ P_zero_dim(&dims); nsplx_adj = NSPLX_MAX; nsply_adj = NSPLY_MAX; if (stepN > stepE) dims.overlap = OVERLAP_SIZE * stepN; else dims.overlap = OVERLAP_SIZE * stepE; P_get_edge(P_BICUBIC, &dims, stepE, stepN); P_set_dim(&dims, stepE, stepN, &nsplx_adj, &nsply_adj); G_verbose_message(_("adjusted EW splines %d"), nsplx_adj); G_verbose_message(_("adjusted NS splines %d"), nsply_adj); /* calculate number of subregions */ edgeE = dims.ew_size - dims.overlap - 2 * dims.edge_v; edgeN = dims.sn_size - dims.overlap - 2 * dims.edge_h; N_extension = original_reg.north - original_reg.south; E_extension = original_reg.east - original_reg.west; nsubregion_col = ceil(E_extension / edgeE) + 0.5; nsubregion_row = ceil(N_extension / edgeN) + 0.5; if (nsubregion_col < 0) nsubregion_col = 0; if (nsubregion_row < 0) nsubregion_row = 0; nsubregions = nsubregion_row * nsubregion_col; elaboration_reg.south = original_reg.north; last_row = FALSE; while (last_row == FALSE) { /* For each row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_ROW); if (elaboration_reg.north > original_reg.north) { /* First row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_ROW); } if (elaboration_reg.south <= original_reg.south) { /* Last row */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_ROW); last_row = TRUE; } nsply = ceil((elaboration_reg.north - elaboration_reg.south) / stepN) + 0.5; /* if (nsply > NSPLY_MAX) { nsply = NSPLY_MAX; } */ G_debug(1, "nsply = %d", nsply); elaboration_reg.east = original_reg.west; last_column = FALSE; while (last_column == FALSE) { /* For each column */ subregion++; if (nsubregions > 1) G_message(_("subregion %d of %d"), subregion, nsubregions); P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, GENERAL_COLUMN); if (elaboration_reg.west < original_reg.west) { /* First column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, FIRST_COLUMN); } if (elaboration_reg.east >= original_reg.east) { /* Last column */ P_set_regions(&elaboration_reg, &general_box, &overlap_box, dims, LAST_COLUMN); last_column = TRUE; } nsplx = ceil((elaboration_reg.east - elaboration_reg.west) / stepE) + 0.5; /* if (nsplx > NSPLX_MAX) { nsplx = NSPLX_MAX; } */ G_debug(1, "nsplx = %d", nsplx); /*Setting the active region */ dim_vect = nsplx * nsply; G_debug(1, "read vector region map"); observ = P_Read_Vector_Region_Map(&In, &elaboration_reg, &npoints, dim_vect, 1); if (npoints > 0) { /* If there is any point falling into elaboration_reg... */ int i, tn; nparameters = nsplx * nsply; /* Mean's calculation */ mean = P_Mean_Calc(&elaboration_reg, observ, npoints); /* Least Squares system */ G_debug(1, _("Allocating memory for bilinear interpolation")); BW = P_get_BandWidth(P_BILINEAR, nsply); /* Bilinear interpolation */ N = G_alloc_matrix(nparameters, BW); /* Normal matrix */ TN = G_alloc_vector(nparameters); /* vector */ parVect_bilin = G_alloc_vector(nparameters); /* Bilinear parameters vector */ obsVect = G_alloc_matrix(npoints + 1, 3); /* Observation vector */ Q = G_alloc_vector(npoints + 1); /* "a priori" var-cov matrix */ lineVect = G_alloc_ivector(npoints + 1); /* Setting obsVect vector & Q matrix */ for (i = 0; i < npoints; i++) { obsVect[i][0] = observ[i].coordX; obsVect[i][1] = observ[i].coordY; obsVect[i][2] = observ[i].coordZ - mean; lineVect[i] = observ[i].lineID; Q[i] = 1; /* Q=I */ } G_free(observ); G_verbose_message(_("Bilinear interpolation")); normalDefBilin(N, TN, Q, obsVect, stepE, stepN, nsplx, nsply, elaboration_reg.west, elaboration_reg.south, npoints, nparameters, BW); nCorrectGrad(N, lambda_B, nsplx, nsply, stepE, stepN); G_math_solver_cholesky_sband(N, parVect_bilin, TN, nparameters, BW); G_free_matrix(N); for (tn = 0; tn < nparameters; tn++) TN[tn] = 0; G_debug(1, _("Allocating memory for bicubic interpolation")); BW = P_get_BandWidth(P_BICUBIC, nsply); N = G_alloc_matrix(nparameters, BW); /* Normal matrix */ parVect_bicub = G_alloc_vector(nparameters); /* Bicubic parameters vector */ G_verbose_message(_("Bicubic interpolation")); normalDefBicubic(N, TN, Q, obsVect, stepE, stepN, nsplx, nsply, elaboration_reg.west, elaboration_reg.south, npoints, nparameters, BW); nCorrectLapl(N, lambda_F, nsplx, nsply, stepE, stepN); G_math_solver_cholesky_sband(N, parVect_bicub, TN, nparameters, BW); G_free_matrix(N); G_free_vector(TN); G_free_vector(Q); G_verbose_message(_("Point classification")); classification(&Out, elaboration_reg, general_box, overlap_box, obsVect, parVect_bilin, parVect_bicub, mean, alpha, grad_H, grad_L, dims.overlap, lineVect, npoints, driver, table_interpolation, table_name); G_free_vector(parVect_bilin); G_free_vector(parVect_bicub); G_free_matrix(obsVect); G_free_ivector(lineVect); } /* IF */ else { G_free(observ); G_warning(_("No data within this subregion. " "Consider changing the spline step.")); } } /*! END WHILE; last_column = TRUE */ } /*! END WHILE; last_row = TRUE */ /* Dropping auxiliar table */ if (npoints > 0) { G_debug(1, _("Dropping <%s>"), table_name); if (P_Drop_Aux_Table(driver, table_name) != DB_OK) G_warning(_("Auxiliar table could not be dropped")); } db_close_database_shutdown_driver(driver); Vect_close(&In); Vect_map_add_dblink(&Out, F_INTERPOLATION, NULL, table_interpolation, "id", db, dvr); Vect_close(&Out); G_done_msg(" "); exit(EXIT_SUCCESS); } /*!END MAIN */
/*--------------------------------------------------------------------------------*/ int main(int argc, char *argv[]) { /* Variables' declarations */ int row, nrows, col, ncols, MaxPoints; int nsubregion_col, nsubregion_row; int subregion = 0, nsubregions = 0; int last_row, last_column; int nlines, nlines_first, line_num; int more; int clas, region = TRUE; double Z_interp; double Thres_j, Thres_d, ew_resol, ns_resol; double minNS, minEW, maxNS, maxEW; const char *mapset; char buf[1024], table_name[GNAME_MAX]; char xname[GNAME_MAX], xmapset[GMAPSET_MAX]; int colorBordo, ripieno, conta, lungPunti, lungHull, xi, c1, c2; double altPiano; extern double **P, **cvxHull, **punti_bordo; /* Struct declarations */ struct Cell_head elaboration_reg, original_reg; struct element_grow **raster_matrix; struct Map_info In, Out, First; struct Option *in_opt, *out_opt, *first_opt, *Thres_j_opt, *Thres_d_opt; struct GModule *module; struct line_pnts *points, *points_first; struct line_cats *Cats, *Cats_first; struct field_info *field; dbDriver *driver; dbString sql; dbTable *table; dbCursor cursor; /*------------------------------------------------------------------------------------------*/ /* Options' declaration */ ; module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("LIDAR")); module->description = _("Building contour determination and Region Growing " "algorithm for determining the building inside"); in_opt = G_define_standard_option(G_OPT_V_INPUT); in_opt->description = _("Input vector (v.lidar.edgedetection output"); out_opt = G_define_standard_option(G_OPT_V_OUTPUT); first_opt = G_define_option(); first_opt->key = "first"; first_opt->type = TYPE_STRING; first_opt->key_desc = "name"; first_opt->required = YES; first_opt->gisprompt = "old,vector,vector"; first_opt->description = _("Name of the first pulse vector map"); Thres_j_opt = G_define_option(); Thres_j_opt->key = "tj"; Thres_j_opt->type = TYPE_DOUBLE; Thres_j_opt->required = NO; Thres_j_opt->description = _("Threshold for cell object frequency in region growing"); Thres_j_opt->answer = "0.2"; Thres_d_opt = G_define_option(); Thres_d_opt->key = "td"; Thres_d_opt->type = TYPE_DOUBLE; Thres_d_opt->required = NO; Thres_d_opt->description = _("Threshold for double pulse in region growing"); Thres_d_opt->answer = "0.6"; /* Parsing */ G_gisinit(argv[0]); if (G_parser(argc, argv)) exit(EXIT_FAILURE); Thres_j = atof(Thres_j_opt->answer); Thres_d = atof(Thres_d_opt->answer); Thres_j += 1; /* Open input vector */ Vect_check_input_output_name(in_opt->answer, out_opt->answer, GV_FATAL_EXIT); if ((mapset = G_find_vector2(in_opt->answer, "")) == NULL) { G_fatal_error(_("Vector map <%s> not found"), in_opt->answer); } /* Setting auxiliar table's name */ if (G_name_is_fully_qualified(in_opt->answer, xname, xmapset)) { sprintf(table_name, "%s_edge_Interpolation", xname); } else sprintf(table_name, "%s_edge_Interpolation", in_opt->answer); Vect_set_open_level(1); /* WITHOUT TOPOLOGY */ if (Vect_open_old(&In, in_opt->answer, mapset) < 1) G_fatal_error(_("Unable to open vector map <%s>"), in_opt->answer); Vect_set_open_level(1); /* WITHOUT TOPOLOGY */ if (Vect_open_old(&First, first_opt->answer, mapset) < 1) G_fatal_error(_("Unable to open vector map <%s>"), first_opt->answer); /* Open output vector */ if (0 > Vect_open_new(&Out, out_opt->answer, WITH_Z)) { Vect_close(&In); Vect_close(&First); exit(EXIT_FAILURE); } /* Copy vector Head File */ Vect_copy_head_data(&In, &Out); Vect_hist_copy(&In, &Out); Vect_hist_command(&Out); /* Starting driver and open db for edgedetection interpolation table */ field = Vect_get_field(&In, F_INTERPOLATION); /*if (field == NULL) G_fatal_error (_("Cannot read field info")); */ driver = db_start_driver_open_database(field->driver, field->database); if (driver == NULL) G_fatal_error(_("No database connection for driver <%s> is defined. Run db.connect."), field->driver); /* is this the right place to open the cursor ??? */ db_init_string(&sql); db_zero_string(&sql); sprintf(buf, "SELECT Interp,ID FROM %s", table_name); G_debug(1, "buf: %s", buf); db_append_string(&sql, buf); if (db_open_select_cursor(driver, &sql, &cursor, DB_SEQUENTIAL) != DB_OK) G_fatal_error(_("Unable to open table <%s>"), table_name); count_obj = 1; /* no topology, get number of lines in input vector */ nlines = 0; points = Vect_new_line_struct(); Cats = Vect_new_cats_struct(); Vect_rewind(&In); while (Vect_read_next_line(&In, points, Cats) > 0) { nlines++; } Vect_rewind(&In); /* no topology, get number of lines in first pulse input vector */ nlines_first = 0; points_first = Vect_new_line_struct(); Cats_first = Vect_new_cats_struct(); Vect_rewind(&First); while (Vect_read_next_line(&First, points_first, Cats_first) > 0) { nlines_first++; } Vect_rewind(&First); /* Setting regions and boxes */ G_debug(1, _("Setting regions and boxes")); G_get_set_window(&original_reg); G_get_set_window(&elaboration_reg); /* Fixing parameters of the elaboration region */ /*! The original_region will be divided into subregions */ ew_resol = original_reg.ew_res; ns_resol = original_reg.ns_res; /* calculate number of subregions */ nsubregion_col = ceil((original_reg.east - original_reg.west) / (LATO * ew_resol)) + 0.5; nsubregion_row = ceil((original_reg.north - original_reg.south) / (LATO * ns_resol)) + 0.5; if (nsubregion_col < 0) nsubregion_col = 0; if (nsubregion_row < 0) nsubregion_row = 0; nsubregions = nsubregion_row * nsubregion_col; /* Subdividing and working with tiles */ elaboration_reg.south = original_reg.north; last_row = FALSE; while (last_row == FALSE) { /* For each strip of LATO rows */ elaboration_reg.north = elaboration_reg.south; if (elaboration_reg.north > original_reg.north) /* First row */ elaboration_reg.north = original_reg.north; elaboration_reg.south = elaboration_reg.north - LATO * ns_resol; if (elaboration_reg.south <= original_reg.south) { /* Last row */ elaboration_reg.south = original_reg.south; last_row = TRUE; } elaboration_reg.east = original_reg.west; last_column = FALSE; while (last_column == FALSE) { /* For each strip of LATO columns */ struct bound_box elaboration_box; subregion++; if (nsubregions > 1) G_message(_("subregion %d of %d"), subregion, nsubregions); elaboration_reg.west = elaboration_reg.east; if (elaboration_reg.west < original_reg.west) /* First column */ elaboration_reg.west = original_reg.west; elaboration_reg.east = elaboration_reg.west + LATO * ew_resol; if (elaboration_reg.east >= original_reg.east) { /* Last column */ elaboration_reg.east = original_reg.east; last_column = TRUE; } /* Setting the active region */ elaboration_reg.ns_res = ns_resol; elaboration_reg.ew_res = ew_resol; nrows = (elaboration_reg.north - elaboration_reg.south) / ns_resol + 0.1; ncols = (elaboration_reg.east - elaboration_reg.west) / ew_resol + 0.1; elaboration_reg.rows = nrows; elaboration_reg.cols = ncols; G_debug(1, _("Rows = %d"), nrows); G_debug(1, _("Columns = %d"), ncols); raster_matrix = structMatrix(0, nrows, 0, ncols); MaxPoints = nrows * ncols; /* Initializing matrix */ for (row = 0; row <= nrows; row++) { for (col = 0; col <= ncols; col++) { raster_matrix[row][col].interp = 0; raster_matrix[row][col].fi = 0; raster_matrix[row][col].bordo = 0; raster_matrix[row][col].dueImp = SINGLE_PULSE; raster_matrix[row][col].orig = 0; raster_matrix[row][col].fo = 0; raster_matrix[row][col].clas = PRE_TERRAIN; raster_matrix[row][col].fc = 0; raster_matrix[row][col].obj = 0; } } G_verbose_message(_("read points in input vector")); Vect_region_box(&elaboration_reg, &elaboration_box); line_num = 0; Vect_rewind(&In); while (Vect_read_next_line(&In, points, Cats) > 0) { line_num++; if ((Vect_point_in_box (points->x[0], points->y[0], points->z[0], &elaboration_box)) && ((points->x[0] != elaboration_reg.west) || (points->x[0] == original_reg.west)) && ((points->y[0] != elaboration_reg.north) || (points->y[0] == original_reg.north))) { row = (int)(Rast_northing_to_row (points->y[0], &elaboration_reg)); col = (int)(Rast_easting_to_col (points->x[0], &elaboration_reg)); Z_interp = 0; /* TODO: make sure the current db_fetch() usage works */ /* why not: */ /* db_init_string(&sql); sprintf(buf, "SELECT Interp,ID FROM %s WHERE ID=%d", table_name, line_num); db_append_string(&sql, buf); if (db_open_select_cursor(driver, &sql, &cursor, DB_SEQUENTIAL) != DB_OK) G_fatal_error(_("Unable to open table <%s>"), table_name); while (db_fetch(&cursor, DB_NEXT, &more) == DB_OK && more) { dbColumn *Z_Interp_col; dbValue *Z_Interp_value; table = db_get_cursor_table(&cursor); Z_Interp_col = db_get_table_column(table, 1); if (db_sqltype_to_Ctype(db_get_column_sqltype(Z_Interp_col)) == DB_C_TYPE_DOUBLE) Z_Interp_value = db_get_column_value(Z_Interp_col); else continue; Z_interp = db_get_value_double(Z_Interp_value); break; } db_close_cursor(&cursor); db_free_string(&sql); */ /* instead of */ while (1) { if (db_fetch(&cursor, DB_NEXT, &more) != DB_OK || !more) break; dbColumn *Z_Interp_col, *ID_col; dbValue *Z_Interp_value, *ID_value; table = db_get_cursor_table(&cursor); ID_col = db_get_table_column(table, 1); if (db_sqltype_to_Ctype(db_get_column_sqltype(ID_col)) == DB_C_TYPE_INT) ID_value = db_get_column_value(ID_col); else continue; if (db_get_value_int(ID_value) == line_num) { Z_Interp_col = db_get_table_column(table, 0); if (db_sqltype_to_Ctype (db_get_column_sqltype(Z_Interp_col)) == DB_C_TYPE_DOUBLE) Z_Interp_value = db_get_column_value(Z_Interp_col); else continue; Z_interp = db_get_value_double(Z_Interp_value); break; } } raster_matrix[row][col].interp += Z_interp; raster_matrix[row][col].fi++; /*if (( clas = Vect_get_line_cat (&In, line_num, F_EDGE_DETECTION_CLASS) ) != UNKNOWN_EDGE) { */ if (Vect_cat_get(Cats, F_EDGE_DETECTION_CLASS, &clas)) { raster_matrix[row][col].clas += clas; raster_matrix[row][col].fc++; } raster_matrix[row][col].orig += points->z[0]; raster_matrix[row][col].fo++; } Vect_reset_cats(Cats); Vect_reset_line(points); } for (row = 0; row <= nrows; row++) { for (col = 0; col <= ncols; col++) { if (raster_matrix[row][col].fc != 0) { raster_matrix[row][col].clas--; raster_matrix[row][col]. clas /= raster_matrix[row][col].fc; } if (raster_matrix[row][col].fi != 0) raster_matrix[row][col]. interp /= raster_matrix[row][col].fi; if (raster_matrix[row][col].fo != 0) raster_matrix[row][col]. orig /= raster_matrix[row][col].fo; } } /* DOUBLE IMPULSE */ Vect_rewind(&First); while (Vect_read_next_line(&First, points_first, Cats_first) > 0) { if ((Vect_point_in_box (points_first->x[0], points_first->y[0], points_first->z[0], &elaboration_box)) && ((points->x[0] != elaboration_reg.west) || (points->x[0] == original_reg.west)) && ((points->y[0] != elaboration_reg.north) || (points->y[0] == original_reg.north))) { row = (int)(Rast_northing_to_row (points_first->y[0], &elaboration_reg)); col = (int)(Rast_easting_to_col (points_first->x[0], &elaboration_reg)); if (fabs (points_first->z[0] - raster_matrix[row][col].orig) >= Thres_d) raster_matrix[row][col].dueImp = DOUBLE_PULSE; } Vect_reset_cats(Cats_first); Vect_reset_line(points_first); } /* REGION GROWING */ if (region == TRUE) { G_verbose_message(_("Region Growing")); punti_bordo = G_alloc_matrix(MaxPoints, 3); P = Pvector(0, MaxPoints); colorBordo = 5; ripieno = 6; for (row = 0; row <= nrows; row++) { G_percent(row, nrows, 2); for (col = 0; col <= ncols; col++) { if ((raster_matrix[row][col].clas >= Thres_j) && (raster_matrix[row][col].clas < colorBordo) && (raster_matrix[row][col].fi != 0) && (raster_matrix[row][col].dueImp == SINGLE_PULSE)) { /* Selecting a connected Object zone */ ripieno++; if (ripieno > 10) ripieno = 6; /* Selecting points on a connected edge */ for (conta = 0; conta < MaxPoints; conta++) { punti_bordo[conta][0] = 0; punti_bordo[conta][1] = 0; punti_bordo[conta][2] = 0; P[conta] = punti_bordo[conta]; /* It only makes indexes to be equal, not coord values!! */ } lungPunti = 0; lungHull = 0; regGrow8(elaboration_reg, raster_matrix, punti_bordo, &lungPunti, row, col, colorBordo, Thres_j, MaxPoints); /* CONVEX-HULL COMPUTATION */ lungHull = ch2d(P, lungPunti); cvxHull = G_alloc_matrix(lungHull, 3); for (xi = 0; xi < lungHull; xi++) { cvxHull[xi][0] = P[xi][0]; cvxHull[xi][1] = P[xi][1]; cvxHull[xi][2] = P[xi][2]; } /* Computes the interpoling plane based only on Object points */ altPiano = pianOriz(punti_bordo, lungPunti, &minNS, &minEW, &maxNS, &maxEW, raster_matrix, colorBordo); for (c1 = minNS; c1 <= maxNS; c1++) { for (c2 = minEW; c2 <= maxEW; c2++) { if (checkHull(c1, c2, cvxHull, lungHull) == 1) { raster_matrix[c1][c2].obj = count_obj; if ((raster_matrix[c1][c2].clas == PRE_TERRAIN) && (raster_matrix[c1][c2].orig >= altPiano) && (lungHull > 3)) raster_matrix[c1][c2].clas = ripieno; } } } G_free_matrix(cvxHull); count_obj++; } } } G_free_matrix(punti_bordo); free_Pvector(P, 0, MaxPoints); } /* WRITING THE OUTPUT VECTOR CATEGORIES */ Vect_rewind(&In); while (Vect_read_next_line(&In, points, Cats) > 0) { /* Read every line for buffering points */ if ((Vect_point_in_box (points->x[0], points->y[0], points->z[0], &elaboration_box)) && ((points->x[0] != elaboration_reg.west) || (points->x[0] == original_reg.west)) && ((points->y[0] != elaboration_reg.north) || (points->y[0] == original_reg.north))) { row = (int)(Rast_northing_to_row (points->y[0], &elaboration_reg)); col = (int)(Rast_easting_to_col (points->x[0], &elaboration_reg)); if (raster_matrix[row][col].clas == PRE_TERRAIN) { if (raster_matrix[row][col].dueImp == SINGLE_PULSE) Vect_cat_set(Cats, F_CLASSIFICATION, TERRAIN_SINGLE); else Vect_cat_set(Cats, F_CLASSIFICATION, TERRAIN_DOUBLE); } else { if (raster_matrix[row][col].dueImp == SINGLE_PULSE) Vect_cat_set(Cats, F_CLASSIFICATION, OBJECT_SINGLE); else Vect_cat_set(Cats, F_CLASSIFICATION, OBJECT_DOUBLE); } Vect_cat_set(Cats, F_COUNTER_OBJ, raster_matrix[row][col].obj); Vect_write_line(&Out, GV_POINT, points, Cats); } Vect_reset_cats(Cats); Vect_reset_line(points); } free_structmatrix(raster_matrix, 0, nrows - 1, 0, ncols - 1); } /*! END WHILE; last_column = TRUE */ } /*! END WHILE; last_row = TRUE */ Vect_close(&In); Vect_close(&First); Vect_close(&Out); db_close_database_shutdown_driver(driver); G_done_msg(" "); exit(EXIT_SUCCESS); }
int main(int argc, char **argv) { int i, nsites, warn_once = 0; int all; long x, y; struct Cell_head window; struct GModule *module; struct { struct Option *input, *tests, *dfield, *layer; } parm; struct { struct Flag *q, *l, *region; } flag; double *w, *z; struct Map_info Map; int line, nlines, npoints; int field; struct line_pnts *Points; struct line_cats *Cats; struct bound_box box; /* Attributes */ int nrecords; int ctype; struct field_info *Fi; dbDriver *Driver; dbCatValArray cvarr; G_gisinit(argv[0]); module = G_define_module(); G_add_keyword(_("vector")); G_add_keyword(_("statistics")); G_add_keyword(_("points")); G_add_keyword(_("point pattern")); module->description = _("Tests for normality for vector points."); parm.input = G_define_standard_option(G_OPT_V_MAP); parm.layer = G_define_standard_option(G_OPT_V_FIELD); parm.tests = G_define_option(); parm.tests->key = "tests"; parm.tests->key_desc = "range"; parm.tests->type = TYPE_STRING; parm.tests->multiple = YES; parm.tests->required = YES; parm.tests->label = _("Lists of tests (1-15)"); parm.tests->description = _("E.g. 1,3-8,13"); parm.dfield = G_define_standard_option(G_OPT_DB_COLUMN); parm.dfield->required = YES; flag.region = G_define_flag(); flag.region->key = 'r'; flag.region->description = _("Use only points in current region"); flag.l = G_define_flag(); flag.l->key = 'l'; flag.l->description = _("Lognormality instead of normality"); if (G_parser(argc, argv)) exit(EXIT_FAILURE); all = flag.region->answer ? 0 : 1; /* Open input */ Vect_set_open_level(2); if (Vect_open_old2(&Map, parm.input->answer, "", parm.layer->answer) < 0) G_fatal_error(_("Unable to open vector map <%s>"), parm.input->answer); field = Vect_get_field_number(&Map, parm.layer->answer); /* Read attributes */ Fi = Vect_get_field(&Map, field); if (Fi == NULL) { G_fatal_error("Database connection not defined for layer %d", field); } Driver = db_start_driver_open_database(Fi->driver, Fi->database); if (Driver == NULL) G_fatal_error(_("Unable to open database <%s> by driver <%s>"), Fi->database, Fi->driver); nrecords = db_select_CatValArray(Driver, Fi->table, Fi->key, parm.dfield->answer, NULL, &cvarr); G_debug(1, "nrecords = %d", nrecords); ctype = cvarr.ctype; if (ctype != DB_C_TYPE_INT && ctype != DB_C_TYPE_DOUBLE) G_fatal_error(_("Only numeric column type supported")); if (nrecords < 0) G_fatal_error(_("Unable to select data from table")); G_verbose_message(_("%d records selected from table"), nrecords); db_close_database_shutdown_driver(Driver); /* Read points */ npoints = Vect_get_num_primitives(&Map, GV_POINT); z = (double *)G_malloc(npoints * sizeof(double)); G_get_window(&window); Vect_region_box(&window, &box); Points = Vect_new_line_struct(); Cats = Vect_new_cats_struct(); nlines = Vect_get_num_lines(&Map); nsites = 0; for (line = 1; line <= nlines; line++) { int type, cat, ret, cval; double dval; G_debug(3, "line = %d", line); type = Vect_read_line(&Map, Points, Cats, line); if (!(type & GV_POINT)) continue; if (!all) { if (!Vect_point_in_box(Points->x[0], Points->y[0], 0.0, &box)) continue; } Vect_cat_get(Cats, 1, &cat); G_debug(3, "cat = %d", cat); /* find actual value */ if (ctype == DB_C_TYPE_INT) { ret = db_CatValArray_get_value_int(&cvarr, cat, &cval); if (ret != DB_OK) { G_warning(_("No record for cat %d"), cat); continue; } dval = cval; } else if (ctype == DB_C_TYPE_DOUBLE) { ret = db_CatValArray_get_value_double(&cvarr, cat, &dval); if (ret != DB_OK) { G_warning(_("No record for cat %d"), cat); continue; } } G_debug(3, "dval = %e", dval); z[nsites] = dval; nsites++; } G_verbose_message(_("Number of points: %d"), nsites); if (nsites <= 0) G_fatal_error(_("No points found")); if (nsites < 4) G_warning(_("Too small sample")); if (flag.l->answer) { warn_once = 0; for (i = 0; i < nsites; ++i) { if (z[i] > 1.0e-10) z[i] = log10(z[i]); else if (!warn_once) { G_warning(_("Negative or very small point values set to -10.0")); z[i] = -10.0; warn_once = 1; } } } for (i = 0; parm.tests->answers[i]; i++) if (!scan_cats(parm.tests->answers[i], &x, &y)) { G_usage(); exit(EXIT_FAILURE); } for (i = 0; parm.tests->answers[i]; i++) { scan_cats(parm.tests->answers[i], &x, &y); while (x <= y) switch (x++) { case 1: /* moments */ fprintf(stdout, _("Moments \\sqrt{b_1} and b_2: ")); w = Cdhc_omnibus_moments(z, nsites); fprintf(stdout, "%g %g\n", w[0], w[1]); break; case 2: /* geary */ fprintf(stdout, _("Geary's a-statistic & an approx. normal: ")); w = Cdhc_geary_test(z, nsites); fprintf(stdout, "%g %g\n", w[0], w[1]); break; case 3: /* extreme deviates */ fprintf(stdout, _("Extreme normal deviates: ")); w = Cdhc_extreme(z, nsites); fprintf(stdout, "%g %g\n", w[0], w[1]); break; case 4: /* D'Agostino */ fprintf(stdout, _("D'Agostino's D & an approx. normal: ")); w = Cdhc_dagostino_d(z, nsites); fprintf(stdout, "%g %g\n", w[0], w[1]); break; case 5: /* Kuiper */ fprintf(stdout, _("Kuiper's V (regular & modified for normality): ")); w = Cdhc_kuipers_v(z, nsites); fprintf(stdout, "%g %g\n", w[1], w[0]); break; case 6: /* Watson */ fprintf(stdout, _("Watson's U^2 (regular & modified for normality): ")); w = Cdhc_watson_u2(z, nsites); fprintf(stdout, "%g %g\n", w[1], w[0]); break; case 7: /* Durbin */ fprintf(stdout, _("Durbin's Exact Test (modified Kolmogorov): ")); w = Cdhc_durbins_exact(z, nsites); fprintf(stdout, "%g\n", w[0]); break; case 8: /* Anderson-Darling */ fprintf(stdout, _("Anderson-Darling's A^2 (regular & modified for normality): ")); w = Cdhc_anderson_darling(z, nsites); fprintf(stdout, "%g %g\n", w[1], w[0]); break; case 9: /* Cramer-Von Mises */ fprintf(stdout, _("Cramer-Von Mises W^2(regular & modified for normality): ")); w = Cdhc_cramer_von_mises(z, nsites); fprintf(stdout, "%g %g\n", w[1], w[0]); break; case 10: /* Kolmogorov-Smirnov */ fprintf(stdout, _("Kolmogorov-Smirnov's D (regular & modified for normality): ")); w = Cdhc_kolmogorov_smirnov(z, nsites); fprintf(stdout, "%g %g\n", w[1], w[0]); break; case 11: /* chi-square */ fprintf(stdout, _("Chi-Square stat (equal probability classes) and d.f.: ")); w = Cdhc_chi_square(z, nsites); fprintf(stdout, "%g %d\n", w[0], (int)w[1]); break; case 12: /* Shapiro-Wilk */ if (nsites > 50) { G_warning(_("Shapiro-Wilk's W cannot be used for n > 50")); if (nsites < 99) G_message(_("Use Weisberg-Binghams's W''")); } else { fprintf(stdout, _("Shapiro-Wilk W: ")); w = Cdhc_shapiro_wilk(z, nsites); fprintf(stdout, "%g\n", w[0]); } break; case 13: /* Weisberg-Bingham */ if (nsites > 99 || nsites < 50) G_warning(_("Weisberg-Bingham's W'' cannot be used for n < 50 or n > 99")); else { fprintf(stdout, _("Weisberg-Bingham's W'': ")); w = Cdhc_weisberg_bingham(z, nsites); fprintf(stdout, "%g\n", w[0]); } break; case 14: /* Royston */ if (nsites > 2000) G_warning(_("Royston only extended Shapiro-Wilk's W up to n = 2000")); else { fprintf(stdout, _("Shapiro-Wilk W'': ")); w = Cdhc_royston(z, nsites); fprintf(stdout, "%g\n", w[0]); } break; case 15: /* Kotz */ fprintf(stdout, _("Kotz' T'_f (Lognormality vs. Normality): ")); w = Cdhc_kotz_families(z, nsites); fprintf(stdout, "%g\n", w[0]); break; default: break; } } exit(EXIT_SUCCESS); }