5. Setup and formating for KF

5.1. Prepare inputs

5.1.1. Format

KFTS-InSAR reads an HDF5 file containing the interferogram stack. It can contain the following dataset, although other format are suported, for N interferograms of (Y, X) pixels linking M acquisitions:

Jmat            Dataset {N, M}     Connectivity matrix [-1,1,0]
bperp           Dataset {N}        Perpendicular baseline array (not used by KFTS yet)
dates           Dataset {M}        Array of ordinal values for SAR acquisition dates
figram          Dataset {N, Y, X}  Array of interferograms
tims            Dataset {M}        Array of SAR acquisition times in decimal years

Note

KFTS is able to read MintPy ifgramStack file without the use of the prepare_input.py routine

5.1.2. Prepare input interferogram stack [optional]

Preprocessing only useful if the HDF5 file discribed above has to be constructed. Functions are optimized for ISCE architecture.

class prepare_input.BuildStack(verbose, ylims, xlims)[source]

Quick class to deal with stack of interferogram

ConnectMatrix(dates)[source]

Gets the connectivity matrix for given set of IFGs. Args:

  • dates -> List of pairs of strings

Returns:
  • Uts -> Unique dates of acquisitions

  • connMat -> Connectivity matrix (1,-1,0)

datenum(datelist)[source]

Converts list of dates in yymmdd/yyyymmdd format to ordinal array.

deramp(data, out, network=True, poly=3, dref=[(0, 0)])[source]

Network deramping of the stack of interferograms. Used when no GPS is available.

Args:
  • data Input stack object (interferograms)

  • out Output stack object (dataset will be overwritten)

Kwargs:
  • network Network deramping or individual deramping

  • poly Polynomial code for deramping

  • dref Coordinates of points delimiting a polygon taken

    as the reference region for estimating a ramp (list of 2-tuple)

findramp(phs, mask)[source]

Estimate the best-fitting ramp parameters for a matrix. Args:

  • phs - Input unwrapped phase matrix.

  • mask - Mask of reliable values.

  • poly - Integer flag for type of correction.

    1 - Constant 3 - Plane 4 - Plane + cross term

Returns:
  • ramppoly - Polynomial corresponding to the rank

referenceIgram(phs)[source]

remove for each single 2D interfero phs

removeramp(phs, ramppoly)[source]

Deramp a matrix with the given ramp polynomial. .. Args:

  • phs Input matrix

  • ramppoly Ramp polynomial

class prepare_input.GetConfig(configfile)[source]

Class to run read config file and keep everything in object

prepare_input.getBaselines(baselineDir, masterDate)[source]

Get baseline values from baseline files.

prepare_input.getPairs(igramsDir)[source]

Get the two dates of each interferogram produced by ISCE, in Igram directory.

prepare_input.makefnames(dirname, dates1, dates2)[source]

Generates paths to the files needed for creating the stack. .. Args:

  • dates1 - Date of master

  • dates2 - Date of slave

prepare_input.save2file(dates1, dates2, outfile, dateList=None, baselineList=None, satnameList=None)[source]

Save baselines to a file. It can be used for GIAnT time series analysis.

5.2. Read and formate data

Processing steps contained in the main KFTS-InSAR routine.

class kf.readinput.SetupKF(h5file, fmt='ISCE', comm=False, mpiarg=0, utime='years', verbose=True, subregion=None, cohTh=None, refyx=None)[source]
Class for reading and modifying interferograms for Kalman filtering.
h5file:

.h5 file path containing

  • time : decimal dates relative to first acquisition (usually in years)

  • dates : absolute data

  • igram : interferograms

  • links : connection between phases to build interfero (M x N), 0 1 and -1

  • bperp : perpendicular baseline between aquisitions (not exploited yet)

Opts:
fmt:

format of H5file input, default is ‘ISCE’ (also ‘RAW’ ‘MintPy’)

comm:

do you use parallel features of mpi4py (True or False, default False)

mpiarg:

precise rank and size of communicator (tuple used if mpi=True)

utime:

time unit as string (years or days)

verbose:

print stuff? (True or False)

subregion:

subregion class instance containing x and y bounds as pixel numbers

cohTh:

minimum coherence used to mask pixels in each interferogram

apply_input_mask(maskfile)[source]
  • maskfile: hdf5 file containing a ‘mask’ dataset of (Y,X) shape

copydata2file(fin)[source]

Create new file to copy and edit data (interferograms) in a safe and memory-saving way before splitting tasks between workers

create_R(rr)[source]
Create covariance matrix of data
rr:

variance in observation (=interferograms) noise (float)

dividepxls(mpi, mpiarg)[source]

Check if MPI used and divide pxls into subsets for different workers

filter_by_coherence(fin, cohTh)[source]
Remove points below a given threshold in all interferograms
fin:

open HDF5 file containing a ‘coherence’ dataset

get_interf_pairs()[source]
Extract indices of phases substracted together to build interfero.
  • imoins : indice of phases substracted to iplus (M)

  • iplus : indices of phases added (M)

mintpy2kfts(fin)[source]

Read and translate information in the hDF5 output of MintPy named ifgramStack.h5 Essentially reformat dates and build matrix mapping interferogram to time

fin:

open HDF5 file

pxl_with_nodata(thres=30, chunks=200, plot=False)[source]
Check for pixels with little data and build mask.
xv,yv:

grid element

Opts:
thres:

minimum number of interferogram required to be available for one pixel

chunks:

chunks of columns to be processed simultaneously

Return:
  • pairs of indices in list

rereference(ref)[source]

Remove the interferoùetric value on a given pixel to all interferograms (optional if interferogram stack already referenced)

ref:

(y,x) pixel index of reference

select_pxl_band(x, y, slope, Xoff1, Xoff2, xmin=0)[source]
Chose subset around fault between two parallel lines (X = slope*Y + off).
slope:

slope of bounds for X as a function of Y (float)

Xoff1:

offset in minimum line (float)

Xoff2:

offset in maximum line (Xoff2 > Xoff1)

Opt:
xmin:

if nonzero in spatial_grid function need to specify its value

NOT ADAPTED FOR MPI USE YET

spatial_grid(xmin=0, xmax=0, ymin=0, ymax=0, truncate=False)[source]

Create spatial grid with possibility of choosing a spatial subset.

Opts:
xmin,xmax,ymin,ymax:

indexes delimiting the study area (integers) default are all pixels (0,Nx,0,Ny)

truncate:

True or False (for quick testing)

Return:
  • meshgrid in x and y (2 2D arrays)

summary()[source]

Print data properties

trunc_time(t_num)[source]
Select subset of dates.
t_num:

the number of dates to keep at the end of array

5.3. Create or open storage file

Processing steps contained in the main KFTS-InSAR routine. Storage files, or outputs are described in Section 2.3.

kf.readinput.initiatefileforKF(statefile, phasefile, L, data, model, store, updtfile=None, comm=False, toverlap=0, tshift=1, t_sep=None)[source]
Open h5py file for kalman filter.
statefile:

file name and location

phasefile:

file name and location

L:

number of parameters

model:

model description in tuple used for timefunction.py

store:

is a tuple of things to store

updtfile:

file name to store additional statistics about KF analysis

comm:

communicator if MPI used (e.g. MPI.COMM_WORLD)

toverlap:

number of overlaping timesteps with past solution (only if restart KF)

tshift:

number of previously estimated phases (may be updated)

kf.readinput.reopenfileforKF(statefile, phasefile, comm=False)[source]

Function to reopen file after closure to test readability