insar class

class csi.insar(name, utmzone=None, ellps='WGS84', verbose=True, lon0=None, lat0=None)
Args:
  • name : Name of the InSAR dataset.

Kwargs:
  • utmzone : UTM zone. (optional, default is 10 (Western US))

  • lon0 : Longitude of the utmzone

  • lat0 : Latitude of the utmzone

  • ellps : ellipsoid (optional, default=’WGS84’)

Returns:
  • None

ModelResolutionDownsampling(faults, threshold, damping, startingsize=10.0, minimumsize=0.5, tolerance=0.1, plot=False)

Downsampling algorythm based on Lohman & Simons, 2005, G3.

Args:
  • faults : List of faults, these need to have a buildGFs routine (ex: for RectangularPatches, it will be Okada).

  • threshold : Resolution threshold, if above threshold, keep dividing.

  • damping : Damping parameter. Damping is enforced through the addition of a identity matrix.

Kwargs:
  • startingsize : Starting size of the downsampling boxes.

  • minimumsize : Minimum window size (km)

  • tolerance : Tolerance on the window size calculation

  • plot : True/False

Returns:
  • None

alongStrikeStep(fault, offset, width, binning, discretized=True, strike='mean', returnBox=False)

Computes the difference in values along 2 lines, parallel to the surface trace of a fault.

Args:
  • fault : fault object with a trace

  • offset : offset between each line and the fault

  • width : width of the averaging box

  • binning : bins width (in km)

Kwargs:
  • disctretized : Use the discretized trace of the fault

  • strike : ‘mean’ (mean strike) or a number between 0 and 2pi

  • returnBox : Return the boxes i which data are taken

buildCd(sigma, lam, function='exp', diagonalVar=False, normalizebystd=False)

Builds the full Covariance matrix from values of sigma and lambda.

If function=’exp’:

\(C_d(i,j) = \sigma^2 e^{-\frac{d[i,j]}{\lambda}}\)

elif function=’gauss’:

\(C_d(i,j) = \sigma^2 e^{-\frac{d_{i,j}^2}{2*\lambda}}\)

Args:
  • sigma : Sigma term of the covariance

  • lam : Caracteristic length of the covariance

Kwargs:
  • function : Can be ‘gauss’ or ‘exp’

  • diagonalVar : Substitute the diagonal by the standard deviation of the measurement squared

  • normalizebystd : Weird option to normalize the covariance matrix by the std deviation

Returns:
  • None

buildDiagCd()

Builds a full Covariance matrix from the uncertainties. The Matrix is just a diagonal matrix.

buildsynth(faults, direction='sd', poly=None, vertical=True, custom=False, computeNormFact=True)

Computes the synthetic data using either the faults and the associated slip distributions or the pressure sources.

Args:
  • faults : List of faults or pressure sources.

Kwargs:
  • direction : Direction of slip to use or None for pressure sources.

  • poly : if a polynomial function has been estimated, build and/or include

  • vertical : always True. Used here for consistency among data types

  • custom : if True, uses the fault.custom and fault.G[data.name][‘custom’] to correct

  • computeNormFact : if False, uses TransformNormalizingFactor set with self.setTransformNormalizingFactor

Returns:
  • None

checkLOS(figure=1, factor=100.0, decim=1)

Plots the LOS vectors in a 3D plot.

Kwargs:
  • figure: Figure number.

  • factor: Increases the size of the vectors.

  • decim : Do not plot all the pixels (takes way too much time)

Returns:
  • None

checkLand(resolution='fine', reject=True)

Checks whether pixels are inland or in water.

Args:
  • resolution: Cartopy resolution. Can be fine, coarse or auto

  • refect: deletes the pixels if True. If False, returns a list of booleans

checkNaNs()

Checks and remove data points that have NaNs in vel, err, lon, lat or los.

checkZeros()

Checks and remove data points that have Zeros in vel, lon or lat

cleanProfile(name, xlim=None, zlim=None)

Cleans a specified profile.

Args:
  • name : name of the profile to work with

Kwargs:
  • xlim : tuple (xmin, xmax). Removes pixels outside of this box

  • zlim : tuple (zmin, zmax). Removes pixels outside of this box

Returns:
  • None

computeCustom(fault)

