Point Kriging#
Ordinary Kriging#
- pyinterpolate.ordinary_kriging(theoretical_model: TheoreticalVariogram, unknown_locations: ndarray | Point | List | Tuple | GeoSeries | GeometryArray | ArrayLike, known_locations: ArrayLike = None, known_values: ArrayLike = None, known_geometries: ArrayLike = None, neighbors_range=None, no_neighbors=4, max_tick=5.0, use_all_neighbors_in_range=False, allow_approximate_solutions=False, progress_bar: bool = True) ndarray[source]
Function predicts value at unknown location with Ordinary Kriging technique.
- Parameters:
- theoretical_modelTheoreticalVariogram
Fitted theoretical variogram model.
- unknown_locationsUnion[ArrayLike, Point]
Points where you want to estimate value
(x, y) <-> (lon, lat).- known_locationsnumpy array
Known locations:
[x, y, value].- known_valuesArrayLike, optional
Observation in the i-th geometry (from
known_geometries). Optional parameter, if not given thenknown_locationsmust be provided.- known_geometriesArrayLike, optional
Array or similar structure with geometries. It must have the same length as
known_values. Optional parameter, if not given thenknown_locationsmust be provided. Point type geometry.- neighbors_rangefloat, default=None
The maximum distance where we search for neighbors. If
Noneis given then the range is selected from the Theoretical Model’srangattribute.- no_neighborsint, default = 4
The number of n-closest neighbors used for interpolation.
- max_tickfloat, default=5.
If searching for neighbors in a specific direction how big should be a tolerance for increasing the search angle (how many degrees more).
- use_all_neighbors_in_rangebool, default = False
True: if the real number of neighbors within theneighbors_rangeis greater than thenumber_of_neighborsthen take all of them anyway.- allow_approximate_solutionsbool, default=False
Allows the approximation of kriging weights based on the OLS algorithm. We don’t recommend set it to
Trueif you don’t know what are you doing. This parameter can be useful when you have clusters in your dataset, that can lead to singular or near-singular matrix creation.- progress_barbool, default=True
Show a progress bar during the interpolation process.
- Returns:
- : numpy array
[predicted value, variance error, longitude (x), latitude (y)]
- Raises:
- RunetimeError
Singular Matrix in the Kriging system.
Examples
>>> import geopandas as gpd >>> import numpy as np >>> import pandas as pd >>> >>> from pyinterpolate import (build_experimental_variogram, ... build_theoretical_variogram, ordinary_kriging) >>> >>> dem = gpd.read_file('dem.gpkg') >>> unknown_locations = gpd.read_file('unknown_locations.gpkg') >>> step_size = 500 >>> max_range = 10000 >>> exp_variogram = build_experimental_variogram( ... values=dem['dem'], ... geometries=dem['geometry'], ... step_size=step_size, ... max_range=max_range ... ) >>> theo_variogram = build_theoretical_variogram(exp_variogram) >>> interp = ordinary_kriging( ... theoretical_model=theo_variogram, ... unknown_locations=unknown_locations['geometry'], ... known_values=dem['dem'], ... known_geometries=dem['geometry'] ... ) >>> print(interp[0]) [7.91222896e+01 9.72740449e+01 2.38012302e+05 5.51466805e+05]
Simple Kriging#
- pyinterpolate.simple_kriging(theoretical_model: TheoreticalVariogram, unknown_locations: ndarray | Point | List | Tuple | GeoSeries | GeometryArray | ArrayLike, process_mean: float, known_locations: ArrayLike = None, known_values: ArrayLike = None, known_geometries: ArrayLike = None, neighbors_range=None, no_neighbors=4, max_tick=5.0, use_all_neighbors_in_range=False, allow_approximate_solutions=False, progress_bar: bool = True) List[source]
Function predicts value at unknown location using Simple Kriging.
- Parameters:
- theoretical_modelTheoreticalVariogram
Fitted theoretical variogram model.
- unknown_locationsUnion[List, Tuple, numpy array]
Point where you want to estimate value
(x, y) <-> (lon, lat).- process_meanfloat
The mean value of a process over a study area. Should be known before processing. That’s why Simple Kriging has a limited number of applications. You must have multiple samples and well-known area to know this parameter.
- known_locationsnumpy array
Known locations:
[x, y, value].- known_valuesArrayLike, optional
Observation in the i-th geometry (from
known_geometries). Optional parameter, if not given thenknown_locationsmust be provided.- known_geometriesArrayLike, optional
Array or similar structure with geometries. It must have the same length as
known_values. Optional parameter, if not given thenknown_locationsmust be provided. Point type geometry.- neighbors_rangefloat, default=None
The maximum distance where we search for neighbors. If
Noneis given then range is selected from the Theoretical Model.- no_neighborsint, default = 4
The number of the n-closest neighbors used for interpolation.
- max_tickfloat, default=5.
If searching for neighbors in a specific direction how large should be tolerance for increasing the search angle (how many degrees more).
- use_all_neighbors_in_rangebool, default = False
True: if the real number of neighbors within theneighbors_rangeis greater than thenumber_of_neighborsparameter then take all of them anyway.- allow_approximate_solutionsbool, default=False
Allows the approximation of kriging weights based on the OLS algorithm. We don’t recommend set it to
Trueif you don’t know what are you doing. This parameter can be useful when you have clusters in your dataset, that can lead to singular or near-singular matrix creation.- progress_barbool, default=True
Show a progress bar during the interpolation process.
- Returns:
- : numpy array
[predicted value, variance error, longitude (x), latitude (y)]
- Raises:
- RunetimeError
Singular matrix in Kriging system.
Examples
Examples#
>>> import geopandas as gpd >>> import numpy as np >>> import pandas as pd >>> >>> from pyinterpolate import (build_experimental_variogram, ... build_theoretical_variogram, simple_kriging) >>> >>> dem = gpd.read_file('dem.gpkg') >>> unknown_locations = gpd.read_file('unknown_locations.gpkg') >>> step_size = 500 >>> max_range = 10000 >>> exp_variogram = build_experimental_variogram( ... values=dem['dem'], ... geometries=dem['geometry'], ... step_size=step_size, ... max_range=max_range ... ) >>> theo_variogram = build_theoretical_variogram(exp_variogram) >>> process_mean = 20.5 >>> interp = simple_kriging( ... theoretical_model=theo_variogram, ... unknown_locations=unknown_locations['geometry'], ... process_mean=process_mean, ... known_values=dem['dem'], ... known_geometries=dem['geometry'] ... ) >>> print(interp[0]) [7.91222896e+01 9.72740449e+01 2.38012302e+05 5.51466805e+05]
Indicator Kriging#
- class pyinterpolate.IndicatorKriging(indicator_variograms: TheoreticalIndicatorVariogram, unknown_locations: ArrayLike, known_locations: ArrayLike = None, known_values: ArrayLike = None, known_geometries: ArrayLike = None, kriging_type: str = 'ok', process_mean: float = None, neighbors_range=None, no_neighbors=4, max_tick=5.0, use_all_neighbors_in_range=False, allow_approximate_solutions=False, get_expected_values=True)[source]
Class performs indicator kriging.
- Parameters:
- indicator_variogramsIndicatorVariograms
Modeled variograms for each threshold.
- unknown_locationsnumpy ndarray
Points where we want to estimate value
(x, y) <-or-> (lon, lat).- known_locationsnumpy array, optional
The known locations:
[x, y, value].- known_valuesArrayLike, optional
Observation in the i-th geometry (from
known_geometries). Optional parameter, if not given thenknown_locationsmust be provided.- known_geometriesArrayLike, optional
Array or similar structure with geometries. It must have the same length as
known_values. Optional parameter, if not given thenknown_locationsmust be provided. Point type geometry.- kriging_typestr, default = ‘ok’
Type of kriging to perform. Possible values: ‘ok’ - ordinary kriging, ‘sk’ - simple kriging.
- process_meanfloat
The mean value of a process over a study area. Should be known before processing. That’s why Simple Kriging has a limited number of applications. You must have multiple samples and well-known area to know this parameter.
- neighbors_rangefloat, default=None
The maximum distance where we search for neighbors. If
Noneis given then range is selected from thetheoretical_modelrangattribute.- no_neighborsint, default = 4
The number of the n-closest neighbors used for interpolation.
- use_all_neighbors_in_rangebool, default = False
True: if the real number of neighbors within theneighbors_rangeis greater than thenumber_of_neighborsparameter then take all of them anyway.- allow_approximate_solutionsbool, default=False
Allows the approximation of kriging weights based on the OLS algorithm. We don’t recommend set it to
Trueif you don’t know what are you doing. This parameter can be useful when you have clusters in your dataset, that can lead to singular or near-singular matrix creation.- get_expected_valuesbool, default=True
If
Truethen expected values and variances are calculated.
- Attributes:
- thresholdsarray
Thresholds used for indicator kriging.
- coordinatesarray
Coordinates of unknown locations.
- indicator_predictionsarray
Indicator kriging predictions for each threshold and each unknown location.
- expected_valuesarray
Expected values derived from
indicator_predictionsfor each unknown location.- variancesarray
Variances derived from
indicator_predictionsfor each unknown location.
Methods
get_indicator_maps()
Returns dictionary with thresholds and indicator maps for each of them.
get_expected_values()
Returns two arrays: one array with coordinates and expected values, and the second with coordinates and variances.
Examples
>>> import numpy as np >>> from pyinterpolate.kriging.point.indicator import IndicatorKriging >>> from pyinterpolate.semivariogram.indicator.indicator import ( ... ExperimentalIndicatorVariogram, ... TheoreticalIndicatorVariogram ... ) >>> >>> >>> dem = np.random.random(size=(1000, 3)) >>> exp_variogram = ExperimentalIndicatorVariogram( ... values=dem[:, -1], ... geometries=dem[:, :-1], ... number_of_thresholds=3, ... step_size=0.1, ... max_range=0.6 ... ) >>> theo_variograms = TheoreticalIndicatorVariogram( ... experimental_indicator_variogram=exp_variogram ... ) >>> theo_variogram.fit() >>> ikriging = IndicatorKriging( ... known_values=dem[:, -1], ... known_geometries=dem[:, :-1], ... indicator_variograms=theo_variograms, ... unknown_locations=np.random.random(size=(100, 2)), ... kriging_type='ok', ... no_neighbors=16 ... ) >>> print(ikriging.indicator_predictions[0]) [0.12942075 0.13930802 0. ] >>> print(ikriging.expected_values[:3]) [0.64983095 0.64983095 0.64983095] >>> print(ikriging.variances[:3]) [0.28827397 0.28539667 0.30135499] >>> imaps = ikriging.get_indicator_maps() >>> print(imaps.keys()) dict_keys([0.35200453863293624, 0.6783974224309782, 0.9929787628823334])
- get_expected_values_maps() Tuple[ndarray, ndarray][source]
Method returns expected values and variances for each threshold.
- Returns:
- expected_values, variancesnumpy ndarray, numpy ndarray
Expected values and variances.
- get_indicator_maps() Dict[source]
Method returns indicator map for each threshold.
- Returns:
- indicator_mapsDict
Indicator map for each threshold.
Universal Kriging#
- class pyinterpolate.UniversalKriging(known_points: ArrayLike = None, known_values: ArrayLike = None, known_geometries: ArrayLike = None, fitted_regression_model: Any = None)[source]
Performs universal kriging.
- Parameters:
- known_pointsArrayLike, optional
Known points and their values
[lon, lat, value]- known_valuesArrayLike, optional
Observation in the i-th geometry (from
known_geometries). Optional parameter, if not given thenknown_locationsmust be provided.- known_geometriesArrayLike, optional
Array or similar structure with geometries. It must have the same length as
known_values. Optional parameter, if not given thenknown_locationsmust be provided. Point type geometry.- fitted_regression_modeloptional
Any kind of regression model with .predict() method that could be used for trend modeling.
- Attributes:
- known_pointsnumpy array
See
observationsparameter.- trend_modelMultivariateRegression or another regression model
The model responsible for trend prediction. It could be a model from the external package, but it must have .predict() method which takes as an input modeling features (coordinates).
- trend_valuesnumpy array
Observation values predicted by regression model.
- bias_valuesnumpy array
The difference between measured observations and
trend_values.- bias_experimental_modelExperimentalVariogram
The experimental semivariogram of bias.
- bias_modelTheoreticalVariogram
Modeled semivariogram of the data bias.
Methods
fit_trend()
Fits observations into linear regression model.
detrend()
Removes trend from observations (calculates bias).
fit_bias()
Fits bias into semivariogram model.
predict()
Predicts values at unknown locations.
plot_experimental_bias_model()
Plots experimental semivariogram of bias.
plot_theoretical_bias_model()
Plots theoretical semivariogram of bias.
plot_trend_surfaces()
Visual comparison of observations, trend, and bias.
Examples
>>> import geopandas as gpd >>> import numpy as np >>> import pandas as pd >>> >>> from pyinterpolate.kriging.point.universal import UniversalKriging >>> >>> >>> dem = gpd.read_file('dem.gpkg') >>> unknown_points = gpd.read_file('unknown_locations.gpkg') >>> uk = UniversalKriging( ... known_values=dem[:, -1], ... known_geometries=dem[:, :-1] ... ) >>> uk.fit_trend() >>> uk.detrend() >>> uk.fit_bias( ... step_size=500, max_range=10000 ... ) >>> predictions = uk.predict(points=unknown_points) >>> print(predictions[0]) # z_hat, x, y [9.72916006e+01 2.38012302e+05 5.51466805e+05]
- detrend(trend_values: ArrayLike = None)[source]
Function removes trend from observations.
- Parameters:
- trend_valuesoptional, numpy array
Optional array with trend values. It may be an output from the external regression model.
- fit_bias(step_size: float, max_range: float, use_all_models=False)[source]
Function fits bias into variogram models.
- Parameters:
- step_sizefloat
Bins width.
- max_rangefloat
Maximum range of analysis.
- use_all_modelsbool, default = False
Use all available semivariogram models (
True), or only the selection of models - linear, spherical, and power model.
- fit_trend()[source]
Function fits data into MultivariateRegression model.
- plot_experimental_bias_model()[source]
Plots experimental variogram of bias.
- plot_theoretical_bias_model()[source]
Plots theoretical semivariogram of bias.
- plot_trend_surfaces()[source]
Plots initial observations, trend surface, and bias surface.
- predict(points: ArrayLike, neighbors_range: float | None = None, no_neighbors: int = 4, use_all_neighbors_in_range=False, allow_approx_solutions=False, number_of_workers: int = 1, show_progress_bar: bool = True)[source]
- Parameters:
- pointsnumpy array
Coordinates with missing values (to estimate results).
- neighbors_rangefloat, default=None
The maximum distance where we search for neighbors. If
Noneis given then range is selected from thetheoretical_modelrangattribute.- no_neighborsint, default = 4
The number of the n-closest neighbors used for interpolation.
- use_all_neighbors_in_rangebool, default = False
True: if the real number of neighbors within theneighbors_rangeis greater than thenumber_of_neighborsparameter then take all of them anyway.- allow_approx_solutionsbool, default=False
Allows the approximation of kriging weights based on the OLS algorithm. We don’t recommend set it to
Trueif you don’t know what are you doing. This parameter can be useful when you have clusters in your dataset, that can lead to singular or near-singular matrix creation.- number_of_workersint, default=1
How many processing units can be used for predictions. Increase it only for a very large number of interpolated points (~10k+).
- show_progress_barbool, default=True
Show progress bar of predictions.
- Returns:
- : numpy array
Predictions
[predicted value, longitude (x), latitude (y)]