pykrige.ok.OrdinaryKriging

class pykrige.ok.OrdinaryKriging(x, y, z, variogram_model='linear', variogram_parameters=None, variogram_function=None, nlags=6, weight=False, anisotropy_scaling=1.0, anisotropy_angle=0.0, verbose=False, enable_plotting=False, enable_statistics=False, coordinates_type='euclidean', exact_values=True, pseudo_inv=False, pseudo_inv_type='pinv')[source]

Convenience class for easy access to 2D Ordinary Kriging.

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) – Specifies 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). This parameter has no effect if coordinate_types is set to ‘geographic’.

  • 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. This parameter has no effect if coordinate_types is set to ‘geographic’.

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

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

  • enable_statistics (bool, optional) – Default is False

  • coordinates_type (str, optional) – One of ‘euclidean’ or ‘geographic’. Determines if the x and y coordinates are interpreted as on a plane (‘euclidean’) or as coordinates on a sphere (‘geographic’). In case of geographic coordinates, x is interpreted as longitude and y as latitude coordinates, both given in degree. Longitudes are expected in [0, 360] and latitudes in [-90, 90]. Default is ‘euclidean’.

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

__init__(x, y, z, variogram_model='linear', variogram_parameters=None, variogram_function=None, nlags=6, weight=False, anisotropy_scaling=1.0, anisotropy_angle=0.0, verbose=False, enable_plotting=False, enable_statistics=False, coordinates_type='euclidean', exact_values=True, pseudo_inv=False, pseudo_inv_type='pinv')[source]

Methods

__init__(x, y, z[, variogram_model, ...])

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.

Attributes

eps

variogram_dict