void logprint_variogram(const VARIOGRAM *v, int verbose) { /* prints contents of VARIOGRAM v to R console */ if (v->id1 < 0 && v->id2 < 0) return; /* never set */ if (v->id1 == v->id2) Rprintf("variogram(%s):\n", name_identifier(v->id1)); else Rprintf("variogram(%s,%s):\n", name_identifier(v->id1), name_identifier(v->id2)); for (int i = 0; i < v->n_models; i++) { Rprintf("# model: %d type: %s sill: %g range: %g\n", i, v_models[v->part[i].model].name_long, v->part[i].sill, v->part[i].range[0]); if (v->part[i].tm_range != NULL) { Rprintf("# range anisotropy, rotation matrix:\n"); for (int j = 0; j < 3; j++) { for (int k = 0; k < 3; k++) Rprintf("%s%8.4f", k == 0 ? "# " : " ", v->part[i].tm_range->tm[j][k]); Rprintf("\n"); } } } Rprintf( "# sum sills %g, max %g, min %g, flat at distance %g\n", v->sum_sills, v->max_val, v->min_val, v->max_range); return; }
void fprint_header_vgm(FILE *f, const DATA *a, const DATA *b, const SAMPLE_VGM *ev) { time_t tp; char *cp = NULL; /* char *pwd; */ fprintf(f, "#gstat %s %s [%s]", GSTAT_OS, VERSION, command_line); /* if (pwd = getenv("PWD")) fprintf(f, "(in %s)", pwd); */ fprintf(f, "\n"); fprintf(f, "#sample %s%s\n", (ev->evt == CROSSVARIOGRAM && ev->pseudo ? "pseudo " : ""), vgm_type_str[ev->evt]); tp = time(&tp); fprintf(f, "#%s", asctime(localtime(&tp))); /* includes \n */ cp = print_data_line(a, &cp); fprintf(f, "#data(%s): %s", name_identifier(a->id), cp); if (a != b) { cp = print_data_line(b, &cp); fprintf(f, " data(%s): %s", name_identifier(b->id), cp); } if (cp != NULL) efree(cp); fprintf(f, "\n"); fprintf(f, "#[1] mean: %g variance: %g", a->mean, a->std * a->std); if (a != b) fprintf(f, " [2] mean: %g variance: %g", b->mean, b->std * a->std); fprintf(f, "\n"); fprintf(f, "#cutoff: %g ", ev->cutoff); if (gl_bounds == NULL) fprintf(f,"%s %g\n","interval width:", ev->iwidth); else fprintf(f, "(fixed boundaries)\n"); if (! ev->is_directional) fprintf(f, "#direction: total "); else fprintf(f,"#alpha <x,y>: %gd +/- %g; beta <alpha,z>: %gd +/- %g%s", gl_alpha, gl_tol_hor, gl_beta, gl_tol_ver, gl_sym_ev ? " (symmetric)" : ""); fprintf(f, "\n"); if (ev->cloud) fprintf(f, "# i j distance %s cloud\n", vgm_type_str[ev->evt]); else fprintf(f, "# from to n_pairs av_dist %s\n", vgm_type_str[ev->evt]); return; }
static void init_qtree(DATA *d) { /* * Initialize the root of the search tree: * Since this is called from qtree_quick_select(), we know allready * quite a lot. This helps choosing sensible values for the * top level bbox origin and size. */ const GRIDMAP *gt = NULL; DATA *simlocs = NULL; int i, mode; BBOX bbox; if (is_simulation(get_method())) { /* * sequential simulation: the simulation path (through DATA or GRIDMAP) * will make up for (most of) the search locations: */ gt = (const GRIDMAP *) NULL; simlocs = get_dataval(); /* in case of simulation one of them will be non-NULL */ } /* get initial estimate for top level bbox: */ if (gt != NULL) bbox = bbox_from_grid(gt, NULL); else if (simlocs != NULL) { bbox = bbox_from_data(simlocs); if (bbox.size <= 0.0) bbox = bbox_from_data(d); } else bbox = bbox_from_data(d); if (bbox.size <= 0.0) bbox = bbox_from_data(get_dataval()); if (bbox.size <= 0.0) ErrMsg(ER_IMPOSVAL, "bbox with zero size: remove neighbourhood settings?"); init_qnode(&(d->qtree_root), d->n_list < gl_split, bbox); /* ML1 */ mode = bbox.mode; for (i = 0; i < d->n_list; i++) qtree_push_point(d, d->list[i]); /* now they won't be rejected */ /* ML2 */ if (DEBUG_DUMP) { printlog("top level search tree statistics for data(%s):\n", name_identifier(d->id)); printlog("bounding box origin ["); if (mode & X_BIT_SET) printlog("%g", d->qtree_root->bb.x); if (mode & Y_BIT_SET) printlog(",%g", d->qtree_root->bb.y); if (mode & Z_BIT_SET) printlog(",%g", d->qtree_root->bb.z); printlog("]; dimension %g\n", d->qtree_root->bb.size); } /* qtree_print(d); */ return; }
void qtree_print(DATA *d) { /* * plot the full tree (2D), in a format that can be read by jgraph, found * at netlib or at http://kenner.cs.utk.edu/~plank/plank/jgraph/jgraph.html */ printlog("newgraph\nxaxis size 3\nyaxis size 3\n"); printlog("title : %s [n = %d]\n", name_identifier(d->id), d->n_list); logprint_qtree(d->qtree_root, 0); return; }
static void read_all_data(DATA **data, DATA *valdata, int n_vars) { int i; DATA *area; init_data_minmax(); area = get_data_area(); for (i = 0; i < n_vars; i++) { if (get_mode() == STRATIFY) printlog("stratum # %d:\n", i + strata_min); printlog("data(%s): ", name_identifier(i)); if (data[i]->id < 0) { message("data(%s) was not specified\n", name_identifier(i)); ErrMsg(ER_SYNTAX, "data specification error"); } read_gstat_data(data[i]); report_data(data[i]); } /* for i */ /* * what to do when area is specified, but no masks or data()? * default prediction to `area'. Create a valdata with one point at * centre of area (for select()); and centre area at (0,0,0) */ if (area && get_n_masks() <= 0 && valdata->id == -1) { valdata->id = ID_OF_VALDATA; valdata->centre = area->centre = 1; } /* * read data() data: */ if (valdata->id > -1) { setup_valdata_X(valdata); if (! valdata->centre) valdata = read_gstat_data(valdata); } /* * read area, if existed */ if (area != NULL && get_method() != POLY) { read_gstat_data(area); /* now, before centring area: */ if (valdata->centre) valdata = get_area_centre(area, valdata); if (area->centre) centre_area(area); printlog("area:%s\n", area->centre ? " (centred around 0)" : ""); report_data(area); if (DEBUG_DATA) print_data_list(area); } /* * read edges, if existed */ if (get_n_edges() > 0) { read_edges(); report_edges(); /* setup_visibility_graph(); */ /*setup_planar_subdivisions();*/ } /* * setup and report data */ if (valdata->id > -1) { printlog("data():%s ", valdata->centre ? " [at area centre]" : ""); report_data(valdata); } for (i = 0; i < n_vars; i++) setup_data_minmax(data[i]); if (valdata->id > -1) setup_data_minmax(valdata); for (i = 0; i < n_vars; i++) calc_polynomials(data[i]); if (valdata->id > -1) calc_polynomials(valdata); if (DEBUG_DATA) { for (i = 0; i < n_vars; i++) print_data_list(data[i]); if (valdata->id > -1) print_data_list(valdata); } }
void check_global_variables(void) { /* * Purpose : check internal variable consistency, add some parameters * Created by : Edzer J. Pebesma * Date : april 13, 1992 * Prerequisites : none * Returns : - * Side effects : none * also check Cauchy-Schwartz unequality on cross/variograms. */ int i, j, nposX, n_merge = 0; METHOD m; VARIOGRAM *v_tmp; /* UK: check if n_masks equals total nr of unbiasedness cond. */ if (gl_nblockdiscr < 2) ErrMsg(ER_RANGE, "nblockdiscr must be >= 2"); if (method == UKR || method == LSLM) { nposX = 0; for (i = 0; i < get_n_vars(); i++) for (j = 0; j < data[i]->n_X; j++) { if (data[i]->colX[j] > 0) nposX++; } } if (method == SPREAD) { for (i = 0; i < get_n_vars(); i++) if (data[i]->sel_rad == DBL_MAX) data[i]->sel_rad *= 0.99; /* force distance calculation */ } if (get_n_beta_set() != 0 && get_n_beta_set() != get_n_vars()) ErrMsg(ER_SYNTAX, "set sk_mean or beta either for all or for no variables"); if (!(method == ISI || method == GSI)) { if (gl_nsim > 1) ErrMsg(ER_IMPOSVAL, "nsim only allowed for simulation"); } if (method == ISI && max_block_dimension(0) > 0.0) ErrMsg(ER_IMPOSVAL, "indicator simulation only for points"); /* * check if both block and area are set */ if (data_area != NULL && (block.x > 0.0 || block.y > 0.0 || block.z > 0.0)) ErrMsg(ER_IMPOSVAL, "both block and area set: choose one"); /* * check for equality of coordinate dimensions: */ for (i = 1; i < get_n_vars(); i++) { if ((data[i]->mode & V_BIT_SET) != (data[0]->mode & V_BIT_SET)) { message("data(%s) and data(%s):\n", name_identifier(0), name_identifier(i)); ErrMsg(ER_IMPOSVAL, "data have different coordinate dimensions"); } } if (valdata->id > 0 && data[0]->dummy == 0 && ((data[0]->mode | (V_BIT_SET | S_BIT_SET)) != (valdata->mode | (V_BIT_SET | S_BIT_SET)))) { message("data() and data(%s):\n", name_identifier(0)); ErrMsg(ER_IMPOSVAL, "data have different coordinate dimensions"); for (i = 0; i < get_n_vars(); i++) { if (data[i]->dummy) { data[i]->mode = (valdata->mode | V_BIT_SET); data[i]->minX = valdata->minX; data[i]->minY = valdata->minY; data[i]->minZ = valdata->minZ; data[i]->maxX = valdata->maxX; data[i]->maxY = valdata->maxY; data[i]->maxZ = valdata->maxZ; set_norm_fns(data[i]); } } } for (i = 0; i < get_n_vars(); i++) { if (data[i]->fname == NULL && !data[i]->dummy) { message("file name for data(%s) not set\n", name_identifier(i)); ErrMsg(ER_NULL, " "); } if (data[i]->id < 0) { message("data(%s) not set\n", name_identifier(i)); ErrMsg(ER_NULL, " "); } if (data[i]->beta && data[i]->beta->size != data[i]->n_X) { pr_warning("beta dimension (%d) should equal n_X (%d)", data[i]->beta->size, data[i]->n_X); ErrMsg(ER_IMPOSVAL, "sizes of beta and X don't match"); } if (data[i]->sel_rad == DBL_MAX && data[i]->oct_max > 0) ErrMsg(ER_IMPOSVAL, "define maximum search radius (rad) for octant search"); if (data[i]->vdist && data[i]->sel_rad == DBL_MAX) ErrMsg(ER_IMPOSVAL, "when using vdist, radius should be set"); if (! data[i]->dummy && ! (data[i]->mode & V_BIT_SET)) { message("no v attribute set for data(%s)\n", name_identifier(data[i]->id)); ErrMsg(ER_NULL, " "); } if (method != SEM && method != COV) { /* check neighbourhood settings */ if (data[i]->sel_rad < 0.0 || data[i]->sel_min < 0 || data[i]->sel_max < 0 || (data[i]->sel_min > data[i]->sel_max)) { message( "invalid neighbourhood selection: radius %g max %d min %d\n", data[i]->sel_rad, data[i]->sel_max, data[i]->sel_min); ErrMsg(ER_IMPOSVAL, " "); } } if (data[i]->id > -1 && (method == OKR || method == SKR || is_simulation(method) || method == UKR)) { if (vgm[LTI(i,i)] == NULL || vgm[LTI(i,i)]->id < 0) { message("variogram(%s) not set\n", name_identifier(i)); ErrMsg(ER_VARNOTSET, "variogram()"); } } n_merge += data[i]->n_merge; } if (n_merge && get_mode() != MULTIVARIABLE) ErrMsg(ER_IMPOSVAL, "merge only works in multivariable mode"); if (mode == SIMPLE && get_method() != UIF) { /* check if it's clean: */ for (i = 0; i < get_n_vars(); i++) for (j = 0; j < i; j++) if (vgm[LTI(i,j)] != NULL && vgm[LTI(i,j)]->id > 0) { message("variogram(%s, %s) %s\n", name_identifier(i), name_identifier(j), "can only be set for ck, cs, uk, sk, ok, sem or cov"); ErrMsg(ER_IMPOSVAL, "variogram()"); } } if ((m = get_default_method()) != get_method()) { if (m == UKR && (get_method() == OKR || get_method() == SKR)) ErrMsg(ER_IMPOSVAL, "\nremove X=... settings for ordinary or simple kriging"); if (m == OKR && get_method() == SKR) ErrMsg(ER_IMPOSVAL, "method: something's terribly wrong!"); if (m == OKR && get_method() == UKR) { message("I would recommend:\n"); message("Do not specify uk if ok is all you'll get\n"); } } if (mode == MULTIVARIABLE && get_method() != UIF && get_method() != SEM && get_method() != COV && n_variograms_set() > 0) check_variography((const VARIOGRAM **) vgm, get_n_vars()); v_tmp = init_variogram(NULL); free_variogram(v_tmp); }
void check_variography(const VARIOGRAM **v, int n_vars) /* * check for intrinsic correlation, linear model of coregionalisation * or else (with warning) Cauchy Swartz */ { int i, j, k, ic = 0, lmc, posdef = 1; MAT **a = NULL; double b; char *reason = NULL; if (n_vars <= 1) return; /* * find out if lmc (linear model of coregionalization) hold: * all models must have equal base models (sequence and range) */ for (i = 1, lmc = 1; lmc && i < get_n_vgms(); i++) { if (v[0]->n_models != v[i]->n_models) { reason = "number of models differ"; lmc = 0; } for (k = 0; lmc && k < v[0]->n_models; k++) { if (v[0]->part[k].model != v[i]->part[k].model) { reason = "model types differ"; lmc = 0; } if (v[0]->part[k].range[0] != v[i]->part[k].range[0]) { reason = "ranges differ"; lmc = 0; } } for (k = 0; lmc && k < v[0]->n_models; k++) if (v[0]->part[k].tm_range != NULL) { if (v[i]->part[k].tm_range == NULL) { reason = "anisotropy for part of models"; lmc = 0; } else if ( v[0]->part[k].tm_range->ratio[0] != v[i]->part[k].tm_range->ratio[0] || v[0]->part[k].tm_range->ratio[1] != v[i]->part[k].tm_range->ratio[1] || v[0]->part[k].tm_range->angle[0] != v[i]->part[k].tm_range->angle[0] || v[0]->part[k].tm_range->angle[1] != v[i]->part[k].tm_range->angle[1] || v[0]->part[k].tm_range->angle[2] != v[i]->part[k].tm_range->angle[2] ) { reason = "anisotropy parameters are not equal"; lmc = 0; } } else if (v[i]->part[k].tm_range != NULL) { reason = "anisotropy for part of models"; lmc = 0; } } if (lmc) { /* * check for ic: */ a = (MAT **) emalloc(v[0]->n_models * sizeof(MAT *)); for (k = 0; k < v[0]->n_models; k++) a[k] = m_get(n_vars, n_vars); for (i = 0; i < n_vars; i++) { for (j = 0; j < n_vars; j++) { /* for all variogram triplets: */ for (k = 0; k < v[0]->n_models; k++) ME(a[k], i, j) = v[LTI(i,j)]->part[k].sill; } } /* for ic: a's must be scaled versions of each other: */ ic = 1; for (k = 1, ic = 1; ic && k < v[0]->n_models; k++) { b = ME(a[0], 0, 0)/ME(a[k], 0, 0); for (i = 0; ic && i < n_vars; i++) for (j = 0; ic && j < n_vars; j++) if (fabs(ME(a[0], i, j) / ME(a[k], i, j) - b) > EPSILON) ic = 0; } /* check posdef matrices */ for (i = 0, lmc = 1, posdef = 1; i < v[0]->n_models; i++) { posdef = is_posdef(a[i]); if (posdef == 0) { reason = "coefficient matrix not positive definite"; if (DEBUG_COV) { printlog("non-positive definite coefficient matrix %d:\n", i); m_logoutput(a[i]); } ic = lmc = 0; } if (! posdef) printlog( "non-positive definite coefficient matrix in structure %d", i+1); } for (k = 0; k < v[0]->n_models; k++) m_free(a[k]); efree(a); if (ic) { printlog("Intrinsic Correlation found. Good.\n"); return; } else if (lmc) { printlog("Linear Model of Coregionalization found. Good.\n"); return; } } /* * lmc does not hold: check on Cauchy Swartz */ pr_warning("No Intrinsic Correlation or Linear Model of Coregionalization found\nReason: %s", reason ? reason : "unknown"); if (gl_nocheck == 0) { pr_warning("[add `set = list(nocheck = 1)' to the gstat() or krige() to ignore the following error]\n"); ErrMsg(ER_IMPOSVAL, "variograms do not satisfy a legal model"); } printlog("Now checking for Cauchy-Schwartz inequalities:\n"); for (i = 0; i < n_vars; i++) for (j = 0; j < i; j++) if (is_valid_cs(v[LTI(i,i)], v[LTI(j,j)], v[LTI(i,j)])) { printlog("variogram(%s,%s) passed Cauchy-Schwartz\n", name_identifier(j), name_identifier(i)); } else pr_warning("Cauchy-Schwartz inequality found for variogram(%s,%s)", name_identifier(j), name_identifier(i) ); return; }