demcompare.dem_tools
This module contains main functions to manipulate DEM raster images.
It represents the primary API to manipulate DEM as xarray dataset in demcompare. Dataset and associated internal functions are described in dataset_tools.py
Module Contents
Classes
Enum type definition for sampling_source parameter |
Functions
|
Reads the input DEM path and parameters and generates |
|
Returns a copy of the input dem. |
|
Writes a Dataset in a tiff file. |
|
Applies pixellic offset to the input dataset's |
|
Creates dem from input array and transform. |
|
Reprojects both DEMs to common grid, common bounds and common georef origin. |
|
Metric computation method |
|
Accumulates the classification layers of each dem |
|
Computes DEM's slope |
|
Compute and optionnally save plots from a dem using pyplot img_show. |
|
Verifies that the input configuration and input dem contain the |
|
Return the curvature of the input dem. |
Attributes
- demcompare.dem_tools.load_dem(path: str, nodata: float = None, band: int = 1, geoid_georef: bool = False, geoid_path: str | None = None, zunit: str = 'm', input_roi: bool | dict | Tuple = False, classification_layers: Dict = None) xarray.Dataset [source]
Reads the input DEM path and parameters and generates the DEM object as xr.Dataset to be handled in demcompare functions.
A DEM can be any raster file opened by rasterio.
- Parameters:
path (str) – path to dem (readable by rasterio)
nodata (float or None) – forcing dem no data value (None by default and if set inside metadata)
band (int) – band to be read in DEM. Default: 1
geoid_georef (bool) – is dem’s georef is geoid
geoid_path (str or None) – optional path to local geoid, default is EGM96
zunit (str) – dem z unit
input_roi (bool, dict or Tuple) – False if dem are to be fully loaded, other options are a dict roi or Tuple
classification_layers (Dict or None) – input classification layers
- Returns:
dem xr.DataSet containing : (see dataset_tools for details)
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:
xr.Dataset
- demcompare.dem_tools.copy_dem(dem: xarray.Dataset) xarray.Dataset [source]
Returns a copy of the input dem.
- Parameters:
dem (xr.Dataset) –
input dem to copy, 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
- Return dem_copy:
copy of the input dem, 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
- Return type:
xr.Dataset
- demcompare.dem_tools.save_dem(dataset: xarray.Dataset, filename: str, nodata: float = DEFAULT_NODATA) xarray.Dataset [source]
Writes a Dataset in a tiff file. If new_array is set, new_array is used as data. Returns written dataset.
- Parameters:
dataset (xr.Dataset) –
xarray.DataSet containing the variables :
image : 2D (row, col) xr.DataArray float32
georef_transform: 1D (trans_len) xr.DataArray
classification_layer_masks : 3D (row, col, indicator) xr.DataArray
filename (str) – output filename
nodata (float) – value of nodata to use
- Returns:
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
- Return type:
xr.Dataset
- demcompare.dem_tools.translate_dem(dataset: xarray.Dataset, x_offset: float | int | numpy.ndarray, y_offset: float | int | numpy.ndarray) xarray.Dataset [source]
Applies pixellic offset to the input dataset’s georeference origin by modifying its transform.
- Parameters:
dataset (xr.Dataset) –
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
x_offset (Union[float, int, ndarray]) – x offset
y_offset (Union[float, int, ndarray]) – y offset
- Returns:
translated 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
- Return type:
xr.Dataset
- demcompare.dem_tools.create_dem(data: numpy.ndarray, transform: numpy.ndarray | rasterio.Affine = None, img_crs: rasterio.crs.CRS | None = None, input_img: str | None = None, bounds: rasterio.coords.BoundingBox = None, classification_layer_masks: Dict | xarray.DataArray = None, nodata: float = None, geoid_georef: bool = False, geoid_path: str | None = None, zunit: str = 'm', source_rasterio: Dict[str, rasterio.DatasetReader] = None) xarray.Dataset [source]
Creates dem from input array and transform.
- The demcompare DEM is an xarray 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
- Parameters:
data (np.ndarray) – image data
transform (np.ndarray or rasterio.Affine) – rasterio georeferencing transformation matrix
input_img (str) – image path
bounds (rasterio.coords.BoundingBox or None) – dem bounds
classification_layer_masks (Dict, xr.DataArray or None) – classification layers
nodata (float or None) – nodata value in the image
geoid_georef (bool) – if dem’s georef is geoid
geoid_path (str or None) – optional path to local geoid, default is EGM96
zunit (str) – unit
source_rasterio (Dict[str,rasterio.DatasetReader] or None) – rasterio dataset reader object
- Img_crs:
image CRS georef
- Returns:
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
- Return type:
xr.Dataset
- class demcompare.dem_tools.SamplingSourceParameter[source]
Bases:
enum.Enum
Enum type definition for sampling_source parameter value are sec or ref to choose in reproject_dems
- demcompare.dem_tools.reproject_dems(sec: xarray.Dataset, ref: xarray.Dataset, initial_shift_x: int | float = 0, initial_shift_y: int | float = 0, sampling_source: str = SamplingSourceParameter.SEC.value) Tuple[xarray.Dataset, xarray.Dataset, Tuple[float, float]] [source]
Reprojects both DEMs to common grid, common bounds and common georef origin.
The common grid, bounds, georef origin are defined by the sampling_source parameter. It defines which is the sampling of the output DEMs.
If sampling_source is “ref”: ref is cropped to the common bounds sec is cropped to the common bounds and resampled
If sampling_source is “sec”: ref is cropped to the common bounds and resampled sec is cropped to the common bounds
- Parameters:
sec (xr.Dataset) –
dem to align xr.DataSet containing :
im : 2D (row, col) xarray.DataArray float32
trans: 1D (trans_len) xarray.DataArray
ref (xr.Dataset) –
ref xr.DataSet containing :
im : 2D (row, col) xarray.DataArray float32
trans: 1D (trans_len) xarray.DataArray
initial_shift_x (Union[int, float]) – optional initial shift x
initial_shift_y (Union[int, float]) – optional initial shift y
sampling_source (str) – ‘ref’ or ‘sec’, the sampling value of the output dems, by defaut “sec”
- Returns:
reproj_cropped_sec xr.DataSet, reproj_cropped_ref xr.DataSet, adapting_factor. The xr.Datasets containing :
im : 2D (row, col) xarray.DataArray float32
trans: 1D (trans_len) xarray.DataArray
- Return type:
xr.Dataset, xr.Dataset, Tuple[float, float]
- demcompare.dem_tools.compute_waveform(dem: xarray.Dataset, output_dir: str = None) Tuple[numpy.ndarray, numpy.ndarray] [source]
Metric computation method
- Parameters:
dem – input data to compute the metric
output_dir (str) – optional output directory
- Returns:
the computed row and col waveforms
- Return type:
Tuple[np.ndarray, np.ndarray]
- demcompare.dem_tools.accumulates_class_layers(ref: xarray.Dataset, sec: xarray.Dataset, diff_ref_sec: xarray.Dataset) xarray.Dataset [source]
Accumulates the classification layers of each dem
- Parameters:
ref – ref dem
sec – sec dem
diff_ref_sec – difference (i.e angular, in altitude, …) dem between ref and sec
- Returns:
the difference dem between ref and sec, with the classification layers
- Return type:
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
- demcompare.dem_tools.compute_dem_slope(dataset: xarray.Dataset, degree: bool = False, add_attribute: bool = True, unit_change: bool = True) xarray.Dataset [source]
Computes DEM’s slope Slope is presented here : http://pro.arcgis.com/ fr/pro-app/tool-reference/spatial-analyst/how-aspect-works.htm
TODO: modify code such as to remove double return at the end
- Parameters:
dataset – dataset
degree (bool) – True if is in degree
add_attribute (bool) – if True, add attribute to the input dataset if False, return the slope
unit_change (bool) – if True change units of the slope if False, no units change
- Returns:
slope
- Return type:
np.ndarray
- demcompare.dem_tools.compute_and_save_image_plots(dem: xarray.Dataset, plot_path: str = None, fig_title: str = None, colorbar_title: str = None, cmap: str = 'terrain', vmin_plot: float = None, vmax_plot: float = None)[source]
Compute and optionnally save plots from a dem using pyplot img_show. see https://matplotlib.org/stable/api/_as_gen/matplotlib.pyplot.imshow.html
- Parameters:
dem (str) – dem object to compute and save image plots
plot_path (str) – path to save the plots if present
fig_title (str) – optional plot figure title
title_colorbar (str) – optional dem path to save the original tif file
cmap (str) – registered colormap name used to map scalar data to colors.
vmin_plot (float) – minimum value of the z axis when plotting
vmax_plot (float) – maximum value of the z axis when plotting
- demcompare.dem_tools.verify_fusion_layers(dem: xarray.Dataset, classif_cfg: Dict, support: str)[source]
Verifies that the input configuration and input dem contain the input necessary layers for the fusion classification.
- Parameters:
dem (str) – name to save the plots
classif_cfg (Dict) – classification layers configuration
support (str) – fusion support, ref or sec
- Returns:
None
- demcompare.dem_tools.compute_curvature_filtering(dem: xarray.Dataset, filter_intensity: float = 0.9, replication: bool = True) xarray.Dataset [source]
Return the curvature of the input dem. First, compute the FFT of the input dem: F(y) = FFT(DEM). Then, apply a filter y^filter_intensity with s=0.9: F(y) = F(y)* y^filter_intensity. Finally, apply the inverse FFT: IFFT(F(y)). We keep the real part (imaginary part = digital noise).
- Parameters:
dem (xr.Dataset) –
dem 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
filter_intensity (float) – parameter of the DEM’s FFT filter intensity. Should be close to 1. Default = 0.9.
replication (bool) – if true, the image is replicated by x4 in order to improve resolution. Default = True.
- Returns:
curvature 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
- Return type:
xr.Dataset