demcompare.coregistration.nuth_kaab_internal
This module contains functions associated to the Nuth and Kaab universal co-registration (Correcting elevation data for glacier change detection 2011).
Based on the work of geoutils project https://github.com/GeoUtils/geoutils/blob/master/geoutils/dem_coregistration.py Authors : Amaury Dehecq, Andrew Tedstone Date : June 2015 License : MIT
Module Contents
Classes
NuthKaab class, allows to perform a Nuth & Kaab coregistration |
- class demcompare.coregistration.nuth_kaab_internal.NuthKaabInternal(cfg: demcompare.internal_typing.ConfigType = None)[source]
Bases:
demcompare.coregistration.coregistration_template.CoregistrationTemplate
NuthKaab class, allows to perform a Nuth & Kaab coregistration from authors above and adapted in demcompare
- fill_conf_and_schema(cfg: demcompare.internal_typing.ConfigType = None) demcompare.internal_typing.ConfigType [source]
Add default values to the dictionary if there are missing elements and define the configuration schema
- Parameters:
cfg (ConfigType) – coregistration configuration
- Return cfg:
coregistration configuration updated
- Return type:
ConfigType
- _coregister_dems_algorithm(sec: xarray.Dataset, ref: xarray.Dataset) Tuple[demcompare.transformation.Transformation, xarray.Dataset, xarray.Dataset] [source]
Coregister_dems_algorithm, computes coregistration transformation and reprojected coregistered DEMs with Nuth et kaab algorithm Plots might be saved if save_optional_outputs is set.
- Parameters:
sec (xarray Dataset) –
sec xr.DataSet containing :
image : 2D (row, col) xr.DataArray float32
georef_transform: 1D (trans_len) xr.DataArray
classification_layer_masks : 3D (row, col, indicator) xr.DataArray
ref (xarray Dataset) –
ref xr.DataSet containing :
image : 2D (row, col) xr.DataArray float32
georef_transform: 1D (trans_len) xr.DataArray
classification_layer_masks : 3D (row, col, indicator) xr.DataArray
- Returns:
- transformation, reproj_coreg_sec xr.DataSet,
reproj_coreg_ref xr.DataSet. The xr.Datasets containing :
image : 2D (row, col) xr.DataArray float32
georef_transform: 1D (trans_len) xr.DataArray
classification_layer_masks : 3D (row, col, indicator) xr.DataArray
- Return type:
Tuple[Transformation, xr.Dataset, xr.Dataset]
- static interpolate_dem_on_grid(interp_dem: numpy.ndarray, xgrid: numpy.ndarray, ygrid: numpy.ndarray) Tuple[scipy.interpolate.RectBivariateSpline, scipy.interpolate.RectBivariateSpline] [source]
interpolate_dem_on_grid, returns the RectBivariateSpline function to do Bivariate spline approximation over a rectangular mesh.
- Parameters:
interp_dem (np.ndarray) – input dem image
xgrid (np.ndarray) – x axis grid
ygrid (np.ndarray) – x axis grid
- Returns:
spline_1, spline_2,
- Return type:
Tuple[RectBivariateSpline, RectBivariateSpline]
- interpolate_classif_layers(dem_classif: xarray.DataArray, xgrid: numpy.ndarray, ygrid: numpy.ndarray, x_offset: float, y_offset: float)[source]
interpolates the classification layers on the input grids with the input offsets.
- Parameters:
dem_classif (xr.Dataarray) – input dem image
xgrid (np.ndarray) – x axis grid
ygrid (np.ndarray) – x axis grid
x_offset (float) – x offset
y_offset (float) – y offset
- Returns:
interpolated classification layers
- Return type:
xr.Dataarray
- static crop_dem_with_offset(dem: numpy.ndarray, x_offset: float, y_offset: float) numpy.ndarray [source]
Crops the input dem with the given offsets.
- Parameters:
dem (np.ndarray) – input dem image
x_offset (float) – x offset
y_offset (float) – y offset
- Returns:
cropped dem
- Return type:
np.ndarray
- crop_classif_layers(dem_classif: xarray.DataArray, x_offset, y_offset) xarray.DataArray [source]
crop_classif_layers crops and updates the input classification layers with the input offsets.
- Parameters:
dem_classif (xr.Dataarray) – classification layers
x_offset (float) – x offset
y_offset (float) – y offset
- Returns:
cropped classification layers
- Return type:
xr.Dataarray
- static _grad2d(dem: numpy.ndarray) Tuple[numpy.ndarray, numpy.ndarray] [source]
Computes input DEM’s slope and aspect
- Parameters:
dem (np.ndarray) – input dem
- Returns:
slope (fast forward style) and aspect
- Return type:
np.ndarray, np.ndarray
- _nuth_kaab_single_iter(dh: numpy.ndarray, slope: numpy.ndarray, aspect: numpy.ndarray, plot_file: str | bool = None) Tuple[float, float, float] [source]
Computes the horizontal shift between 2 DEMs using the method presented in Nuth & Kaab 2011
- Parameters:
dh (np.ndarray) – elevation difference sec - ref
slope (np.ndarray) – slope for the same locations as the dh
aspect (np.ndarray) – aspect for the same locations as the dh
plot_file (str or bool) – file to where store plot. Set to None if plot is to be printed. Set to False for no plot at all.
- Returns:
east, north, c
- Return type:
float, float, float
- _filter_target(aspect: numpy.ndarray, target: numpy.ndarray) Tuple[numpy.ndarray, numpy.ndarray] [source]
Filter target slice outliers of an input array by a 3*sigma filtering to improve Nuth et kaab fit.
- Parameters:
aspect (np.ndarray) – elevation difference sec - ref
target (np.ndarray) – slope for the same locations as the dh
- Returns:
slice_filt_median, target_filt
- Return type:
Tuple[List]
- _save_fit_plots(plot_file: str, aspect: numpy.ndarray, target: numpy.ndarray, target_filt: numpy.ndarray, slice_filt_median: numpy.ndarray, yfit: numpy.ndarray) None [source]
Compute and save the Nuth et Kaab fit plots of each iteration
- Parameters:
plot_file (str) – file to where store plot
aspect (np.ndarray) – terrain aspect
target (np.ndarray) – nuth et kaab target
target_filt (np.ndarray) – filtered target
slice_filt_median (np.ndarray) – filtered target median
yfit (np.ndarray) – fit polynome of target_filt
- Returns:
None