TriangularPatches class

class csi.TriangularPatches(name, utmzone=None, ellps='WGS84', verbose=True, lon0=None, lat0=None)

Classes implementing a fault made of triangular patches. Inherits from Fault

Args:
  • name : Name of the fault.

Kwargs:
  • utmzone : UTM zone (optional, default=None)

  • lon0 : Longitude of the center of the UTM zone

  • lat0 : Latitude of the center of the UTM zone

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

  • verbose : Speak to me (default=True)

AverageAlongStrikeOffsets(name, insars, filename, discretized=True, smooth=None)

!Untested in a looong time…!

If the profiles have the lon lat vectors as the fault, This routines averages it and write it to an output file.

ExtractAlongStrikeAllDepths(filename=None, discret=0.5)

Extracts the Along Strike Variations of the creep at all depths for the discretized fault trace.

Kwargs:
  • filename : Name of the output file

  • discret : Fault discretization

Returns:
  • None

ExtractAlongStrikeVariations(depth=0.5, origin=None, filename=None, orientation=0.0)

Extract the Along Strike Variations of the creep at a given depth

Kwargs:
  • depth : Depth at which we extract the along strike variations of slip.

  • origin : Computes a distance from origin. Give [lon, lat].

  • filename: Saves to a file.

  • orientation: defines the direction of positive distances.

Returns:
  • None

ExtractAlongStrikeVariationsOnDiscretizedFault(depth=0.5, filename=None, discret=0.5)

! Untested in a looong time !

Extracts the Along Strike variations of the slip at a given depth, resampled along the discretized fault trace.

Kwargs:
  • depth : Depth at which we extract the along strike variations of slip.

  • discret : Discretization length

  • filename : Saves to a file.

Returns:
  • None

addpatch(patch, slip=[0, 0, 0])

Adds a patch to the list.

Args:
  • patch : Geometry of the patch to add (km, not lon lat)

Kwargs:
  • slip : List of the strike, dip and tensile slip.

Returns:
  • None

addpatches(patches)

Adds patches to the list.

Args:
  • patches : List of patch geometries

Returns:
  • None

buildAdjacencyMap(verbose=True)

For each triangle, find the indices of the adjacent (edgewise) triangles.

Kwargs:
  • verbose

Returns:
  • None

buildLaplacian(verbose=True, method=None, irregular=False)

Build a discrete Laplacian smoothing matrix.

Kwargs:
  • verbose : Speak to me

  • method : Not used, here for consistency purposes

  • irregular : Not used, here for consistency purposes

Returns:
  • Laplacian : 2D array

computeArea()

Computes the area of all triangles.

Returns:
  • None

computeSlipDirection(scale=1.0, factor=1.0, ellipse=False, nsigma=1.0, neg_depth=False)

Computes the segment indicating the slip direction.

Kwargs:
  • scale : can be a real number or a string in ‘total’, ‘strikeslip’, ‘dipslip’ or ‘tensile’

  • factor : Multiply by a factor

  • ellipse : Compute the ellipse

  • nsigma : How many times sigma for the ellipse

Returns:
  • None

computetotalslip()

Computes the total slip and stores it self.totalslip

deletepatch(patch, checkVertices=True, checkSlip=False)

Deletes a patch.

Args:
  • patch : index of the patch to remove.

Kwargs:
  • checkVertices : Make sure vertice array corresponds to patch corners

  • checkSlip : Check that slip vector corresponds to patch corners

Returns:
  • None

deletepatches(tutu, checkVertices=True, checkSlip=True)

Deletes a list of patches.

Args:
  • tutu : List of indices

Kwargs:
  • checkVertices : Check and delete if patches are concerned.

  • checkSlip : Check and delete is slip is concerned.

Returns:
  • None

deletevertex(iVertex, checkPatch=True, checkSlip=True)

Delete a Vertex. If some patches are composed of this vertex and checkPatch is True, deletes the patches.

Args:
  • iVertex : index of the vertex to delete

Kwargs:
  • checkPatch : Check and delete if patches are concerned.

  • checkSlip : Check and delete is slip is concerned.

Returns:
  • None

deletevertices(iVertices, checkPatch=True, checkSlip=True)

Deletes some vertices. If some patches are composed of these vertices and checkPatch is True, deletes the patches.

Args:
  • iVertices : List of vertices to delete.

Kwargs:
  • checkPatch : Check and delete if patches are concerned.

  • checkSlip : Check and delete if slip terms are concerned.

