pykrige.uk.UniversalKriging

class pykrige.uk.UniversalKriging(x, y, z, variogram_model='linear', variogram_parameters=None, variogram_function=None, nlags=6, weight=False, anisotropy_scaling=1.0, anisotropy_angle=0.0, drift_terms=None, point_drift=None, external_drift=None, external_drift_x=None, external_drift_y=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

Provides greater control over 2D kriging by utilizing drift terms.

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

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

  • z (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 (float, optional) – Scalar stretching value to take into account anisotropy. 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, if anisotropy_angle is not 0).

  • anisotropy_angle (float, optional) – CCW angle (in degrees) by which to rotate coordinate system in order to take into account anisotropy. Default is 0 (no rotation). Note that the coordinate system is rotated.

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

  • point_drift (array_like, optional) – Array-like object that contains the coordinates and strengths of the point-logarithmic drift terms. Array shape must be (N, 3), where N is the number of point drift terms. First column (index 0) must contain x-coordinates, second column (index 1) must contain y-coordinates, and third column (index 2) must contain the strengths of each point term. Strengths are relative, so only the relation of the values to each other matters. Note that the code will appropriately deal with point-logarithmic terms that are at the same coordinates as an evaluation point or data point, but Python will still kick out a warning message that an ln(0) has been encountered. If the problem involves anisotropy, the well coordinates will be adjusted and the drift values will be calculated in the adjusted data frame.

  • external_drift (array_like, optional) – Gridded data used for the external Z scalar drift term. Must be shape (M, N), where M is in the y-direction and N is in the x-direction. Grid spacing does not need to be constant. If grid spacing is not constant, must specify the grid cell sizes. If the problem involves anisotropy, the external drift values are extracted based on the pre-adjusted coordinates (i.e., the original coordinate system).

  • external_drift_x (array_like, optional) – X-coordinates for gridded external Z-scalar data. Must be shape (M,) or (M, 1), where M is the number of grid cells in the x-direction. The coordinate is treated as the center of the cell.

  • external_drift_y (array_like, optional) – Y-coordinates for gridded external Z-scalar data. Must be shape (N,) or (N, 1), where N is the number of grid cells in the y-direction. The coordinate is treated as the center of the cell.

  • 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 two spatial coordinates and must return a single value for each coordinate pair. It must be set up to be called with only two arguments, first an array of x values and second an array of y values. If the problem involves anisotropy, the drift values are calculated in the adjusted data frame.

  • verbose (bool, 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 variogram model with the actual binned data.

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

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).

get_variogram_points()

Returns both the lags and the variogram function evaluated at each of them.

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()

Allows user to switch plot display on/off.

switch_verbose()

Allows user to switch code talk-back on/off.

update_variogram_model(variogram_model[, ...])

Allows user to update variogram type and/or variogram model parameters.

display_variogram_model()[source]

Displays variogram model with the actual binned data.

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

Calculates a kriged grid and the associated variance. Includes drift terms.

Parameters:
  • style (str) – Specifies how to treat input kriging points. Specifying ‘grid’ treats xpoints and ypoints as two arrays of x and y coordinates that define a rectangular grid. Specifying ‘points’ treats xpoints and ypoints as two arrays that provide coordinate pairs at which to solve the kriging system. Specifying ‘masked’ treats xpoints and ypoints as two arrays of x and y 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 MxN 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 MxN grid. If style is specified as ‘points’, y-coordinates of specific points at which to solve kriging system. Note that in this case, xpoints and ypoints must have the same dimensions (i.e., M = N).

  • mask (boolean array, shape (M, N), optional) – Specifies the points in the rectangular grid defined by xpoints and ypoints 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 should will not be solved at the point.

  • backend (str, 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’. Note that Cython backend is not supported for UK.

  • 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 UniversalKriging 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 (M, N), where M is the number of y grid-points and N is the number of x grid-points, or shape (M, ) or (N, 1), where M is the number of points at which to evaluate the kriging system.

Returns:

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

  • sigmasq (ndarray, shape (M, 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.

get_statistics()[source]

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

get_variogram_points()[source]

Returns both the lags and the variogram function evaluated at each of them.

The evaluation of the variogram function and the lags are produced internally. This method is convenient when the user wants to access to the lags and the resulting variogram (according to the model provided) for further analysis.

Returns:

lags (array) - the lags at which the variogram was evaluated variogram (array) - the variogram function evaluated at the lags

Return type:

(tuple) tuple containing

plot_epsilon_residuals()[source]

Plots the epsilon residuals for the variogram fit.

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]

Allows user to switch plot display on/off. Takes no arguments.

switch_verbose()[source]

Allows user to switch code talk-back on/off. Takes no arguments.

update_variogram_model(variogram_model, variogram_parameters=None, variogram_function=None, nlags=6, weight=False, anisotropy_scaling=1.0, anisotropy_angle=0.0)[source]

Allows user to update variogram type and/or variogram model parameters.

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. Defualt is 6.

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

  • anisotropy_scaling (float, optional) – Scalar stretching value to take into account anisotropy. Default is 1 (effectively no stretching). Scaling is applied in the y-direction.

  • anisotropy_angle (float, optional) – CCW angle (in degrees) by which to rotate coordinate system in order to take into account anisotropy. Default is 0 (no rotation).

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>}