pykrige.uk3d.UniversalKriging3D

class pykrige.uk3d.UniversalKriging3D(x, y, z, val, variogram_model='linear', variogram_parameters=None, variogram_function=None, nlags=6, weight=False, anisotropy_scaling_y=1.0, anisotropy_scaling_z=1.0, anisotropy_angle_x=0.0, anisotropy_angle_y=0.0, anisotropy_angle_z=0.0, drift_terms=None, specified_drift=None, functional_drift=None, verbose=False, enable_plotting=False, exact_values=True, pseudo_inv=False, pseudo_inv_type='pinv')[source]

Bases: object

Three-dimensional universal kriging.

Parameters:
  • x (array_like) – X-coordinates of data points.

  • y (array_like) – Y-coordinates of data points.

  • z (array_like) – Z-coordinates of data points.

  • val (array_like) – Values at data points.

  • variogram_model (str or GSTools CovModel, optional) – Specified which variogram model to use; may be one of the following: linear, power, gaussian, spherical, exponential, hole-effect. Default is linear variogram model. To utilize a custom variogram model, specify ‘custom’; you must also provide variogram_parameters and variogram_function. Note that the hole-effect model is only technically correct for one-dimensional problems. You can also use a GSTools CovModel.

  • variogram_parameters (list or dict, optional) –

    Parameters that define the specified variogram model. If not provided, parameters will be automatically calculated using a “soft” L1 norm minimization scheme. For variogram model parameters provided in a dict, the required dict keys vary according to the specified variogram model:

    # linear
        {'slope': slope, 'nugget': nugget}
    # power
        {'scale': scale, 'exponent': exponent, 'nugget': nugget}
    # gaussian, spherical, exponential and hole-effect:
        {'sill': s, 'range': r, 'nugget': n}
        # OR
        {'psill': p, 'range': r, 'nugget': n}
    

    Note that either the full sill or the partial sill (psill = sill - nugget) can be specified in the dict. For variogram model parameters provided in a list, the entries must be as follows:

    # linear
        [slope, nugget]
    # power
        [scale, exponent, nugget]
    # gaussian, spherical, exponential and hole-effect:
        [sill, range, nugget]
    

    Note that the full sill (NOT the partial sill) must be specified in the list format. For a custom variogram model, the parameters are required, as custom variogram models will not automatically be fit to the data. Furthermore, the parameters must be specified in list format, in the order in which they are used in the callable function (see variogram_function for more information). The code does not check that the provided list contains the appropriate number of parameters for the custom variogram model, so an incorrect parameter list in such a case will probably trigger an esoteric exception someplace deep in the code. NOTE that, while the list format expects the full sill, the code itself works internally with the partial sill.

  • variogram_function (callable, optional) – A callable function that must be provided if variogram_model is specified as ‘custom’. The function must take only two arguments: first, a list of parameters for the variogram model; second, the distances at which to calculate the variogram model. The list provided in variogram_parameters will be passed to the function as the first argument.

  • nlags (int, optional) – Number of averaging bins for the semivariogram. Default is 6.

  • weight (bool, optional) – Flag that specifies if semivariance at smaller lags should be weighted more heavily when automatically calculating variogram model. The routine is currently hard-coded such that the weights are calculated from a logistic function, so weights at small lags are ~1 and weights at the longest lags are ~0; the center of the logistic weighting is hard-coded to be at 70% of the distance from the shortest lag to the largest lag. Setting this parameter to True indicates that weights will be applied. Default is False. (Kitanidis suggests that the values at smaller lags are more important in fitting a variogram model, so the option is provided to enable such weighting.)

  • anisotropy_scaling_y (float, optional) – Scalar stretching value to take into account anisotropy in the y direction. Default is 1 (effectively no stretching). Scaling is applied in the y direction in the rotated data frame (i.e., after adjusting for the anisotropy_angle_x/y/z, if anisotropy_angle_x/y/z is/are not 0).

  • anisotropy_scaling_z (float, optional) – Scalar stretching value to take into account anisotropy in the z direction. Default is 1 (effectively no stretching). Scaling is applied in the z direction in the rotated data frame (i.e., after adjusting for the anisotropy_angle_x/y/z, if anisotropy_angle_x/y/z is/are not 0).

  • anisotropy_angle_x (float, optional) – CCW angle (in degrees) by which to rotate coordinate system about the x axis in order to take into account anisotropy. Default is 0 (no rotation). Note that the coordinate system is rotated. X rotation is applied first, then y rotation, then z rotation. Scaling is applied after rotation.

  • anisotropy_angle_y (float, optional) – CCW angle (in degrees) by which to rotate coordinate system about the y axis in order to take into account anisotropy. Default is 0 (no rotation). Note that the coordinate system is rotated. X rotation is applied first, then y rotation, then z rotation. Scaling is applied after rotation.

  • anisotropy_angle_z (float, optional) – CCW angle (in degrees) by which to rotate coordinate system about the z axis in order to take into account anisotropy. Default is 0 (no rotation). Note that the coordinate system is rotated. X rotation is applied first, then y rotation, then z rotation. Scaling is applied after rotation.

  • drift_terms (list of strings, optional) – List of drift terms to include in three-dimensional universal kriging. Supported drift terms are currently ‘regional_linear’, ‘specified’, and ‘functional’.

  • specified_drift (list of array-like objects, optional) – List of arrays that contain the drift values at data points. The arrays must be shape (N,) or (N, 1), where N is the number of data points. Any number of specified-drift terms may be used.

  • functional_drift (list of callable objects, optional) – List of callable functions that will be used to evaluate drift terms. The function must be a function of only the three spatial coordinates and must return a single value for each coordinate triplet. It must be set up to be called with only three arguments, first an array of x values, the second an array of y values, and the third an array of z values. If the problem involves anisotropy, the drift values are calculated in the adjusted data frame.

  • verbose (boolean, optional) – Enables program text output to monitor kriging process. Default is False (off).

  • enable_plotting (boolean, optional) – Enables plotting to display variogram. Default is False (off).

  • exact_values (bool, optional) – If True, interpolation provides input values at input locations. If False, interpolation accounts for variance/nugget within input values at input locations and does not behave as an exact-interpolator [2]. Note that this only has an effect if there is variance/nugget present within the input data since it is interpreted as measurement error. If the nugget is zero, the kriged field will behave as an exact interpolator.

  • pseudo_inv (bool, optional) – Whether the kriging system is solved with the pseudo inverted kriging matrix. If True, this leads to more numerical stability and redundant points are averaged. But it can take more time. Default: False

  • pseudo_inv_type (str, optional) –

    Here you can select the algorithm to compute the pseudo-inverse matrix:

    • ”pinv”: use pinv from scipy which uses lstsq

    • ”pinvh”: use pinvh from scipy which uses eigen-values

    Default: “pinv”