Returns:
  • None

distanceMatrix(distance='center', lim=None)

Returns a matrix of the distances between patches.

Kwargs:
  • distance : distance estimation mode

    • center : distance between the centers of the patches.

    • no other method is implemented for now.

  • lim : if not None, list of two float, the first one is the distance above which d=lim[1].

Returns:
  • distances : Array of floats

distancePatchToPatch(patch1, patch2, distance='center', lim=None)

Measures the distance between two patches.

Args:
  • patch1 : first patch or its index

  • patch2 : second patch or its index

Kwargs:
  • distance : distance estimation mode

    • center : distance between the centers of the patches.

    • no other method is implemented for now.

  • lim : if not None, list of two float, the first one is the distance above which d=lim[1].

Returns:
  • distace : float

distanceVertexToVertex(vertex1, vertex2, distance='center', lim=None)

Measures the distance between two vertexes.

Args:
  • vertex1 : first patch or its index

  • vertex2 : second patch or its index

Kwargs:
  • lim : if not None, list of two float, the first one is the distance above which d=lim[1].

  • distance : Useless argument only here for compatibility reasons

Returns:
  • distance : float

distancesMatrix(distance='center', lim=None)

Returns two matrices of the distances between patches. One for the horizontal dimensions, the other for verticals

Kwargs:
  • distance : distance estimation mode

    • center : distance between the centers of the patches.

    • no other method is implemented for now.

  • lim : if not None, list of two float, the first one is the distance above which d=lim[1].

Returns:
  • distances : Array of floats

findAsperities(function, slip='strikeslip', verbose=True)

Finds the number, size and location of asperities that are identified by the given function.

Args:
  • function : Function that takes an array the size of the number of patches and returns an array of bolean the same size. Trues are within the asperity.

Kwargs:
  • slip : Which slip vector do you want to apply the function to

  • verbose : Talk to me?

Returns:
  • Asperities

getDepths(center=True)

Returns an array of depths.

Kwargs:
  • center : If True, returns the center of the patches

getDips()

Returns an array of dips.

getEllipse(patch, ellipseCenter=None, Npoints=10, factor=1.0, nsigma=1.0)

Compute the ellipse error given Cm for a given patch

Args:
  • patch : Which patch to consider

Kwargs:
  • center : center of the ellipse

  • Npoints : number of points on the ellipse

  • factor : scaling factor

  • nsigma : will design a nsigma*sigma error ellipse

Returns:
  • Ellipse : Array containing the ellipse

getStrikes()

Returns an array of strikes.

getSubSourcesFault(verbose=True)

Returns a TriangularPatches fault object with each triangle corresponding to the subsources used for plotting.

Kwargs:
  • verbose : Talk to me (default: True)

Returns:
  • fault : Returns a triangularpatches instance

getcenter(p)

Get the center of one triangular patch.

Args:
  • p : Patch geometry.

Returns:
  • x,y,z : floats

getcenters(coordinates='xy')

Get the center of the patches.

Kwargs:
  • coordinates : xy or ll

Returns:
  • centers: list of triplets

getpatchgeometry(patch, center=False, retNormal=False, checkindex=True)

Returns the patch geometry as needed for triangleDisp.

Args:
  • patch : index of the wanted patch or patch

Kwargs:
  • center : if true, returns the coordinates of the center of the patch. if False, returns the first corner

  • checkindex : Checks the index of the patch

  • retNormal : If True gives, also the normal vector to the patch

Returns:
  • x, y, z, width, length, strike, dip, (normal)

homogeneizeStrike(direction=1, sign=1.0)

Rotates the vertices in {Faces} so that the normals are all pointing in a common direction. The {direction} is checked by the axis (0=x, 1=y, 2=z).

Kwargs:
  • direction : Direction of the normal to check.

  • sign : +1 or -1

mapFault2Fault(Map, fault)

User provides a Mapping function np.array((len(self.patch), len(fault.patch))) and a fault and the slip from the argument fault is mapped into self.slip.

Function just does self.slip[:,0] = np.dot(Map,fault.slip)

mapUnder2Above(deepfault)

This routine is very very particular. It only works with 2 vertical faults. It Builds the mapping function from one fault to another, when these are vertical. These two faults must have the same surface trace. If the deep fault has more than one raw of patches, it might go wrong and give some unexpected results.

Args:
  • deepfault : Deep section of the fault.

patchArea(patch)

Returns the area of one patch.

