fit_tools

pyddem.fit_tools provides tools to derive filter and interpolate DEM stacks into elevation time series

pyddem.fit_tools.constrain_var_slope_corr(ds, ds_arr, ds_corr, t_vals, uncert, fn_ref_dem=None, nproc=1)[source]

Given a data cube of stacked DEMs, a data cube of stereo-correlations and a per-DEM vector of co-registration RMSE: estimate a median slope and compute a data cube of elevation errors. (error dependance on the slope and stereo-correlation is defined independently) #TODO: allow for a user-defined function of error?

Parameters:
  • ds – xarray Dataset of datacube
  • ds_arr – Data cube of elevations
  • ds_corr – Data cube of stereo-correlations
  • t_vals – Time vector of data cube
  • uncert – Vector of co-registration RMSEs of the stacked DEMs
  • fn_ref_dem – Filename for input reference DEM
  • nproc – Number of cores for multiprocessing [1]
Returns:

Data cube of errors, Raster of median slope values

pyddem.fit_tools.cube_to_stack(ds, out_cube, y0, nice_fit_t, outfile, slope_arr=None, ci=True, clobber=False, filt_bool=False)[source]

Write data cube to netcdf with new time axis, confidence interval, slope, etc… #TODO: this could be simplified quite a bit…

Parameters:
  • ds – xarray Dataset of data cube
  • out_cube – Data cube of time series of elevations
  • y0 – Origin of time vector
  • nice_fit_t – Time vector relative to the origin
  • outfile – Filename of output netcdf
  • slope_arr – Slope raster to write as concomitant DataArray
  • ci – Boolean to look for confidence interval in the out_cube
  • clobber – Boolean to replace output file
  • filt_bool – Boolean to force a boolean datacube as output (e.g., to the mask data cube of the multi-step filtering)
Returns:

pyddem.fit_tools.dask_time_filter_ref(ds, ds_arr, ref_dem, t_vals, ref_date, max_dhdt=[-50, 50], base_thresh=100.0, nproc=1)[source]

Dask wrapper of temporal filtering of stacked DEMs with a reference DEM: removes all elevations outside a positive/negative linear elevation trend from the reference, with a base threshold.

Parameters:
  • ds_arr – Data cube of elevations
  • ds_arr – Data cube of elevations
  • ref_dem – GeoImg of reference DEM
  • t_vals – Time vector of data cube
  • ref_date – Date of reference DEM
  • max_dhdt – Maximum positive/negative elevation change rate per year from reference elevations [-50 m,50 m]
  • base_thresh – Base vertical threshold for the temporal filter (avoid elevations close to reference) [100 m]
  • nproc – Number of cores for multiprocessing [1]
Returns:

Filtered data cube

pyddem.fit_tools.estimate_var(dh, mask_cube, arr_mask=None, cube_mask=None, nsamp=100000)[source]

Estimation of elevation variance from data cube of difference to reference elevations

Parameters:
  • dh – Data cube of difference between stack of DEMs and reference DEM
  • mask_cube – Mask where sampling must occur
  • arr_mask – Add a raster mask (no time dimension)
  • cube_mask – Add a data cube mask
  • nsamp – Number of samples to draw
Returns:

NMAD, number of samples, STD

pyddem.fit_tools.estimate_vgm(ds, ds_arr, inc_mask=None, exc_mask=None, rast_mask=None, nsamp=10000, tstep=0.25, lag_cutoff=None, min_obs=8, nproc=1, pack_size=50)[source]

Aggregate 1D temporal variograms for multiple stack pixels: random sampling within a inclusion mask

Parameters:
  • ds – xarray Dataset of data cube
  • ds_arr – Data cube of elevations
  • inc_mask – Filename of inclusion shapefile
  • exc_mask – Filename of exclusion shapefile
  • rast_mask – Additional mask (boolean raster)
  • nsamp – Number of pixels sampled in time
  • tstep – Time interval for aggregating pairwise time lags
  • lag_cutoff – Maximum time lag to sample
  • min_obs – Minimum of valid observation per pixel for sampling
  • nproc – Number of cores for multiprocessing [1]
