@@ -159,7 +159,6 @@ class AggregatedVariogram:
159159 Weighted distances between all blocks:
160160 ``{block id : [distances to other blocks]}``
161161
162-
163162 Methods
164163 -------
165164 regularize()
@@ -176,6 +175,59 @@ class AggregatedVariogram:
176175 the Presence of Irregular Geographical Units, Mathematical Geology 40(1),
177176 101-128, 2008
178177
178+ Examples
179+ --------
180+ >>> import os
181+ >>> import geopandas as gpd
182+ >>> from pyinterpolate import (
183+ ... AggregatedVariogram,
184+ ... Blocks,
185+ ... PointSupport
186+ ... )
187+ >>>
188+ >>>
189+ >>> FILENAME = 'cancer_data.gpkg'
190+ >>> LAYER_NAME = 'areas'
191+ >>> DS = gpd.read_file(FILENAME, layer=LAYER_NAME)
192+ >>> AREA_VALUES = 'rate'
193+ >>> AREA_INDEX = 'FIPS'
194+ >>> AREA_GEOMETRY = 'geometry'
195+ >>> PS_LAYER_NAME = 'points'
196+ >>> PS_VALUES = 'POP10'
197+ >>> PS_GEOMETRY = 'geometry'
198+ >>> PS = gpd.read_file(FILENAME, layer=PS_LAYER_NAME)
199+ >>>
200+ >>> CANCER_DATA = {
201+ ... 'ds': DS,
202+ ... 'index_column_name': AREA_INDEX,
203+ ... 'value_column_name': AREA_VALUES,
204+ ... 'geometry_column_name': AREA_GEOMETRY
205+ ... }
206+ >>> POINT_SUPPORT_DATA = {
207+ ... 'ps': PS,
208+ ... 'value_column_name': PS_VALUES,
209+ ... 'geometry_column_name': PS_GEOMETRY
210+ ... }
211+ >>> BLOCKS = Blocks(**CANCER_DATA)
212+ >>> indexes = BLOCKS.block_indexes
213+ >>>
214+ >>> PS = PointSupport(
215+ ... points=POINT_SUPPORT_DATA['ps'],
216+ ... ps_blocks=BLOCKS,
217+ ... points_value_column=POINT_SUPPORT_DATA['value_column_name'],
218+ ... points_geometry_column=POINT_SUPPORT_DATA['geometry_column_name']
219+ ... )
220+ >>> STEP_SIZE = 20000
221+ >>> MAX_RANGE = 300000
222+ >>> ag = AggregatedVariogram(
223+ ... blocks=BLOCKS,
224+ ... point_support=PS,
225+ ... step_size=STEP_SIZE,
226+ ... max_range=MAX_RANGE
227+ ... )
228+ >>> reg_variogram = ag.regularize()
229+ >>> print(reg_variogram[0])
230+ [20000. 53.13549729]
179231 """
180232
181233 def __init__ (self ,
@@ -368,7 +420,7 @@ def regularize(self,
368420 Returns
369421 -------
370422 regularized_model : numpy array
371- ``[lag, semivariance, number of point pairs, number of blocks included ]``
423+ ``[lag, semivariance]``
372424
373425 Notes
374426 -----
@@ -673,6 +725,59 @@ def regularize(blocks: Blocks,
673725 [1] Goovaerts P., Kriging and Semivariogram Deconvolution in
674726 the Presence of Irregular Geographical
675727 Units, Mathematical Geology 40(1), 101-128, 2008
728+
729+ Examples
730+ --------
731+ >>> import os
732+ >>> import geopandas as gpd
733+ >>> from pyinterpolate import (
734+ ... regularize,
735+ ... Blocks,
736+ ... PointSupport,
737+ ... )
738+ >>>
739+ >>>
740+ >>> FILENAME = 'cancer_data.gpkg'
741+ >>> LAYER_NAME = 'areas'
742+ >>> DS = gpd.read_file(FILENAME, layer=LAYER_NAME)
743+ >>> AREA_VALUES = 'rate'
744+ >>> AREA_INDEX = 'FIPS'
745+ >>> AREA_GEOMETRY = 'geometry'
746+ >>> PS_LAYER_NAME = 'points'
747+ >>> PS_VALUES = 'POP10'
748+ >>> PS_GEOMETRY = 'geometry'
749+ >>> PS = gpd.read_file(FILENAME, layer=PS_LAYER_NAME)
750+ >>>
751+ >>> CANCER_DATA = {
752+ ... 'ds': DS,
753+ ... 'index_column_name': AREA_INDEX,
754+ ... 'value_column_name': AREA_VALUES,
755+ ... 'geometry_column_name': AREA_GEOMETRY
756+ ... }
757+ >>> POINT_SUPPORT_DATA = {
758+ ... 'ps': PS,
759+ ... 'value_column_name': PS_VALUES,
760+ ... 'geometry_column_name': PS_GEOMETRY
761+ ... }
762+ >>> BLOCKS = Blocks(**CANCER_DATA)
763+ >>> indexes = BLOCKS.block_indexes
764+ >>>
765+ >>> PS = PointSupport(
766+ ... points=POINT_SUPPORT_DATA['ps'],
767+ ... ps_blocks=BLOCKS,
768+ ... points_value_column=POINT_SUPPORT_DATA['value_column_name'],
769+ ... points_geometry_column=POINT_SUPPORT_DATA['geometry_column_name']
770+ ... )
771+ >>> STEP_SIZE = 20000
772+ >>> MAX_RANGE = 300000
773+ >>> reg_variogram = regularize(
774+ ... blocks=BLOCKS,
775+ ... point_support=PS,
776+ ... step_size=STEP_SIZE,
777+ ... max_range=MAX_RANGE
778+ ... )
779+ >>> print(reg_variogram[0])
780+ [20000. 53.13549729]
676781 """
677782
678783 agg_var = AggregatedVariogram (
0 commit comments