Args:
  • patch : one item of the patch list.

Returns:
  • Area : float

patches2triangles(fault, numberOfTriangles=4)

Takes a fault with rectangular patches and splits them into triangles to initialize self.

Args:
  • fault : instance of rectangular patches.

Kwargs:
  • numberOfTriangles : Split each patch in 2 or 4 (default) triangle

Returns:
  • None

plot(figure=134, slip='total', equiv=False, show=True, Map=True, Fault=True, norm=None, linewidth=1.0, plot_on_2d=True, shadedtopo=None, view=None, alpha=1.0, colorbar=True, cbaxis=[0.1, 0.2, 0.1, 0.02], cborientation='horizontal', cblabel='', drawCoastlines=True, expand=0.2, savefig=False, figsize=(None, None))

Plot the available elements of the fault.

Kwargs:
  • figure : Number of the figure.

  • slip : What slip to plot

  • equiv : useless. For consitency between fault objects

  • show : Show me

  • Norm : colorbar min and max values

  • linewidth : Line width in points

  • plot_on_2d : Plot on a map as well

  • drawCoastline : Self-explanatory argument…

  • expandExpand the map by {expand} degree around the edges

    of the fault.

  • savefig : Save figures as eps.

Returns:
  • None

plotMayavi(neg_depth=True, value_to_plot='total', colormap='jet', reverseSign=False)

! OBSOLETE BUT KEPT HERE TO BE TESTED IN THE FUTURE ! Plot 3D representation of fault using MayaVi.

Args:
  • neg_depth : Flag to specify if patch depths are negative or positive

  • value_to_plot : What to plot on patches

  • colormap : Colormap for patches

  • reverseSign : Flag to reverse sign of value_to_plot

! OBSOLETE BUT KEPT HERE TO BE TESTED IN THE FUTURE !

pointRotation3D(iPatch, iPoint, theta, p_axis1, p_axis2)

Rotate a point with an arbitrary axis (fault tip) Used in rotatePatch

Args:
  • iPatch: index of the patch to be rotated

  • iPoint: index of the patch corner (point) to be rotated

  • theta : angle of rotation in degrees

  • p_axis1 : first point of axis (ex: one side of a fault)

  • p_axis2 : second point to define the axis (ex: the other side of a fault)

Returns:
  • rotated point

Reference: ‘Rotate A Point About An Arbitrary Axis (3D)’ - Paul Bourke

readGocadPatches(filename, neg_depth=False, utm=False, factor_xy=1.0, factor_depth=1.0, verbose=False)

Load a triangulated surface from a Gocad formatted file. Vertices must be in geographical coordinates.

Args:
  • filename: tsurf file to read

Kwargs:
  • neg_depth: if true, use negative depth

  • utm: if true, input file is given as utm coordinates (if false -> lon/lat)

  • factor_xy: if utm==True, multiplication factor for x and y

  • factor_depth: multiplication factor for z

  • verbose: Speak to me

Returns:
  • None

readPatchesFromFile(filename, readpatchindex=True, donotreadslip=False, gmtslip=True, inputCoordinates='lonlat')

Reads patches from a GMT formatted file.

Args:
  • filename : Name of the file

Kwargs:
  • inputCoordinates : Default is ‘lonlat’. Can be ‘utm’

  • readpatchindex : Default True.

  • donotreadslip : Default is False. If True, does not read the slip

  • gmtslip : A -Zxxx is in the header of each patch

  • inputCoordinates : Default is ‘lonlat’, can be ‘xyz’

Returns:
  • None

readPatchesFromTriFile(filename, readpatchindex=False, donotreadslip=True, inputCoordinates='lonlat')

Reads patches from a file with triangles having all three vertices on each line modified from the GMT file reader

Args:
  • filename : Name of the file

Kwargs:
  • inputCoordinates : Default is ‘lonlat’. Can be ‘utm’

  • readpatchindex : Default False.

  • donotreadslip : Default is True. If True, does not read the slip

Returns:
  • None

refineMesh()

Cuts all the patches in 4, based on the mid-point of each triangle and builds a new fault from that.

Returns:
  • None

replacePatch(patch, iPatch)

Replaces one patch by the given geometry.

Args:
  • patch : Patch geometry.

  • iPatch : index of the patch to replace.

Returns:
  • None

rotatePatch(iPatch, theta, p_axis1, p_axis2, verbose=False)

Rotate a patch with an arbitrary axis (fault tip) Used by fault class uncertainties