Returns:

Lags, Variance, Variance dispersion

pyddem.fit_tools.fit_stack(fn_stack, fit_extent=None, fn_ref_dem=None, ref_dem_date=None, filt_ref='min_max', time_filt_thresh=[-30, 30], inc_mask=None, exc_mask=None, nproc=1, method='gpr', opt_gpr=False, kernel=None, filt_ls=False, conf_filt_ls=0.99, tlim=None, tstep=0.25, outfile='fit.nc', write_filt=False, clobber=False, merge_date=False, dask_parallel=False)[source]

Given a netcdf stack of DEMs, perform filtering and temporal fitting of elevation with uncertainty propagation

Parameters:
  • fn_stack – Filename for input netcdf file
  • fit_extent – extent over which to limit fit, as [xmin, xmax, ymin, ymax]
  • fn_ref_dem – Filename for input reference DEM
  • ref_dem_date – Date of reference DEM
  • filt_ref – Type of filtering. One of ‘min_max’, ‘time’, or ‘both’. Requires a reference DEM.
  • time_filt_thresh – Maximum dh/dt from reference DEM for time filtering
  • inc_mask – Optional inclusion mask. Pixels outside of the mask will not be fit.
  • exc_mask – Optional exclusion mask. Pixels inside of the mask will not be fit.
  • nproc – Number of cores for multiprocessing [1]
  • method – Fitting method, currently supported: Gaussian Process Regression “gpr”, Ordinary Least Squares “ols” and Weighted Least Squares “wls” [“gpr”]
  • opt_gpr – Run learning optimization in the GPR fitting [False]
  • kernel – Kernel
  • filt_ls – Filter least square with a first fit [False]
  • conf_filt_ls – Confidence interval to filter least square fit [99%]
  • tlim – time range [t_min t_max] over which to fit. Default is full range of stack.
  • tstep – Temporal step for fitted stack [0.25 year]
  • outfile – Path to outfile
  • write_filt – Write filtered stack to file
  • clobber – Overwrite existing output files
  • merge_date – Merge any DEMs which have the same acquisition date
  • dask_parallel – run with dask parallel tools
Returns:

pyddem.fit_tools.get_dem_date(ds, ds_filt, t, outname)[source]

Extract a DEM from an elevation time series, with elevation error and time lag to the closest observation

Parameters:
  • ds – xarray Dataset of elevation time series
  • ds_filt – xarray dataset of boolean data cube of valid observation used to generate the time series

;param t: Date of extraction :outname: Filename for output rasters

Returns:
pyddem.fit_tools.get_var_by_corr_slope_bins(ds, ds_arr, arr_slope, bins_slope, cube_corr, bins_corr, outfile, inc_mask=None, exc_mask=None, nproc=1)[source]

Sampling method for elevation measurement error from stacks of DEMs: binning of external variables and sampling of variance

Parameters:
  • ds – xarray Dataset of data cube
  • ds_arr – Data cube of elevations
  • arr_slope – Slope raster
  • bins_slope – Bins to sample categories of slope
  • cube_corr – Data cube of stereo-correlations
  • bins_corr – Bins to sample categories of stereo-correlations
  • outfile – Filename for output files
  • inc_mask – Filename of inclusion shapefile for the sampling of variance
  • exc_mask – Filename of exclusion shapefile for the sampling of variance
  • nproc – Number of cores for multiprocessing [1]
Returns:

pyddem.fit_tools.gpr(data, t_vals, uncert, t_pred, opt=False, kernel=None, not_stable=True, detrend_ls=False, loop_detrend=False)[source]

Gaussian Process regression wrapper

Parameters:
  • data – Data cube of elevations
  • t_vals – Time vector of data cube
  • uncert – Data cube of errors
  • t_pred – Time vector to predict
  • opt – Boolean, run kernel optimization
  • kernel – Provide kernel object
  • not_stable – Mask of unstable terrain (to apply different kernels)
  • detrend_ls – Boolean, remove linear trend before applying GP
  • loop_detrend – Boolean, remove linear trend at each iteration before applying GP
