# # Author: Heresh Fattahi # Copyright 2017 # # Modified by V. Brancato (10.12.2019) # Including offset filtering with no SNR masking # import isce import isceobj from osgeo import gdal from scipy import ndimage import numpy as np import os from astropy.convolution import convolve def mask_filterNoSNR(denseOffsetFile,filterSize,outName): # Masking the offsets with a data-based approach # Open the offsets ds = gdal.Open(denseOffsetFile+'.vrt',gdal.GA_ReadOnly) off_az = ds.GetRasterBand(1).ReadAsArray() off_rg = ds.GetRasterBand(2).ReadAsArray() ds = None # Remove values reported as missing data (no value data from ampcor) off_rg[np.where(off_rg < -9999)]=0 off_az[np.where(off_az < -9999)]=0 # Store the offsets in a complex variable off = off_rg + 1j*off_az # Mask the offset based on MAD mask = off_masking(off,filterSize,thre=3) xoff_masked = np.ma.array(off.real,mask=mask) yoff_masked = np.ma.array(off.imag,mask=mask) # Delete not used variables mask = None off = None # Remove residual noisy spots with a median filter on the range offmap xoff_masked.mask = xoff_masked.mask | \ (ndimage.median_filter(xoff_masked.filled(fill_value=0),3) == 0) | \ (ndimage.median_filter(yoff_masked.filled(fill_value=0),3) == 0) # Fill the range offset map iteratively with smoothed values data = xoff_masked.data data[xoff_masked.mask]=np.nan off_rg_filled = fill_with_smoothed(data,filterSize) # Apply the median filter on the offset off_rg_filled = ndimage.median_filter(off_rg_filled,filterSize) # Save the filtered offsets length, width = off_rg_filled.shape # writing the masked and filtered offsets to a file print ('writing masked and filtered offsets to: ', outName) ##Write array to offsetfile off_rg_filled.tofile(outName) # write the xml file img = isceobj.createImage() img.setFilename(outName) img.setWidth(width) img.setAccessMode('READ') img.bands = 1 img.dataType = 'FLOAT' img.scheme = 'BIP' img.renderHdr() return def off_masking(off,filterSize,thre=2): vram = ndimage.median_filter(off.real, filterSize) vazm = ndimage.median_filter(off.imag, filterSize) mask = (np.abs(off.real-vram) > thre) | (np.abs(off.imag-vazm) > thre) | (off.imag == 0) | (off.real == 0) return mask def fill(data, invalid=None): """ Replace the value of invalid 'data' cells (indicated by 'invalid') by the value of the nearest valid data cell Input: data: numpy array of any dimension invalid: a binary array of same shape as 'data'. data value are replaced where invalid is True If None (default), use: invalid = np.isnan(data) Output: Return a filled array. """ if invalid is None: invalid = np.isnan(data) ind = ndimage.distance_transform_edt(invalid, return_distances=False, return_indices=True) return data[tuple(ind)] def fill_with_smoothed(off,filterSize): off_2filt=np.copy(off) kernel = np.ones((filterSize,filterSize),np.float32)/(filterSize*filterSize) loop = 0 cnt2=1 while (cnt2 !=0 & loop<100): loop += 1 idx2= np.isnan(off_2filt) cnt2 = np.sum(np.count_nonzero(np.isnan(off_2filt))) print(cnt2) if cnt2 != 0: off_filt= convolve(off_2filt,kernel,boundary='extend',nan_treatment='interpolate') off_2filt[idx2]=off_filt[idx2] idx3 = np.where(off_filt == 0) off_2filt[idx3]=np.nan off_filt=None return off_2filt def mask_filter(denseOffsetFile, snrFile, band, snrThreshold, filterSize, outName): #masking and Filtering ##Read in the offset file ds = gdal.Open(denseOffsetFile + '.vrt', gdal.GA_ReadOnly) Offset = ds.GetRasterBand(band).ReadAsArray() ds = None ##Read in the SNR file ds = gdal.Open(snrFile + '.vrt', gdal.GA_ReadOnly) snr = ds.GetRasterBand(1).ReadAsArray() ds = None # Masking the dense offsets based on SNR print ('masking the dense offsets with SNR threshold: ', snrThreshold) Offset[snr