Args:
  • iPatch: index of the patch to be rotated

  • theta : angle of rotation in degrees

  • p_axis1 : first point of axis (ex: one side of a fault)

  • p_axis2 : second point to define the axis (ex: the other side of a fault)

Returns:
  • rotated patch

selectPatches(minlon, maxlon, minlat, maxlat, mindep, maxdep)

Removes patches that are outside of a 3D box.

Args:
  • minlon : west longitude

  • maxlon : east longitude

  • minlat : south latitude

  • maxlat : north latitude

  • mindep : Minimum depth

  • maxdep : Maximum depth

Returns:
  • None

setMu(model_file)

Set shear modulus values for seismic moment calculation from an edks model_file

Args:
  • model_file: EDKS .model file

Returns:
  • None

setVerticesFromPatches()

Takes the patches and constructs a list of Vertices and Faces

Returns:
  • None

setdepth()

Set depth patch attributes

Returns:
  • None

slip2dis(data, patch, slip=None)

Computes the surface displacement for a given patch at the data location using a homogeneous half-space.

Args:
  • data : data object from gps or insar.

  • patch : number of the patch that slips

Kwargs:
  • slip : if a number is given, that is the amount of slip along strike. If three numbers are given, that is the amount of slip along strike, along dip and opening. if None, values from self.slip are taken.

Returns:
  • ss_dis : Surface displacements due to strike slip

  • ds_dis : Surface displacements due to dip slip

  • ts_dis : Surface displacements due to tensile opening

splitPatch(patch)

Splits a patch into 4 patches, based on the mid-point of each side.

Args:
  • patch : item of the patch list.

Returns:
  • t1, t2, t3, t4 : 4 patches

surfacesimulation(box=None, disk=None, err=None, npoints=10, lonlat=None, slipVec=None)

Takes the slip vector and computes the surface displacement that corresponds on a regular grid.

Kwargs:
  • box : A list of [minlon, maxlon, minlat, maxlat].

  • disk : list of [xcenter, ycenter, radius, n]

  • lonlat : Arrays of lat and lon. [lon, lat]

  • err : Compute random errors and scale them by {err}

  • slipVec : Replace slip by what is in slipVec

translatePatch(iPatch, tr_vector)

Translate a patch Used by class uncertainties

Args:
  • iPatch: index of the patch to be rotated

  • tr_vector: array, translation vector in 3D

Returns:
  • None

vertices2ll()

Converts all the vertices into lonlat coordinates.

Returns:
  • None

writeGocadPatches(filename, utm=False)

Write a triangulated Gocad surface file.

Args:
  • filename : output file name

Kwargs:
  • utm : Write in utm coordinates if True

Returns:
  • None

writePatches2File(filename, add_slip=None, scale=1.0, stdh5=None, decim=1)

Writes the patch corners in a file that can be used in psxyz.

Args:
  • filename : Name of the file.

Kwargs:
  • add_slipPut the slip as a value for the color.

    Can be None, strikeslip, dipslip, total, coupling

  • scale : Multiply the slip value by a factor.

  • patch : Can be ‘normal’ or ‘equiv’

  • stdh5 : Get the standard deviation from a h5 file

  • decim : Decimate the h5 file

Returns:
  • None

writeSlipCenter2File(filename, add_slip=None, scale=1.0, neg_depth=False)

Write a psxyz, surface or greenspline compatible file with the center of each patch and magnitude slip. Scale can be a real number or a string in ‘total’, ‘strikeslip’, ‘dipslip’ or ‘tensile’

Args:
  • filename : Name of the output file

Kwargs:
  • add_slipPut the slip as a value at the center

    Can be None, strikeslip, dipslip, total, coupling

  • scale : Multiply the slip value by a factor.

  • neg_depth : if True, depth is a negative nmber

Returns:
  • None

writeSlipDirection2File(filename, scale=1.0, factor=1.0, neg_depth=False, ellipse=False, nsigma=1.0)

Write a psxyz compatible file to draw lines starting from the center of each patch, indicating the direction of slip. Scale can be a real number or a string in ‘total’, ‘strikeslip’, ‘dipslip’ or ‘tensile’

Args:
  • filename : Name of the output file

Kwargs:
  • scale : Scale of the line

  • factor : Multiply slip by a factor

  • neg_depth : if True, depth is a negative nmber

  • ellipse : Write the ellipse

  • nsigma : Nxsigma for the ellipse

Returns:
  • None