Returns:

Data cube of fitted elevations, Data cube of propagated 1-sigma elevation errors, Data cube mask of valid data

pyddem.fit_tools.isel_maskout(ds, inc_mask)[source]

Mask out values out of shapefile and minimize spatial extent, remove void time steps (to decrease computing time)

Parameters:ds – xarray Dataset of datacube
Returns:ds: reduced Dataset
pyddem.fit_tools.isel_merge_dupl_dates(ds)[source]

Merge duplicate dates of a xarray dataset by taking the mean of similar dates (e.g., overlapping same-date ASTER DEMs)

Parameters:ds – xarray Dataset of datacube
Returns:ds: merged Dataset
pyddem.fit_tools.ls(subarr, t_vals, err, weigh, filt_ls=False, conf_filt=0.99)[source]

Ordinary or Weighted Least-Squares wrapper

Parameters:
  • subarr – Data cube of elevations
  • t_vals – Time vector of data cube
  • err – Data cube of errors
  • weigh – Boolean, use weighted least-squares
  • filt_ls – Boolean, use a first CI filtering before fitting again
  • conf_filt – Confidence interval of the filtering fit
Returns:

Raster concatenation (LS slope, LS intercept, slope 1-sigma error, number of samples, date of first sample, date of last sample), Data cube mask of valid data

pyddem.fit_tools.make_dh_animation(ds, fn_shp=None, month_a_year=None, rates=False, figsize=(10, 10), t0=None, t1=None, dh_max=20, var='z', cmap='RdYlBu', xlbl='easting (km)', ylbl='northing (km)')[source]

Generate a GIF from elevation time series

Parameters:
  • ds – xarray Dataset of elevation time series
  • month_a_year – Numerical month to keep only one month per year in the animation
  • rates – Display rates instead of cumulative change
  • figsize – Tuple of figure size
  • t0 – Starting date
  • t1 – End date
  • dh_max – Max scale for colorbar
  • var – Variable to display in the dataset
  • cmap – Colormap
  • xlbl – xlabel for animation
  • ylbl – ylabel for animation
Returns:

Figure, List of images and annotations for animation

pyddem.fit_tools.manual_refine_sampl_temporal_vgm(fn_stack, fn_ref_dem, out_dir, filt_ref='both', max_dhdt=[-50, 50], ref_dem_date=<MagicMock name='mock()' id='140317251097040'>, inc_mask=None, gla_mask=None, nproc=1)[source]

Full sampling method of temporal variogram from stacks of DEMs: pre-filtering of data, binning of external variables and sampling of temporal variograms

Parameters:
  • fn_stack – Filename of netcdf stack of DEMs
  • fn_ref_dem – Filename of reference DEM
  • out_dir – Output directory
  • filt_ref – Filtering method
  • max_dhdt – Maximum positive/negative elevation change rate per year from reference elevations [-50 m,50 m]
  • ref_dem_date – Date of reference DEM
  • inc_mask – Filename of inclusion shapefile for the stack
  • gla_mask – Filename of inclusion shapefile for the sampling of variogram
  • nproc – Number of cores for multiprocessing [1]
Returns:

pyddem.fit_tools.ols_matrix(X, Y, conf_interv=0.68, conf_slope=0.68)[source]

Ordinary Least-Squares (matrix, for optimized processing time)

Parameters:
  • X – Data cube of x
  • Y – Data cube of y
  • conf_interv – Confidence interval for y output
  • conf_slope – Confidence interval for WLS slope output
Returns:

Slope, Intercept, Slope CI, Y lower CB, Y upper CB

pyddem.fit_tools.prefilter_stack(ds, ds_arr, fn_ref_dem, t_vals, filt_ref='min_max', ref_dem_date=None, max_dhdt=[-50, 50], nproc=1)[source]

Given a data cube of stacked DEMs and a reference DEM: condition a first-order spatial and temporal filtering of outliers.