References

[1]

P.K. Kitanidis, Introduction to Geostatistcs: Applications in Hydrogeology, (Cambridge University Press, 1997) 272 p.

[2]

N. Cressie, Statistics for spatial data, (Wiley Series in Probability and Statistics, 1993) 137 p.

Methods

display_variogram_model()

Displays semivariogram and variogram model.

execute(style, xpoints, ypoints, zpoints[, ...])

Calculates a kriged grid and the associated variance.

get_epsilon_residuals()

Returns the epsilon residuals for the variogram fit.

get_statistics()

Returns the Q1, Q2, and cR statistics for the variogram fit (in that order).

plot_epsilon_residuals()

Plots the epsilon residuals for the variogram fit.

print_statistics()

Prints out the Q1, Q2, and cR statistics for the variogram fit.

switch_plotting()

Enables/disable variogram plot display.

switch_verbose()

Enables/disables program text output.

update_variogram_model(variogram_model[, ...])

Changes the variogram model and variogram parameters for the kriging system.

display_variogram_model()[source]

Displays semivariogram and variogram model.

execute(style, xpoints, ypoints, zpoints, mask=None, backend='vectorized', specified_drift_arrays=None)[source]

Calculates a kriged grid and the associated variance.

This is now the method that performs the main kriging calculation. Note that currently measurements (i.e., z values) are considered ‘exact’. This means that, when a specified coordinate for interpolation is exactly the same as one of the data points, the variogram evaluated at the point is forced to be zero. Also, the diagonal of the kriging matrix is also always forced to be zero. In forcing the variogram evaluated at data points to be zero, we are effectively saying that there is no variance at that point (no uncertainty, so the value is ‘exact’).

In the future, the code may include an extra ‘exact_values’ boolean flag that can be adjusted to specify whether to treat the measurements as ‘exact’. Setting the flag to false would indicate that the variogram should not be forced to be zero at zero distance (i.e., when evaluated at data points). Instead, the uncertainty in the point will be equal to the nugget. This would mean that the diagonal of the kriging matrix would be set to the nugget instead of to zero.