Computes the displacements associated with the custom green’s functions.

Args:
  • fault : A fault instance

Returns:
  • None

computePoly(fault, computeNormFact=True)

Computes the orbital bias estimated in fault

Args:
  • fault : Fault object that has a polysol structure.

Kwargs:
  • computeNormFact : Recompute the norm factors

Returns:
  • None

computeTransformNormalizingFactor()

Compute orbit normalizing factors and store them in insar object.

Returns:
  • None

curve2prof(xl, yl, width, widthDir)

Routine returning the profile along a curve. !!!! Not tested in a long time !!!!

Args:
  • xl : List of the x coordinates of the line.

  • yl : List of the y coordinates of the line.

  • width : Width of the zone around the line.

  • widthDir : Direction to of the width.

Returns:
  • None

deletePixels(u)

Delete the pixels indicated by index in u.

Args:
  • u : array of indexes

Returns:
  • None

distance2point(lon, lat)

Returns the distance of all pixels to a point.

Args:
  • lon : Longitude of a point

  • lat : Latitude of a point

Returns:
  • array

distancePixel2Pixel(i, j)

Returns the distance in km between two pixels.

Args:
  • i : index of a pixel

  • h : index of a pixel

Returns:
  • float

extractAroundGPS(inp, distance, doprojection=True, verbose=False)

Returns a gps object with values projected along the LOS around the gps stations included in gps. In addition, it projects the gps displacements along the LOS

Args:
  • inp : gps or gpstimeseries object

  • distance : distance to consider around the stations

Kwargs:
  • doprojection : Projects the gps enu disp into the los as well

Retunrs:
  • gps instance

get2DstrainEst(computeNormFact=True)

Returns the matrix to estimate the 2d aerial strain tensor. First column is the Epsilon_xx component Second column is the Epsilon_xy component Fourth column is the Epsilon_yy component.

Kwargs:
  • computeNormFact : Recompute the normalization factor

Returns:
  • None

getDistance2Faults(faults)

Returns the minimum distance to a fault for each pixel

Args:
  • faults : a fault or a list of faults.

Return:
  • distance : array

getMisfit()

Computes the Summed Misfit of the data and if synthetics are computed, the RMS of the residuals

Returns:
  • float, float

getOffsetFault(profile, fault, distance, threshold=25, discretized=True, verbose=False, useerrors=True)

Computes the offset next to a fault along a profile. This is a copy from a script written by Manon Dalaison during her PhD in 2020 (infamous year).

Args:
  • profile : Name of the profile

  • fault : fault object

  • distance : Distance to take on each side of the fault (min, max)

Kwargs:
  • threshold : minimum number of data points on each side to move on

Returns:
  • offset, uncertainty, los

getPolyEstimator(ptype, computeNormFact=True)

Returns the Estimator for the polynomial form to estimate in the InSAR data.

Args:
  • ptype: integer
    • 1: constant offset to the data

    • 3: constant and linear function of x and y

    • 4: constant, linear term and cross term.

Kwargs:
  • computeNormFact : bool

Returns:
  • None

getRMS()

Computes the RMS of the data and if synthetics are computed, the RMS of the residuals

Returns:
  • float, float

getStepFromProfile(profname, fault, leftdistance, rightdistance, discretized=False, method='mean')

Returns the offset across a fault from a profile.

Args:

profname : Name of the profile fault : fault object leftdistance: tuple of distances between which to average displacement to the left rightdistance: tuple of distances between which to average displacement to the right

Returns:

step, std, los

getTransformEstimator(trans, computeNormFact=True)

Returns the Estimator for the transformation to estimate in the InSAR data.

Args:
  • transTransformation type
    • 1: constant offset to the data

    • 3: constant and linear function of x and y

    • 4: constant, linear term and cross term.

    • strain: Estimates an aerial strain tensor

Kwargs:
  • computeNormFact : Recompute the normalization factor

Returns:
  • None

getVariance()

Computes the Variance of the data and if synthetics are computed, the RMS of the residuals

Returns:
  • float, float

getprofile(name, loncenter, latcenter, length, azimuth, width)

Project the SAR velocities onto a profile. Works on the lat/lon coordinates system. Profile is stored in a dictionary called {self}.profile

