Block and Poisson Kriging#
Centroid-based Poisson Kriging#
- pyinterpolate.centroid_poisson_kriging(semivariogram_model: TheoreticalVariogram, point_support: PointSupport, unknown_block_index: str | Hashable, number_of_neighbors: int, neighbors_range: float = None, is_weighted_by_point_support=True, raise_when_negative_prediction=True, raise_when_negative_error=True, allow_lsa=False, negative_prediction_to_zero=False) Dict[source]
Function performs centroid-based Poisson Kriging of blocks (areal) data.
- Parameters:
- semivariogram_modelTheoreticalVariogram
A fitted variogram.
- point_supportPointSupport
The point support of polygons.
- unknown_block_indexHashable
The id of the block with the unknown value.
- number_of_neighborsint
The minimum number of neighbours that can potentially affect the unknown block.
- neighbors_rangefloat, optional
The maximum range for neighbors search. If not provided then it is read from the semivariogram model.
- is_weighted_by_point_supportbool, default = True
Are distances between blocks weighted by the point support?
- raise_when_negative_predictionbool, default=True
Raise error when prediction is negative.
- raise_when_negative_errorbool, default=True
Raise error when prediction error is negative.
- allow_lsabool, 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.- negative_prediction_to_zerobool, default=False
While prediction is negative, then set it to 0.
- Returns:
- resultsDict
Block index, prediction, error:
{"block id", "zhat", "sig"}
- Raises:
- ValueError
Raised when prediction or prediction error are negative.
Examples
>>> import os >>> import geopandas as gpd >>> from pyinterpolate import ( ... centroid_poisson_kriging, ... Blocks, ... ExperimentalVariogram, ... PointSupport, ... TheoreticalVariogram ... ) >>> >>> >>> FILENAME = 'cancer_data.gpkg' >>> LAYER_NAME = 'areas' >>> DS = gpd.read_file(FILENAME, layer=LAYER_NAME) >>> AREA_VALUES = 'rate' >>> AREA_INDEX = 'FIPS' >>> AREA_GEOMETRY = 'geometry' >>> PS_LAYER_NAME = 'points' >>> PS_VALUES = 'POP10' >>> PS_GEOMETRY = 'geometry' >>> PS = gpd.read_file(FILENAME, layer=PS_LAYER_NAME) >>> >>> CANCER_DATA = { ... 'ds': DS, ... 'index_column_name': AREA_INDEX, ... 'value_column_name': AREA_VALUES, ... 'geometry_column_name': AREA_GEOMETRY ... } >>> POINT_SUPPORT_DATA = { ... 'ps': PS, ... 'value_column_name': PS_VALUES, ... 'geometry_column_name': PS_GEOMETRY ... } >>> BLOCKS = Blocks(**CANCER_DATA) >>> indexes = BLOCKS.block_indexes >>> >>> PS = PointSupport( ... points=POINT_SUPPORT_DATA['ps'], ... ps_blocks=BLOCKS, ... points_value_column=POINT_SUPPORT_DATA['value_column_name'], ... points_geometry_column=POINT_SUPPORT_DATA['geometry_column_name'] ... ) >>> >>> EXPERIMENTAL = ExperimentalVariogram( ... ds=BLOCKS.representative_points_array(), ... step_size=40000, ... max_range=300001 ... ) >>> >>> THEO = TheoreticalVariogram() >>> THEO.autofit( ... experimental_variogram=EXPERIMENTAL, ... sill=150 ... ) >>> centroid_pk = centroid_poisson_kriging( ... semivariogram_model=THEO, ... point_support=PS, ... unknown_block_index=34017, ... number_of_neighbors=8 ... ) >>> print(centroid_pk) {'block_id': 34017, 'zhat': 134.13113195517153, 'sig': 11.795520367595483}
Area-to-area Poisson Kriging#
- pyinterpolate.area_to_area_pk(semivariogram_model: TheoreticalVariogram, point_support: PointSupport, unknown_block_index: str | Hashable, number_of_neighbors: int, neighbors_range: float = None, raise_when_negative_prediction=True, raise_when_negative_error=True, negative_prediction_to_zero=False) dict[source]
Function predicts areal value in an unknown location based on the area-to-area Poisson Kriging
- Parameters:
- semivariogram_modelTheoreticalVariogram
A fitted variogram.
- point_supportPointSupport
The point support of polygons.
- unknown_block_indexHashable
The id of the block with the unknown value.
- number_of_neighborsint
The minimum number of neighbours that can potentially affect the unknown block.
- neighbors_rangefloat, optional
The maximum range for neighbors search. If not provided then it is read from the semivariogram model.
- raise_when_negative_predictionbool, default=True
Raise error when prediction is negative.
- raise_when_negative_errorbool, default=True
Raise error when prediction error is negative.
- negative_prediction_to_zerobool, default=False
While prediction is negative, then set it to 0.
- Returns:
- resultsDict
Block index, prediction, error:
{"block_id", "zhat", "sig"}
- Raises:
- ValueError
Raised when prediction or prediction error are negative.
Examples
>>> import os >>> import geopandas as gpd >>> from pyinterpolate import ( ... area_to_area_pk, ... Blocks, ... ExperimentalVariogram, ... PointSupport, ... TheoreticalVariogram ... ) >>> >>> >>> FILENAME = 'cancer_data.gpkg' >>> LAYER_NAME = 'areas' >>> DS = gpd.read_file(FILENAME, layer=LAYER_NAME) >>> AREA_VALUES = 'rate' >>> AREA_INDEX = 'FIPS' >>> AREA_GEOMETRY = 'geometry' >>> PS_LAYER_NAME = 'points' >>> PS_VALUES = 'POP10' >>> PS_GEOMETRY = 'geometry' >>> PS = gpd.read_file(FILENAME, layer=PS_LAYER_NAME) >>> >>> CANCER_DATA = { ... 'ds': DS, ... 'index_column_name': AREA_INDEX, ... 'value_column_name': AREA_VALUES, ... 'geometry_column_name': AREA_GEOMETRY ... } >>> POINT_SUPPORT_DATA = { ... 'ps': PS, ... 'value_column_name': PS_VALUES, ... 'geometry_column_name': PS_GEOMETRY ... } >>> BLOCKS = Blocks(**CANCER_DATA) >>> indexes = BLOCKS.block_indexes >>> >>> PS = PointSupport( ... points=POINT_SUPPORT_DATA['ps'], ... ps_blocks=BLOCKS, ... points_value_column=POINT_SUPPORT_DATA['value_column_name'], ... points_geometry_column=POINT_SUPPORT_DATA['geometry_column_name'] ... ) >>> >>> EXPERIMENTAL = ExperimentalVariogram( ... ds=BLOCKS.representative_points_array(), ... step_size=40000, ... max_range=300001 ... ) >>> >>> THEO = TheoreticalVariogram() >>> THEO.autofit( ... experimental_variogram=EXPERIMENTAL, ... sill=150 ... ) >>> ata_pk = area_to_area_pk( ... semivariogram_model=THEO, ... point_support=PS, ... unknown_block_index=34017, ... number_of_neighbors=8 ... ) >>> print(ata_pk) {'block_id': 34017, 'zhat': 133.89984110553695, 'sig': 5.9614557671535335}
Area-to-point Poisson Kriging#
- pyinterpolate.area_to_point_pk(semivariogram_model: TheoreticalVariogram, point_support: PointSupport, unknown_block_index: str | Hashable, number_of_neighbors: int, neighbors_range: float = None, raise_when_negative_prediction=True, raise_when_negative_error=True, err_to_nan=True, negative_prediction_to_zero=False)[source]
Function predicts point-support value in the unknown location based on the area-to-point Poisson Kriging
- Parameters:
- semivariogram_modelTheoreticalVariogram
A fitted variogram.
- point_supportPointSupport
The point support of polygons.
- unknown_block_indexHashable
The id of the block with the unknown value.
- number_of_neighborsint
The minimum number of neighbours that can potentially affect the unknown block.
- neighbors_rangefloat, optional
The maximum range for neighbors search. If not provided then it is read from the semivariogram model.
- raise_when_negative_predictionbool, default=True
Raise error when prediction is negative.
- raise_when_negative_errorbool, default=True
Raise error when prediction error is negative.
- err_to_nanbool, default=True
When point interpolation returns
ValueErrorthen set prediction or variance error toNaN.- negative_prediction_to_zerobool, default=False
While prediction is negative, then set it to 0.
- Returns:
- resultsDict[numpy array]
`{"block_id": [["x", "y", "zhat", "sig"], ...]}
- Raises:
- ValueError
Prediction or prediction error are negative.
Examples
>>> import os >>> import geopandas as gpd >>> from pyinterpolate import ( ... area_to_point_pk, ... Blocks, ... ExperimentalVariogram, ... PointSupport, ... TheoreticalVariogram ... ) >>> >>> >>> FILENAME = 'cancer_data.gpkg' >>> LAYER_NAME = 'areas' >>> DS = gpd.read_file(FILENAME, layer=LAYER_NAME) >>> AREA_VALUES = 'rate' >>> AREA_INDEX = 'FIPS' >>> AREA_GEOMETRY = 'geometry' >>> PS_LAYER_NAME = 'points' >>> PS_VALUES = 'POP10' >>> PS_GEOMETRY = 'geometry' >>> PS = gpd.read_file(FILENAME, layer=PS_LAYER_NAME) >>> >>> CANCER_DATA = { ... 'ds': DS, ... 'index_column_name': AREA_INDEX, ... 'value_column_name': AREA_VALUES, ... 'geometry_column_name': AREA_GEOMETRY ... } >>> POINT_SUPPORT_DATA = { ... 'ps': PS, ... 'value_column_name': PS_VALUES, ... 'geometry_column_name': PS_GEOMETRY ... } >>> BLOCKS = Blocks(**CANCER_DATA) >>> indexes = BLOCKS.block_indexes >>> >>> PS = PointSupport( ... points=POINT_SUPPORT_DATA['ps'], ... ps_blocks=BLOCKS, ... points_value_column=POINT_SUPPORT_DATA['value_column_name'], ... points_geometry_column=POINT_SUPPORT_DATA['geometry_column_name'] ... ) >>> >>> EXPERIMENTAL = ExperimentalVariogram( ... ds=BLOCKS.representative_points_array(), ... step_size=40000, ... max_range=300001 ... ) >>> >>> THEO = TheoreticalVariogram() >>> THEO.autofit( ... experimental_variogram=EXPERIMENTAL, ... sill=150 ... ) >>> atp_pk = area_to_point_pk( ... semivariogram_model=THEO, ... point_support=PS, ... unknown_block_index=34017, ... number_of_neighbors=8 ... ) >>> print(atp_pk) {34017: array([[1.82287362e+06, 4.06124500e+05, 1.34010668e+02, 5.47866897e+00]])}