int push_to_merge_table(DATA *d, int to_var, int col_this_X, int col_other_X) { int i; DATA **data; data = get_gstat_data(); if (to_var >= d->id) { /* should not occur by construction */ pr_warning("use push_to_merge_table only backwards (%d >= %d)", to_var, d->id); return 1; } if (col_this_X >= d->n_X || col_other_X >= data[to_var]->n_X) { pr_warning("merge error: variable out of range"); return 1; } if (d->beta || data[to_var]->beta) { pr_warning("cannot merge to or from fixed (known) parameters"); return 1; } for (i = 0; i < d->n_merge; i++) { if (col_this_X == d->mtbl[i].col_this_X) { pr_warning("merge error: cannot merge column twice"); return 1; } } d->n_merge++; d->mtbl = (MERGE_TABLE *) erealloc(d->mtbl, d->n_merge * sizeof (MERGE_TABLE)); d->mtbl[d->n_merge - 1].to_var = to_var; d->mtbl[d->n_merge - 1].col_this_X = col_this_X; d->mtbl[d->n_merge - 1].col_other_X = col_other_X; return 0; }
int fit_variogram(VARIOGRAM *v) { DATA **d = NULL; int i = 0; long n = 0; if (v->ev->refit == 0) return 0; if (v->ev->fit != NO_FIT && v->ev->fit != MIVQUE_FIT) { if (! v->ev->cloud) { while (n == 0 && i < v->ev->n_est) n += v->ev->nh[i++]; /* check if estimates exist */ if (n == 0) /* bad luck */ return 1; } else if (v->ev->n_est == 0) return 1; } if (v->ev->fit != NO_FIT) { if (v->ev->map) ErrMsg(ER_IMPOSVAL, "cannot fit model to variogram map"); for (i = 0; i < v->n_models; i++) if (v->part[i].sill == 0.0 && v->part[i].fit_sill != 0) v->part[i].sill = 1.0; /* avoid lot'o trouble */ } if (v->ev->fit == WLS_FIT_MOD && !is_variogram(v)) pr_warning("this fit method is not recommended for covariograms"); v->ev->direction.x = sin(gl_alpha * PI / 180.0) * cos(gl_beta * PI / 180.0); v->ev->direction.y = cos(gl_alpha * PI / 180.0) * cos(gl_beta * PI / 180.0); v->ev->direction.z = sin(gl_beta * PI / 180.0); switch (v->ev->fit) { case NO_FIT: break; case OLS_FIT: /* BREAKTHROUGH: */ case WLS_FIT: /* BREAKTHROUGH: */ case WLS_FIT_MOD: case WLS_NHH: wls_fit(v); break; case MIVQUE_FIT: if (v->id1 != v->id2) return 1; d = get_gstat_data(); reml_sills(d[v->id1], v); break; default: Rprintf("%d\n", v->ev->fit); ErrMsg(ER_IMPOSVAL, "fit_vgm(): value for fit not recognized or not implemented"); /* no default: force compile warning on missing option! */ } return 0; }
static void exit_predictions(PRED_AT what) { int i; if (gl_nsim > 1 && gl_lhs) lhs(get_gstat_data(), get_n_vars(), get_mode() == STRATIFY); switch (what) { case AT_POINTS: write_points(NULL, NULL, NULL, NULL, 0); if (gl_nsim > 1) save_simulations_to_ascii(o_filename); break; case AT_GRIDMAP: if (gl_nsim > 1) { if (DEBUG_DUMP) printlog("\nWriting results to files..."); save_simulations_to_maps(masks[0]); if (DEBUG_DUMP) printlog("done"); } else { for (i = 0; i < get_n_outfile(); i++) { if (get_outfile_namei(i)) { map_sign(outmap[i], what_is_outfile(i)); (outmap[i])->write(outmap[i]); map_free(outmap[i]); } } } for (i = 0; i < get_n_masks(); i++) map_free(masks[i]); efree(masks); efree(outmap); break; } print_orvc(); efree(est); } /* exit_predictions() */
int ossfim(int argc, char *argv[]) { int c, n = 25, dx = 9, dy = 9, i, j, plot_vgm = 0; double b = 1, B = 10, s = 1, S = 10, blocksize, samplespacing, est[2], **table; DATA **d = NULL; DPOINT *block = NULL, where; char *vgm_str = "1Exp(10)", *map_name = NULL; VARIOGRAM *vgm; while ((c = getopt(argc, argv, "n:m:B:b:S:s:V:v:x:y:")) != EOF) { switch (c) { case 'n': if (read_int(optarg, &n) || n <= 0) ErrMsg(ER_ARGOPT, "n"); break; case 'b': if (read_double(optarg, &b) || b < 0) ErrMsg(ER_ARGOPT, "b"); break; case 'B': if (read_double(optarg, &B) || B <= 0) ErrMsg(ER_ARGOPT, "B"); break; case 's': if (read_double(optarg, &s) || s <= 0) ErrMsg(ER_ARGOPT, "s"); break; case 'S': if (read_double(optarg, &S) || S <= 0) ErrMsg(ER_ARGOPT, "S"); break; case 'x': if (read_int(optarg, &dx) || dx <= 0) ErrMsg(ER_ARGOPT, "x"); break; case 'y': if (read_int(optarg, &dy) || dy <= 0) ErrMsg(ER_ARGOPT, "y"); break; case 'v': vgm_str = optarg; break; case 'V': plot_vgm = 1; vgm_str = optarg; break; case 'm': map_name = optarg; break; default: ErrClo(optopt); break; } } which_identifier("dummy grid"); d = get_gstat_data(); init_one_data(d[0]); d[0]->id = 0; d[0]->n_list = d[0]->n_max = 0; d[0]->mode = X_BIT_SET | Y_BIT_SET | V_BIT_SET; set_norm_fns(d[0]); vgm = get_vgm(0); if (read_variogram(vgm, vgm_str)) ErrMsg(ER_SYNTAX, vgm_str); vgm->ev->evt = SEMIVARIOGRAM; vgm->id1 = vgm->id2 = d[0]->id; block = get_block_p(); block->z = 0.0; block->x = block->y = -1.0; est[0] = 0.0; est[1] = -1.0; where.x = where.y = where.z = 0.0; where.X = (double *) emalloc(sizeof(double)); where.X[0] = 1.0; if (plot_vgm) return fprint_gnuplot_variogram(stdout, vgm, "", GIF, 0); table = (double **) emalloc((dy + 1) * sizeof(double *)); for (i = 0; i <= dy; i++) table[i] = (double *) emalloc((dx + 1) * sizeof(double)); /* do it: */ for (i = 0; i <= dx; i++) { /* sample spacing loop */ samplespacing = s + (i / (1.0 * dx)) * (S - s); generate_grid(d[0], samplespacing, n); select_at(d[0], &where); for (j = 0; j <= dy; j++) { /* block sizes loop */ reset_block_discr(); vgm_init_block_values(vgm); blocksize = b + (j / (1.0 * dy)) * (B - b); block->x = block->y = blocksize; if (blocksize == 0.0) SET_POINT(&where); else SET_BLOCK(&where); gls(d, 1, GLS_BLUP, &where, est); if (map_name) table[i][j] = sqrt(est[1]); else printlog("%g %g %g\n", samplespacing, blocksize, sqrt(est[1])); } } if (map_name) ossfim2map(table, map_name, s, S, b, B, dx, dy); return 0; }
int main(int argc, char *argv[]) { #else int gstat_main(int argc, char *argv[]) { #endif DATA **data = NULL, *valdata = NULL; /* * initialise some global variables: */ atexit(close_gstat_log_file); init_userio(1); init_global_variables(); argv0 = argv[0]; /* * register command line arguments on command_line: */ command_line = store_argv(argc, argv); parse_gstatrc(); /* * INPUT: command line options; */ parse_options(argc, argv); /* exits on -e options */ /* * start with program heading: */ printlog("%s: %s version %s\n", GSTAT_NAME, GSTAT_OS, VERSION); printlog("%s\n", GSTAT_CR); gstat_start(); /* * INPUT: Parse command files: */ if (optind == argc) { /* there's no command file name left */ if (get_method() != UIF) { /* -i or -m command line option */ /* no arguments */ printlog("Updates, manuals and source code: %s\n", GSTAT_HOME); printlog("%s\n", USAGE); ErrMsg(ER_NOCMD, ""); } else { start_ui(); exit(0); } } else { /* we have a command file to be read */ for ( ; optind < argc; optind++) { command_file_name = argv[optind]; parse_file(command_file_name); if (logfile_name != NULL) set_gstat_log_file(efopen(logfile_name, "w")); /* * get global variables locally: */ data = get_gstat_data(); valdata = get_dataval(); set_seed(gl_seed); /* * check variable settings and next * INPUT: read data values from file: */ read_all_data(data, valdata, get_n_vars()); if (get_method() == NSP) /* Still no answer to this: choose default */ set_method(get_default_method()); set_mode(); check_global_variables(); setup_meschach_error_handler(); if (DEBUG_DUMP) dump_all(); if (get_method() != NSP) printlog("[%s]\n", method_string(get_method())); if (check_only) set_method(NSP); /* * start calculations && OUTPUT routines: */ switch (get_method()) { case UIF: start_ui(); break; case NSP: break; case COV: case SEM: do_variogram(get_n_vars(), get_method()); break; case POLY: setup_poly_method(); /*FALLTHROUGH*/ default: if (gl_xvalid) /* validation/cross validation */ cross_valid(data); else predict_all(data); /* or prediction: */ break; } /* switch get_method() */ remove_all(); /* free all data etc. */ init_global_variables(); /* re-init for next round */ } } if (DEBUG_DUMP) atexit(print_file_record); if (get_method() != UIF) elapsed(); /* * file closing & data freeing time: */ if (plotfile != NULL) efclose(plotfile); exit(0); } /* end of main() */
/* * calculate sample variogram from data * calculates variogram, crossvariogram, covariogram or crosscovariogram * from sample data. Data are obtained from the central data base (glvars) * using get_gstat_data(), and the variogram requested is that of data id * v->id1 and v->id2 -- a direct (co) variogram when id1 == id2, a cross * (co) variogram when id1 != id2. * * if v->fname is set and (one of) id1 or id2 is a dummy data, the * actual sample variogram is not calculated but rather read from the * file v->fname. This is done to enable separate sample variogram * calculation (in batch or on a fast remote computer) and model fitting * (e.g. on the desk top). * * returns: non-zero if writing sample variogram to file failed. */ int calc_variogram(VARIOGRAM *v /* pointer to VARIOGRAM structure */, const char *fname /* pointer to output file name, or NULL if no output has to be written to file */ ) { DATA **d = NULL, *d1 = NULL, *d2 = NULL; FILE *f = NULL; assert(v); d = get_gstat_data(); d1 = d[v->id1]; d2 = d[v->id2]; if (v->fname && (d1->dummy || d2->dummy)) { if ((v->ev = load_ev(v->ev, v->fname)) == NULL) ErrMsg(ER_READ, "could not read sample variogram"); v->ev->cloud = 0; v->ev->recalc = 0; return 0; } if (d1->sel == NULL) select_at(d1, NULL); /* global selection (sel = list) */ if (d2->sel == NULL) select_at(d2, NULL); if (v->ev->evt == CROSSVARIOGRAM && v->ev->is_asym == -1) { /* v's first time */ if (coordinates_are_equal(d[v->id1], d[v->id2])) v->ev->pseudo = 0; else v->ev->pseudo = 1; if (gl_sym_ev == 0) v->ev->is_asym = v->ev->pseudo; /* pseudo: always, else: only if set */ else v->ev->is_asym = 0; } if (gl_zero_est == ZERO_DEFAULT) { /* choose a suitable default */ if (is_covariogram(v)) v->ev->zero = ZERO_SPECIAL; else { /* v is variogram */ if (v->ev->pseudo) v->ev->zero = ZERO_SPECIAL; else v->ev->zero = ZERO_INCLUDE; } } else v->ev->zero = zero_int2enum(gl_zero_est); assert(v->ev->zero != ZERO_DEFAULT); fill_cutoff_width(d1, v); if (v->ev->map && v->fname == NULL) return -1; v->ev->cloud = (v->ev->iwidth <= 0.0); if (v->ev->cloud && (d[v->id1]->n_list >= MAX_NH || d[v->id2]->n_list >= MAX_NH)) pr_warning("observation numbers in cloud will be wrong"); set_direction_values(gl_alpha, gl_beta, gl_tol_hor, gl_tol_ver); v->ev->is_directional = is_directional(v); if (v->ev->recalc) { switch (v->ev->evt) { case SEMIVARIOGRAM: semivariogram(d[v->id1], v->ev); break; case CROSSVARIOGRAM: cross_variogram(d[v->id1], d[v->id2], v->ev); break; case COVARIOGRAM: v->ev->is_asym = gl_sym_ev; covariogram(d[v->id1], v->ev); break; case CROSSCOVARIOGRAM: cross_covariogram(d[v->id1], d[v->id2], v->ev); break; case NOTSPECIFIED: default: assert(0); /* aborts */ break; } } if (v->ev->map) ev2map(v); else if (fname != NULL) { f = efopen(fname, "w"); fprint_header_vgm(f, d[v->id1], d[v->id2], v->ev); fprint_sample_vgm(f, v->ev); } if (f && f != stdout) return efclose(f); return 0; }
static void init_predictions(PRED_AT w) { int i; DPOINT *bp; DATA **d = NULL; #ifdef WITH_SPIRAL DATA_GRIDMAP *grid; #endif est = (double *) emalloc(get_n_outfile() * sizeof(double)); bp = get_block_p(); d = get_gstat_data(); switch (w) { case AT_POINTS: if (o_filename == NULL) ErrMsg(ER_VARNOTSET, "please specify output file"); write_points(o_filename, val_data, NULL, NULL, get_n_outfile()); if (bp->x == -1.0) { /* set default */ bp->x = bp->y = 1.0; pr_warning("default block size set to: dx=1, dy=1"); } break; case AT_GRIDMAP: /* open mask files: */ get_maskX(NULL, NULL, 0, 0); /* re-initializes static arrays */ masks = (GRIDMAP **) emalloc(get_n_masks() * sizeof(GRIDMAP *)); for (i = 0; i < get_n_masks(); i++) masks[i] = check_open(get_mask_name(i), i); /* read as float */ if (n_pred_locs > 0) strata_min = floor(masks[0]->cellmin); outmap = (GRIDMAP **) emalloc(get_n_outfile() * sizeof(GRIDMAP *)); printlog("initializing maps "); for (i = 0; i < get_n_outfile(); i++) { if (get_outfile_namei(i) != NULL) { printlog("."); /* creating maps ..... */ if (get_method() == ISI) masks[0]->celltype = CT_UINT8; outmap[i] = map_dup(get_outfile_namei(i), masks[0]); } else outmap[i] = NULL; } printlog("\n"); if (bp->x == -1.0) { /* set default to map cellsize */ bp->x = masks[0]->cellsizex; bp->y = masks[0]->cellsizey; pr_warning("default block size set to dx=%g, dy=%g", bp->x, bp->y); } for (i = 0; i < get_n_vars(); i++) { if (d[i]->dummy) { d[i]->minX = masks[0]->x_ul + 0.5 * masks[0]->cellsizex; d[i]->maxX = masks[0]->x_ul + masks[0]->cellsizex * (masks[0]->cols - 0.5); d[i]->maxY = masks[0]->y_ul - 0.5 * masks[0]->cellsizey; d[i]->minY = masks[0]->y_ul - masks[0]->cellsizey * (masks[0]->rows - 0.5); d[i]->minZ = d[i]->maxZ = 0.0; } if (d[i]->togrid) datagrid_rebuild(d[i], 1); } break; } if (gl_nsim > 1) init_simulations(d); if (is_simulation(get_method()) && get_n_beta_set() != get_n_vars()) setup_beta(d, get_n_vars(), gl_nsim); } /* init_predictions() */