Args:
  • name : Name of the profile.

  • loncenter : Profile origin along longitude.

  • latcenter : Profile origin along latitude.

  • length : Length of profile.

  • azimuth : Azimuth in degrees.

  • width : Width of the profile.

Returns:
  • None

getprofileAlongCurve(name, lon, lat, width, widthDir)

Project the SAR velocities onto a profile. Works on the lat/lon coordinates system. Has not been tested in a long time…

Args:
  • name : Name of the profile.

  • lon : Longitude of the Line around which we do the profile

  • lat : Latitude of the Line around which we do the profile

  • width : Width of the zone around the line.

  • widthDir : Direction to of the width.

Returns:
  • None

incaz2los(incidence, azimuth, origin='onefloat', dtype=<class 'numpy.float32'>)

From the incidence and the heading, defines the LOS vector.

Args:
  • incidence : Incidence angle.

  • azimuth : Azimuth angle of the LOS

Kwargs:
  • originWhat are these numbers
    • onefloat : One number

    • grd : grd files

    • binary : Binary files

    • binaryfloat : Arrays of float

  • dtype : Data type (default is np.float32)

Returns:
  • None

inchd2los(incidence, heading, origin='onefloat')

From the incidence and the heading, defines the LOS vector.

Args:
  • incidence : Incidence angle.

  • heading : Heading angle.

Kwargs:
  • originWhat are these numbers
    • onefloat : One number

    • grd : grd files

    • binary : Binary files

    • binaryfloat : Arrays of float

Returns:
  • None

intersectProfileFault(name, fault, discretized=False)

Gets the distance between the fault/profile intersection and the profile center.

Args:
  • name : name of the profile.

  • fault : fault object from verticalfault.

Returns:
  • float

keepPixels(u)

Keep the pixels indexed u and ditch the other ones

Args:
  • u : array of indexes

Returns:
  • None

mergeInsar(sar)

Combine the existing data set with another insar object.

Args:
  • sar : Instance of the insar class

Returns:
  • None

pickleProfiles(names, filename)

Add on by M. Dalaison (Jan 2019) To save in one file several profiles drawn in one image under the form of a list of python dictionaries

Args:
  • names : profile name or list of profile names (list of strings)

  • filename : for storage of the list of profile object

plot(faults=None, figure=None, gps=None, norm=None, data='data', show=True, Map=True, Fault=True, lognorm=False, drawCoastlines=True, expand=0.2, edgewidth=1, figsize=None, markersize=1.0, plotType='scatter', cmap='jet', alpha=1.0, box=None, titleyoffset=1.1, landcolor='lightgrey', seacolor=None, shadedtopo=None, title=False, los=None, colorbar=True, cbaxis=[0.1, 0.2, 0.1, 0.02], cborientation='horizontal', cblabel='')

Plot the data set, together with a fault, if asked.

Kwargs:
  • faults : list of fault objects.

  • figure : number of the figure.

  • gps : list of gps objects.

  • norm : colorbar limits

  • data : ‘data’, ‘synth’ or ‘res’

  • show : bool. Show on screen?

  • drawCoastlines : bool. default is True

  • expand : default expand around the limits covered by the data

  • edgewidth : width of the edges of the decimation process patches

  • plotType : ‘decim’, ‘scatter’ or ‘flat’

  • figsize : tuple of figure sizes

  • box : Lon/lat box [lonmin, lonmax, latmin, latmax]

Returns:
  • None

plotprofile(name, figsize=(10, 10), fault=None, norm=None, synth=False, cbaxis=[0.1, 0.1, 0.1, 0.01], alpha=0.3, plotType='scatter', drawCoastlines=True)

Plot profile.

Args:
  • name : Name of the profile.

Kwargs:
  • legendscale: Length of the legend arrow.

  • fault : Fault object

  • norm : Colorscale limits

  • synth : Plot synthetics (True/False).

Returns:
  • None

read_from_ascii(filename, factor=1.0, step=0.0, header=0)

Read the InSAR data from an ascii file.

Args:
  • filename : Name of the input file. Format is Lon, Lat, data, uncertainty, los E, los N, los U.

Kwargs:
  • factor : Factor to multiply the LOS velocity.

  • step : Add a value to the velocity.

  • header : Size of the header.

Returns:
  • None