Parameters:
  • ds – xarray Dataset of datacube
  • ds_arr – Data cube of elevations
  • t_vals – Time vector of data cube
  • filt_ref – Method of filtering: spatial “min_max”, temporal “time” or “both” [“min_max”]
  • fn_ref_dem – Filename for input reference DEM
  • ref_dem_date – Date of reference DEM
  • max_dhdt – Maximum positive/negative elevation change rate per year from reference elevations [-50 m,50 m]
  • nproc – Number of cores for multiprocessing [1]
Returns:

Filtered data cube

pyddem.fit_tools.robust_wls_ref_filter_stack(ds, ds_arr, err_arr, t_vals, fn_ref_dem, ref_dem_date=<MagicMock name='mock()' id='140317250602432'>, max_dhdt=[-50, 50], nproc=1, cutoff_kern_size=1000, max_deltat_ref=2.0, base_thresh=100.0)[source]

Given a data cube of stacked DEMs, of elevation errors and a reference DEM: condition a spatial and temporal filtering of outliers based on WLS linear trends.

Parameters:
  • ds – xarray Dataset of data cube
  • ds_arr – Data cube of elevations
  • err_arr – Data cube of elevation errors
  • t_vals – Time vector of data cube
  • fn_ref_dem – Filename for input reference DEM
  • ref_dem_date – Date of reference DEM
  • max_dhdt – Maximum positive/negative elevation change rate per year from reference elevations [-50 m,50 m]
  • nproc – Number of cores for multiprocessing [1]
  • cutoff_kern_size – Radius size for kernel spatial filtering [1000 m]
  • max_deltat_ref – Maximum time lag between acquisition of reference DEM and date provided [2 yr]
  • base_thresh – Base vertical threshold for the temporal filter (avoid elevations close to reference) [100 m]
Returns:

Filtered data cube

pyddem.fit_tools.spat_filter_ref(ds_arr, ref_dem, cutoff_kern_size=500, cutoff_thr=20.0, nproc=1)[source]

Spatial filtering of stacked DEMs with a reference DEM: removes all elevations smaller than the minimum or larger than the maximum reference elevation found within a circle of certain size, with a base threshold.

Parameters:
  • ds_arr – Data cube of elevations
  • ref_dem – GeoImg of reference DEM
  • cutoff_kern_size – Radius size for kernel spatial filtering [500 m]
  • cutoff_thr – Base vertical threshold for the spatial filter (avoid elevations close to reference)
Returns:

ds_arr: Filtered data cube

pyddem.fit_tools.time_filter_ref(ds_arr, ref_arr, t_vals, ref_date, max_dhdt=[-50, 50], base_thresh=100.0)[source]
Temporal filtering of stacked DEMs with a reference DEM: removes all elevations outside a positive/negative linear
elevation trend from the reference, with a base threshold.
Parameters:
  • ds_arr – Data cube of elevations
  • ref_arr – Raster of reference DEM
  • t_vals – Time vector of data cube
  • ref_date – Date of reference DEM
  • max_dhdt – Maximum positive/negative elevation change rate per year from reference elevations [-50 m,50 m]
  • base_thresh – Base vertical threshold for the temporal filter (avoid elevations close to reference) [100 m]
Returns:

ds_arr: Filtered data cube

pyddem.fit_tools.vgm_1d(argsin)[source]

1D variogram sampling

Parameters:argsin – tuple of x values, y values, lag cutoff value and interval of x sampling (for easier parallel proc)
Returns:semivariance
pyddem.fit_tools.vgm_1d_med(argsin)[source]

1D median variogram sampling

Parameters:argsin – tuple of x values, y values, lag cutoff value and interval of x sampling (for easier parallel proc)
Returns:semivariance
pyddem.fit_tools.wls_matrix(x, y, w, conf_interv=0.68, conf_slope=0.68)[source]

Weighted Least-Squares (matrix, for optimized processing time)

Parameters:
  • x – Data cube of x
  • y – Data cube of y
  • w – Data cube of weights
  • conf_interv – Confidence interval for y output
  • conf_slope – Confidence interval for WLS slope output
Returns:

Slope, Intercept, Slope CI, Y lower CB, Y upper CB