Parameters:
  • style (str) – Specifies how to treat input kriging points. Specifying ‘grid’ treats xpoints, ypoints, and zpoints as arrays of x, y, and z coordinates that define a rectangular grid. Specifying ‘points’ treats xpoints, ypoints, and zpoints as arrays that provide coordinates at which to solve the kriging system. Specifying ‘masked’ treats xpoints, ypoints, and zpoints as arrays of x, y, and z coordinates that define a rectangular grid and uses mask to only evaluate specific points in the grid.

  • xpoints (array_like, shape (N,) or (N, 1)) – If style is specific as ‘grid’ or ‘masked’, x-coordinates of LxMxN grid. If style is specified as ‘points’, x-coordinates of specific points at which to solve kriging system.

  • ypoints (array_like, shape (M,) or (M, 1)) – If style is specified as ‘grid’ or ‘masked’, y-coordinates of LxMxN grid. If style is specified as ‘points’, y-coordinates of specific points at which to solve kriging system. Note that in this case, xpoints, ypoints, and zpoints must have the same dimensions (i.e., L = M = N).

  • zpoints (array_like, shape (L,) or (L, 1)) – If style is specified as ‘grid’ or ‘masked’, z-coordinates of LxMxN grid. If style is specified as ‘points’, z-coordinates of specific points at which to solve kriging system. Note that in this case, xpoints, ypoints, and zpoints must have the same dimensions (i.e., L = M = N).

  • mask (boolean array, shape (L, M, N), optional) – Specifies the points in the rectangular grid defined by xpoints, ypoints, zpoints that are to be excluded in the kriging calculations. Must be provided if style is specified as ‘masked’. False indicates that the point should not be masked, so the kriging system will be solved at the point. True indicates that the point should be masked, so the kriging system will not be solved at the point.

  • backend (string, optional) – Specifies which approach to use in kriging. Specifying ‘vectorized’ will solve the entire kriging problem at once in a vectorized operation. This approach is faster but also can consume a significant amount of memory for large grids and/or large datasets. Specifying ‘loop’ will loop through each point at which the kriging system is to be solved. This approach is slower but also less memory-intensive. Default is ‘vectorized’.

  • specified_drift_arrays (list of array-like objects, optional) – Specifies the drift values at the points at which the kriging system is to be evaluated. Required if ‘specified’ drift provided in the list of drift terms when instantiating the UniversalKriging3D class. Must be a list of arrays in the same order as the list provided when instantiating the kriging object. Array(s) must be the same dimension as the specified grid or have the same number of points as the specified points; i.e., the arrays either must be shape (L, M, N), where L is the number of z grid-points, M is the number of y grid-points, and N is the number of x grid-points, or shape (N,) or (N, 1), where N is the number of points at which to evaluate the kriging system.

Returns:

  • kvalues (ndarray, shape (L, M, N) or (N,) or (N, 1)) – Interpolated values of specified grid or at the specified set of points. If style was specified as ‘masked’, kvalues will be a numpy masked array.

  • sigmasq (ndarray, shape (L, M, N) or (N,) or (N, 1)) – Variance at specified grid points or at the specified set of points. If style was specified as ‘masked’, sigmasq will be a numpy masked array.

get_epsilon_residuals()[source]

Returns the epsilon residuals for the variogram fit. No arguments.

get_statistics()[source]

Returns the Q1, Q2, and cR statistics for the variogram fit (in that order). No arguments.

plot_epsilon_residuals()[source]

Plots the epsilon residuals for the variogram fit. No arguments.

print_statistics()[source]

Prints out the Q1, Q2, and cR statistics for the variogram fit. NOTE that ideally Q1 is close to zero, Q2 is close to 1, and cR is as small as possible.

switch_plotting()[source]

Enables/disable variogram plot display. No arguments.

switch_verbose()[source]

Enables/disables program text output. No arguments.

update_variogram_model(variogram_model, variogram_parameters=None, variogram_function=None, nlags=6, weight=False, anisotropy_scaling_y=1.0, anisotropy_scaling_z=1.0, anisotropy_angle_x=0.0, anisotropy_angle_y=0.0, anisotropy_angle_z=0.0)[source]

Changes the variogram model and variogram parameters for the kriging system.

Parameters:
  • variogram_model (str or GSTools CovModel) –

    May be any of the variogram models listed above. May also be ‘custom’, in which case variogram_parameters and variogram_function must be specified. You can also use a GSTools CovModel.

  • variogram_parameters (list or dict, optional) – List or dict of variogram model parameters, as explained above. If not provided, a best fit model will be calculated as described above.

  • variogram_function (callable, optional) – A callable function that must be provided if variogram_model is specified as ‘custom’. See above for more information.

  • nlags (int, optional)) – Number of averaging bins for the semivariogram. Default is 6.

  • weight (boolean, optional) – Flag that specifies if semivariance at smaller lags should be weighted more heavily when automatically calculating variogram model. See above for more information. True indicates that weights will be applied. Default is False.

  • anisotropy_scaling_y (float, optional) – Scalar stretching value to take into account anisotropy in y-direction. Default is 1 (effectively no stretching). See above for more information.

  • anisotropy_scaling_z (float, optional) – Scalar stretching value to take into account anisotropy in z-direction. Default is 1 (effectively no stretching). See above for more information.

  • anisotropy_angle_x (float, optional) – Angle (in degrees) by which to rotate coordinate system about the x axis in order to take into account anisotropy. Default is 0 (no rotation). See above for more information.

  • anisotropy_angle_y (float, optional) – Angle (in degrees) by which to rotate coordinate system about the y axis in order to take into account anisotropy. Default is 0 (no rotation). See above for more information.

  • anisotropy_angle_z (float, optional) – Angle (in degrees) by which to rotate coordinate system about the z axis in order to take into account anisotropy. Default is 0 (no rotation). See above for more information.

UNBIAS = True
eps = 1e-10
variogram_dict = {'exponential': <function exponential_variogram_model>, 'gaussian': <function gaussian_variogram_model>, 'hole-effect': <function hole_effect_variogram_model>, 'linear': <function linear_variogram_model>, 'power': <function power_variogram_model>, 'spherical': <function spherical_variogram_model>}