read_from_ascii_simple(filename, factor=1.0, step=0.0, header=0, los=None)

Read the InSAR data from an ascii file with 3 cols.

Args:
  • filename : Name of the input file (format is Lon, Lat. data)

Kwargs:
  • factor : Factor to multiply the LOS velocity.

  • step : Add a value to the velocity.

  • header : Size of the header.

  • los : LOS unit vector (3 column array)

Returns:
  • None

read_from_binary(data, lon, lat, err=None, factor=1.0, step=0.0, incidence=None, heading=None, azimuth=None, los=None, dtype=<class 'numpy.float32'>, remove_nan=True, downsample=1, remove_zeros=True)

Read from binary file or from array.

Args:
  • data : binary array containing the data or binary file

  • lon : binary arrau containing the longitude or binary file

  • lat : binary array containing the latitude or binary file

Kwargs:
  • err : Uncertainty (array)

  • factor : multiplication factor (default is 1.0)

  • step : constant added to the data (default is 0.0)

  • incidence : incidence angle (degree)

  • heading : heading angle (degree)

  • azimuth : Azimuth angle (degree)

  • los : LOS unit vector 3 component array (3-column array)

  • dtype : data type (default is np.float32 if data is a file)

  • remove_nan : True/False

  • downsample : default is 1 (take one pixel out of those)

  • remove_zeros : True/False

Return:
  • None

read_from_grd(filename, factor=1.0, step=0.0, incidence=None, heading=None, los=None, keepnans=False)

Reads velocity map from a grd file.

Args:
  • filename : Name of the input file

Kwargs:
  • factor : scale by a factor

  • step : add a value.

  • incidence : incidence angle (degree)

  • heading : heading angle (degree)

  • los : LOS unit vector (3 column array)

  • keepnans : True/False

Returns:
  • None

read_from_mat(filename, factor=1.0, step=0.0, incidence=35.88, heading=-13.115)

Reads velocity map from a mat file.

Args:
  • filename : Name of the input matlab file

Kwargs:
  • factor : scale by a factor.

  • step : add a step.

  • incidence : incidence angle (degree)

  • heading : heading angle (degree)

Returns:
  • None

read_from_varres(filename, factor=1.0, step=0.0, header=2, cov=False)

Read the InSAR LOS rates from the VarRes output.

Args:
  • filename : Name of the input file. Two files are opened filename.txt and filename.rsp.

Kwargs:
  • factor : Factor to multiply the LOS velocity.

  • step : Add a value to the velocity.

  • header : Size of the header.

  • cov : Read an additional covariance file (binary np.float32, Nd*Nd elements).

Returns:
  • None

reference2area(lon, lat, radius)

References the data to an area. Selects the area and sets all dates to zero at this area.

Args:
  • lon : longitude of the center of the area

  • lat : latitude of the center of the area

  • radius : Radius of the area

Returns:
  • None

referenceProfile(name, xmin, xmax, method='mean')

Removes the mean value of points between xmin and xmax. Optionally, can remove the linear best fit between these 2 values

Args:
  • name : Name of the profile

  • xmin : minimum value along X-axis to consider

  • xmax : Maximum value along X-axis to consider

Kwargs:
  • method : ‘mean’ or ‘linear’

Returns:
  • None

reject_pixel(u)

Reject pixels.

Args:
  • u : Index of the pixel to reject.

Returns:
  • None

reject_pixels_fault(dis, faults)

Rejects the pixels that are {dis} km close to the fault.

Args:
  • dis : Threshold distance.

  • faults : list of fault objects.

Returns:
  • None

removePoly(fault, verbose=False, custom=False, computeNormFact=True)

Removes a polynomial from the parameters that are in a fault.

Args:
  • fault : a fault instance

Kwargs:
  • verbose : Show us stuff

  • custom : include custom green’s functions

  • computeNormFact : recompute norm factors

Returns:
  • None

removeSynth(faults, direction='sd', poly=None, vertical=True, custom=False, computeNormFact=True)

Removes the synthetics using the faults and the slip distributions that are in there.

Args:
  • faults : List of faults.

