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; }
/* * 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; }
VARIOGRAM *reml_sills(DATA *data, VARIOGRAM *vp) { int i, j, k; MAT **Vk = NULL, *X = MNULL; VEC *Y = VNULL, *init = VNULL; DPOINT *dpa, *dpb; double dx, dy = 0.0, dz = 0.0, dzero2; if (data == NULL || vp == NULL) ErrMsg(ER_NULL, "reml()"); select_at(data, (DPOINT *) NULL); if (vp->n_models <= 0) ErrMsg(ER_VARNOTSET, "reml: please define initial variogram model"); /* * create Y, X, Vk's only once: */ Y = get_y(&data, Y, 1); X = get_X(&data, X, 1); Vk = (MAT **) emalloc(vp->n_models * sizeof(MAT)); init = v_resize(init, vp->n_models); for (i = 0; i < vp->n_models; i++) { init->ve[i] = vp->part[i].sill; /* remember init. values for updating */ vp->part[i].sill = 1; Vk[i] = m_resize(MNULL, X->m, X->m); } dzero2 = gl_zero * gl_zero; for (i = 0; i < data->n_list; i++) { for (j = 0; j < vp->n_models; j++) /* fill diagonals */ Vk[j]->me[i][i] = Covariance(vp->part[j], 0.0, 0.0, 0.0); for (j = 0; j < i; j++) { /* off-diagonal elements: */ dpa = data->list[i]; dpb = data->list[j]; /* * if different points coincide on a locations, shift them, * or the covariance matrix will become singular */ dx = dpa->x - dpb->x; dy = dpa->y - dpb->y; dz = dpa->z - dpb->z; if (data->pp_norm2(dpa, dpb) < dzero2) { if (data->mode & X_BIT_SET) dx = (dx >= 0 ? gl_zero : -gl_zero); if (data->mode & Y_BIT_SET) dy = (dy >= 0 ? gl_zero : -gl_zero); if (data->mode & Z_BIT_SET) dz = (dz >= 0 ? gl_zero : -gl_zero); } for (k = 0; k < vp->n_models; k++) Vk[k]->me[i][j] = Vk[k]->me[j][i] = Covariance(vp->part[k], dx, dy, dz); } } if (reml(Y, X, Vk, vp->n_models, gl_iter, gl_fit_limit, init)) vp->ev->refit = 0; else /* on convergence */ pr_warning("no convergence while fitting variogram"); for (i = 0; i < vp->n_models; i++) vp->part[i].sill = init->ve[i]; update_variogram(vp); if (DEBUG_VGMFIT) logprint_variogram(vp, 1); for (i = 0; i < vp->n_models; i++) m_free(Vk[i]); efree(Vk); m_free(X); v_free(Y); v_free(init); return vp; }
void predict_all(DATA **data) { int i = 0, random_path = 0; DPOINT *here = NULL, *where = NULL; PRED_AT at_what; unsigned int row, col; n_done = 0; val_data = get_dataval(); if (val_data->id > -1) { at_what = AT_POINTS; n_pred_locs = val_data->n_list; if (val_data->colns) strata_min = val_data->minstratum; } else if (get_n_masks() > 0) { at_what = AT_GRIDMAP; here = (DPOINT *) emalloc(sizeof(DPOINT)); here->u.stratum = -2; /* only NON-MV cells */ if (max_block_dimension(0) > 0.0) SET_BLOCK(here); else SET_POINT(here); } else /* what else ? */ return; if (at_what == AT_GRIDMAP && get_n_outfile() == 0) { pr_warning("no output maps defined"); return; } init_predictions(at_what); if (at_what == AT_GRIDMAP && !data[0]->dummy) { if (data[0]->maxX < masks[0]->x_ul || data[0]->minX > (masks[0]->x_ul + masks[0]->cols * masks[0]->cellsizex) || data[0]->minY > masks[0]->y_ul || data[0]->maxY < (masks[0]->y_ul - masks[0]->rows * masks[0]->cellsizey)) { pr_warning("ALL data are outside the map boundaries"); printlog("data x[%g,%g], y[%g,%g]; map x[%g,%g], y[%g,%g]\n", data[0]->minX, data[0]->maxX, data[0]->minY, data[0]->maxY, masks[0]->x_ul, masks[0]->x_ul + masks[0]->cols * masks[0]->cellsizex, masks[0]->y_ul - masks[0]->rows * masks[0]->cellsizey, masks[0]->y_ul ); } else if (map_xy2rowcol(masks[0], data[0]->minX, data[0]->minY, &row, &col) || map_xy2rowcol(masks[0], data[0]->maxX, data[0]->minY, &row, &col) || map_xy2rowcol(masks[0], data[0]->minX, data[0]->maxY, &row, &col) || map_xy2rowcol(masks[0], data[0]->maxX, data[0]->maxY, &row, &col)) pr_warning("at least some data are outside the map boundaries"); /* this is not a sufficient test! */ } if (gl_rp) /* Yes, by default */ random_path = is_simulation(get_method()); row = col = 0; while ((where = next_location(here, at_what, random_path, &row, &col, data)) != NULL) { for (i = 0; i < get_n_outfile(); i++) set_mv_double(&(est[i])); /* initialize estimates */ if (where->u.stratum >= 0) { if (get_mode() != STRATIFY) { for (i = 0; i < get_n_vars(); i++) select_at(data[i], where); } else if (where->u.stratum < get_n_vars()) select_at(data[where->u.stratum], where); get_est(data, get_method(), where, est); } /* printf("%g %g\n", est[0], est[1]); */ write_output(est, at_what, where, row, col); } exit_predictions(at_what); if (here != NULL) efree(here); print_progress(100, 100); }