4. Kalman Filter

4.1. Theory

See Dalaison & Jolivet, 2020 A Kalman Filter Time Series Analysis method for InSAR, Journal of Geophysical Research - Solid Earth, doi:10.1029/2019JB019150.

4.2. Implementation

4.2.1. Introduction

Below is the class computing the time series analysis pixel by pixel. Once the class is initialized, the main steps of the procedure are as follow:

  • First you start a new time series analysis (start_new) OR you restart an existing one (restart_from_file)

  • Then the function kf loops over (new) timesteps and compute iteratively
    1. the forecast (predict)

    2. the analysis (update)

  • Finally the output for the given pixel is stored (write_output)

4.2.2. Functions

class kf.KF_class.Kalman(data, fctmod, j=0, i=0, verbose=True)[source]

Class for a Kalman filter for an InSAR time series analysis Initialize the object

  • dataobject

    observations/measurements in class from readimput_mpi.py

  • fctmodobject

    functional model in class from timefunction.py

  • j,iintegers

    indexes for 2-D image used for storage

  • verboseboolean

    print stuffs

check_fit(X, P, eps_interf=10)[source]

Check quality of fit of phases if verbose activated. Compute residual weighted by its Covariance for analysed state and print warning if pb

  • eps_interffloat

    accepted difference between computed and real interferograms.

create_H_R_and_D(k, indxs)[source]

Produce the measurement vector (D), the measurement matrix (H), and the measurement covariance matrix (R) at a specific timestep (0≤ k <N) –> if len(D)=n for this timestep, then H will be (n x (L+k+1)) and R (n x n)

  • kinteger

    itteration number

  • indxinteger

    indexes (with respect to t0) of phases in self.m[L:]

create_Q(m_err, phi_err, add_err, M)[source]

Create process covariance Q from uncertainty on model (m_err) and interferograms (phi_err) at kth time.

  • m_errfloat or an array of length L

    model uncertainty

  • phi_errfloat

    systematic error on phases (should be zero)

  • add_errfloat

    systematic error on last forecast (square of std of mismodeling)

  • Minteger

    the state vector length

detect_event(k, kmod, m_all)[source]

IN PROGRESS TESTED ON SYNTHETIC DATA Add model parameter for unexpected events not in model

  • kinteger

    iteration

  • kmodinteger

    minimum k at which modification can be applied

expend_m_P(L, n, PL)[source]

Open state vector and covariance (m and P) to add building parameters

  • Linteger

    index at which we open and insert new parameters in m and P

  • ninteger

    number of parameters to add

  • PLfloat

    apriori variance of the new parameters

get_model_from_num_of_param(N)[source]

Truncate model if the number of parameters in the input (N) is smaller than the maximum number of parameters as predicted by the initial functional model (self.Lref). N (or equivalently self.L) may increase with latter iterations

innovation(Xf, Y)[source]

Compute residual or innovation vector Innovation for phases is not informative. After a few steps, reflects noise of data around model.

kf(m_err, phi_err, add_err, t_sep=6, plots=True, cm='jet', ax1=0, ax2=0)[source]

Run kalman filter combining other functions of class (i.e. MAIN)

  • m_errarray

    systematic error on model (should be 0)

  • phi_errfloat

    systematic error on interferograms (should be 0)

  • t_sepinteger

    maximum time separation between interferograms, fix the minimum number of phases that must be kept in the state vector. Constrain the maximum length of the state vector

  • plotsboolean, optional

    WARNING - activate only if one instance of KF (=one pixel), then subsequent parameters must be specified

    • ax1pyplot axis

      in which plot evolution of parameters

    • ax2pyplot axis

      in which plot evolution of predicted value and model

    • cmstring or colormap

      the colormap of reference later discretised

plot_gain(k, ax, cmap)[source]

Plot gain for each parameter over time

plot_model(k, ax, cmap)[source]

Plot resulting model

plot_params(k, ax, cmap)[source]

Plot each parameter over time with its uncertainty in subplots of size len(ax)

predict(X, P, A, Q)[source]

Forecast step

  • Xarray (N)

    The mean state estimate of the previous step ( k −1).

  • Parray (N, N)

    The state covariance of previous step ( k −1).

  • Aarray (N+1, N)

    The transition matrix.

  • Qarray (N+1, N+1)

    The process noise covariance matrix.

reduce_sizes_m_P(k)[source]

Remove phases in m if not used to build interferograms and has converged For ulterior long baseline interferograms, phase and associated standard deviation can be recovered but state Covariance terms are lost (too heavy to store)

  • kinteger

    number of iteration

restart_from_file(fin, pasttime, indxs, dtmax=None)[source]

Extract initial condition from OPENED infile (fin) which stores previously computed mk and Pk for all pixels including pixel[i,j]

  • finobject

    opened H5 file containing formely computed state

  • pasttimearray

    already loaded time array in fin

  • indxsarray

    already loaded index array

start_new(m0, P0)[source]

Start from skratches

  • m01D array (N)

    Initial Model. The length of the vector will determine how many element of the model will be kept (in the order given in the model vector)

  • P02D array (N,N)

    Initial Covariance

title_labels(ax1)[source]

Add axes label and titles for subplots from plot_params and plot gain functions

update(Xf, Pf)[source]
  • Xfarray

    forecast mean of the state

  • Pfarray

    forecast covariance of the state

write_output(outstates, outphase, outupdate=None)[source]

Save outputs of kalman filter for next run

  • outstatesh5py file

    Open h5file for state storage

  • outphaseh5py file

    Open h5file for phase storage

  • outupdateh5py file

    Open h5file for gain and innovation (Optional)