Kwargs:
  • direction : Direction of slip to use.

  • poly : if a polynomial function has been estimated, build and/or include

  • vertical : always True - used here for consistency among data types

  • custom : if True, uses the fault.custom and fault.G[data.name][‘custom’] to correct

  • computeNormFact : if False, uses TransformNormalizingFactor set with self.setTransformNormalizingFactor

Returns:
  • None

removeTransformation(fault, verbose=False, custom=False)

Wrapper of removePoly to ensure consistency between data sets.

Args:
  • fault : a fault instance

Kwargs:
  • verbose : talk to us

  • custom : Remove custom GFs

Returns:
  • None

returnAverageNearPoint(lon, lat, distance)

Returns the phase value, the los and the errors averaged over a distance from a point (lon, lat).

Args:
  • lon : longitude of the point

  • lat : latitude of the point

  • distance : distance around the point

Returns:
  • float, float, tuple

select_pixels(minlon, maxlon, minlat, maxlat)

Select the pixels in a box defined by min and max, lat and lon.

Args:
  • minlon : Minimum longitude.

  • maxlon : Maximum longitude.

  • minlat : Minimum latitude.

  • maxlat : Maximum latitude.

Retunrs:
  • None

setGFsInFault(fault, G, vertical=True)

From a dictionary of Green’s functions, sets these correctly into the fault object fault for future computation.

Args:
  • fault : Instance of Fault

  • G : Dictionary with 3 entries ‘strikeslip’, ‘dipslip’ and ‘tensile’. These can be a matrix or None.

Kwargs:
  • vertical : Set here for consistency with other data objects, but will always be set to True, whatever you do.

Returns:
  • None

setTransformNormalizingFactor(x0, y0, normX, normY, base)

Set orbit normalizing factors in insar object.

Args:
  • x0 : Reference point x coordinate

  • y0 : Reference point y coordinate

  • normX : Normalization distance along x

  • normY : Normalization distance along y

  • base : baseline

Returns:
  • None

smoothProfile(name, window, method='mean')

Computes smoothed profile but running a mean or median window on the profile.

Args:
  • name : Name of the profile to work on

  • window : Width of the window (km)

Kwargs:
  • method : ‘mean’ or ‘median’

Returns:
  • None

write2file(fname, data='data', outDir='./', err=False)

Write to an ascii file

Args:
  • fname : Filename

Kwargs:
  • data : can be ‘data’, ‘synth’ or ‘resid’

  • outDir : output Directory

Returns:
  • None

write2grd(fname, oversample=1, data='data', interp=100, cmd='surface', tension=None, useGMT=False, verbose=False, outDir='./')

Write to a grd file.

Args:
  • fname : Filename

Kwargs:
  • oversample : Oversampling factor.

  • data : can be ‘data’ or ‘synth’.

  • interp : Number of points along lon and lat.

  • cmd : command used for the conversion( i.e., surface or xyz2grd)

  • tension : Tension in the gmt command

  • useGMT : use surface or xyz2grd or the direct wrapper to netcdf

  • verbose : Talk to us

  • outDir : output directory

Returns:
  • None

writeAlongStrikeOffsets2File(name, filename)

Write the variations of the offset along strike in a file.

Args:
  • name : name of the profile to work on

  • filename : output file name

Returns:
  • None

writeDecim2file(filename, data='data', outDir='./')

Writes the decimation scheme to a file plottable by GMT psxy command.

Args:
  • filename : Name of the output file (ascii file)

Kwargs:
  • data : Add the value with a -Z option for each rectangle. Can be ‘data’, ‘synth’, ‘res’, ‘transformation’

  • outDir : output directory

Returns:
  • None

writeDownsampled2File(prefix, rsp=False)

Writes the downsampled image data to a file. The file will be called prefix.txt. If rsp is True, then it writes a file called prefix.rsp containing the boxes of the downsampling. If prefix has white spaces, those are replaced by “_”.

Args:
  • prefix : Prefix of the output file

Kwargs:
  • rsp : Write the rsp file?

Returns:
  • None

writeEDKSdata()

*Obsolete*

writeProfile2File(name, filename, fault=None)

Writes the profile named ‘name’ to the ascii file filename.

Args:
  • name : Name of the profile to work with

  • filename : Output file name

Kwargs:
  • fault : Instance of fault. Can be a list of faults. Adds the intersection with the fault in the header of the file.

Returns:
  • None