rsqsim_api.fault.segment
RSQSim fault segment and related classes.
Provides DisplacementArray for storing surface-displacement
observations, RsqSimSegment for a single triangulated fault
segment with geometry, slip-rate properties, and discretisation
utilities, RsqSimFault as a container for one or more
segments, and OpenQuakeSegment for OpenQuake polygon
representations.
Attributes
Classes
Container for 3-component surface-displacement observations at a set of stations. |
|
A single triangulated fault segment for use with RSQSim. |
|
Container for one or more |
|
Lightweight wrapper around a list of polygons for OpenQuake fault representation. |
Module Contents
- class rsqsim_api.fault.segment.DisplacementArray(x_array, y_array, z_array=None, e_array=None, n_array=None, v_array=None)[source]
Container for 3-component surface-displacement observations at a set of stations.
- Parameters:
x_array (numpy.ndarray) – Easting (NZTM x) coordinates of the stations.
y_array (numpy.ndarray) – Northing (NZTM y) coordinates.
z_array (numpy.ndarray or None, optional) – Elevation coordinates. If
None, assumed to be zero.e_array (numpy.ndarray or None, optional) – East-component displacement.
n_array (numpy.ndarray or None, optional) – North-component displacement.
v_array (numpy.ndarray or None, optional) – Vertical-component displacement.
- class rsqsim_api.fault.segment.RsqSimSegment(segment_number, patch_type='triangle', fault_name=None)[source]
A single triangulated fault segment for use with RSQSim.
Stores triangular patch geometry, slip-rate components, adjacency and Laplacian matrices, and surface-trace information. Provides classmethods for constructing a segment from various file formats (triangle arrays, tsurf, DXF, STL, VTK, pandas DataFrames) and methods for exporting to mesh formats and RSQSim input files.
- Parameters:
segment_number (int)
patch_type (str)
fault_name (str)
- patch_outlines[source]
List of
RsqSimTriangularPatchobjects.- Type:
list or None
- vertices[source]
Array of unique vertex coordinates, shape (n_verts, 3).
- Type:
numpy.ndarray or None
- triangles[source]
Triangle connectivity, shape (n_tri, 3), referencing
verticesby index.- Type:
numpy.ndarray or None
- property bounds[source]
Axis-aligned bounding box of all vertices in the XY plane.
- Returns:
[x_min, y_min, x_max, y_max]in NZTM metres.- Return type:
numpy.ndarray of shape (4,)
- generate_triangles()[source]
Build
vertices,triangles, andedge_linesfrom the patch outlines.Deduplicates vertices across all patches, then constructs the indexed triangle and edge-line arrays.
- find_triangles_from_vertex_index(vertex_index)[source]
Return the indices of triangles that contain a given vertex index.
- classmethod from_triangles(triangles, segment_number=0, patch_numbers=None, fault_name=None, strike_slip=None, dip_slip=None, rake=None, total_slip=None, min_patch_area=1.0)[source]
Create a segment from a flat triangle-vertex array.
Populates the segment with
RsqSimTriangularPatchobjects. Strike-slip and dip-slip components can be supplied directly or computed fromtotal_slipandrake. Degenerate patches with area less thanmin_patch_areaare removed.- Parameters:
triangles (array-like of shape (n, 9)) – Flat array of triangle vertex coordinates; each row is
[x1,y1,z1, x2,y2,z2, x3,y3,z3].segment_number (int, optional) – Segment identifier. Defaults to 0.
patch_numbers (array-like or None, optional) – Global patch numbers. If
None, assigned sequentially from 0.fault_name (str or None, optional) – Name for the segment.
strike_slip (float or None, optional) – Strike-slip rate (m/s).
dip_slip (float or None, optional) – Dip-slip rate (m/s).
rake (float, numpy.ndarray, or None, optional) – Rake angle(s) in degrees. Required when
total_slipis given.total_slip (numpy.ndarray or None, optional) – Total slip rate (m/s) per patch; used with
raketo compute the ss/ds components.min_patch_area (float, optional) – Patches smaller than this area (m²) are discarded. Defaults to 1.0.
- Returns:
Populated segment object.
- Return type:
- classmethod from_tsurface(tsurface_file, segment_number=0, patch_numbers=None, fault_name=None, strike_slip=None, dip_slip=None)[source]
Create a segment from a tsurf (.ts) mesh file.
- Parameters:
tsurface_file (str) – Path to the tsurf file.
segment_number (int, optional) – Segment identifier. Defaults to 0.
patch_numbers (array-like or None, optional) – Global patch numbers.
fault_name (str or None, optional) – Name for the segment.
strike_slip (float or None, optional) – Strike-slip rate (m/s).
dip_slip (float or None, optional) – Dip-slip rate (m/s).
- Return type:
- classmethod from_dxf(dxf_file, segment_number=0, patch_numbers=None, fault_name=None, strike_slip=None, dip_slip=None)[source]
Create a segment from a DXF mesh file.
- Parameters:
dxf_file (str) – Path to the DXF file.
segment_number (int, optional) – Segment identifier. Defaults to 0.
patch_numbers (array-like or None, optional) – Global patch numbers.
fault_name (str or None, optional) – Name for the segment.
strike_slip (float or None, optional) – Strike-slip rate (m/s).
dip_slip (float or None, optional) – Dip-slip rate (m/s).
- Return type:
- classmethod from_pandas(dataframe, segment_number, patch_numbers, fault_name=None, strike_slip=None, dip_slip=None, read_rake=True, read_slip_rate=True, transform_from_utm=False)[source]
Create a segment from a pandas DataFrame with triangle and slip-rate columns.
- Parameters:
dataframe (pandas.DataFrame) – DataFrame with the first 9 columns being vertex coordinates (x1,y1,z1,…,x3,y3,z3), plus optional
"rake"and"slip_rate"columns.segment_number (int) – Segment identifier.
patch_numbers (array-like) – Global patch numbers, one per row.
fault_name (str or None, optional) – Name for the segment.
strike_slip (float or None, optional) – Strike-slip rate (m/s); overridden when
read_rakeisTrue.dip_slip (float or None, optional) – Dip-slip rate (m/s); overridden when
read_rakeisTrue.read_rake (bool, optional) – If
True(default), read rake from the"rake"column and compute ss/ds from total slip × rake.read_slip_rate (bool, optional) – If
True(default), read slip rate from the"slip_rate"column.transform_from_utm (bool, optional) – If
True, transform vertex coordinates from UTM zone 59S (EPSG:32759) to NZTM (EPSG:2193). Defaults toFalse.
- Return type:
- classmethod from_pickle(dataframe, segment_number, patch_numbers, fault_name=None)[source]
Create a segment from a pickled patch DataFrame.
- Parameters:
dataframe (pandas.DataFrame) – DataFrame produced by
RsqSimMultiFault.pickle_model(), with columns[vertices, normal_vector, down_dip_vector, dip, along_strike_vector, centre, area, dip_slip, strike_slip, fault_num].segment_number (int) – Segment identifier.
patch_numbers (array-like) – Global patch numbers.
fault_name (str or None, optional) – Name for the segment.
- Return type:
- classmethod from_stl(stl_file, segment_number=0, patch_numbers=None, fault_name=None, strike_slip=None, dip_slip=None, rake=None, total_slip=None)[source]
Create a segment from an STL mesh file.
- Parameters:
stl_file (str) – Path to the STL file.
segment_number (int, optional) – Segment identifier. Defaults to 0.
patch_numbers (array-like or None, optional) – Global patch numbers.
fault_name (str or None, optional) – Name for the segment.
strike_slip (float or None, optional) – Slip-rate components (m/s).
dip_slip (float or None, optional) – Slip-rate components (m/s).
rake (float or None, optional) – Rake angle in degrees.
total_slip (numpy.ndarray or None, optional) – Total slip rates (m/s) per patch.
- Return type:
- classmethod from_vtk(vtk_file, segment_number=0, patch_numbers=None, fault_name=None)[source]
Create a segment from a VTK mesh file containing slip and rake data.
- Parameters:
vtk_file (str) – Path to the VTK file.
segment_number (int, optional) – Segment identifier. Defaults to 0.
patch_numbers (array-like or None, optional) – Global patch numbers.
fault_name (str or None, optional) – Name for the segment.
- Return type:
- build_adjacency_map(verbose=False)[source]
Build a triangle adjacency map (triangles sharing an edge).
For each triangle, finds all other triangles that share exactly two vertices (i.e. a full edge). Populates
adjacency_mapas a list of lists.- Parameters:
verbose (bool, optional) – If
True, print triangle indices as they are appended. Defaults toFalse.
- write_neighbour_file(filename)[source]
Write a text file listing the neighbouring triangle indices for each triangle.
- Parameters:
filename (str) – Output file path. Each line lists the space-separated neighbour indices for one triangle.
- build_laplacian_matrix(double=True, verbose=True)[source]
Build a distance-weighted discrete Laplacian smoothing matrix.
Uses inter-patch centre distances to form the scale-dependent Laplacian operator (Desbrun et al., 1999). Normalises by the maximum absolute diagonal entry.
- Parameters:
double (bool, optional) – If
True(default), horizontally stack two copies of the matrix (for use in slip inversions with separate ss/ds columns). IfFalse, store an N×N matrix inlaplacian_sing.verbose (bool, optional) – If
True(default), print progress messages.
- find_top_vertex_indices(depth_tolerance=100.0, complicated_faults=False)[source]
Return indices of the shallowest vertices within a depth tolerance.
- Parameters:
depth_tolerance (float or int, optional) – Vertical tolerance in metres below the shallowest vertex. Defaults to 100.
complicated_faults (bool, optional) – If
True, use edge-vertex logic to handle non-planar or folded fault geometries. Defaults toFalse.
- Returns:
Vertex indices (into
vertices) of the top edge.- Return type:
numpy.ndarray of int
- find_vertex_indices(depth_tolerance=100.0, complicated_faults=False)[source]
- Parameters:
depth_tolerance (float | int)
complicated_faults (bool)
- find_top_vertices(depth_tolerance=100, complicated_faults=False)[source]
Return 3-D coordinates of the shallowest vertices.
- Parameters:
depth_tolerance (float or int, optional) – Vertical tolerance (m) below the shallowest vertex. Defaults to 100.
complicated_faults (bool, optional) – Passed to
find_top_vertex_indices().
- Return type:
numpy.ndarray of shape (n, 3)
- find_top_edges(depth_tolerance=100, complicated_faults=False)[source]
Return the edge-line pairs (vertex index pairs) forming the top edge of the fault.
- Parameters:
depth_tolerance (float or int, optional) – Vertical tolerance (m) below the shallowest vertex. Defaults to 100.
complicated_faults (bool, optional) – Passed to
find_top_vertex_indices().
- Return type:
numpy.ndarray of shape (n_edges, 2)
- find_top_patch_numbers(depth_tolerance=100)[source]
Return the patch numbers whose triangles include at least one top-edge vertex.
- Parameters:
depth_tolerance (float or int, optional) – Vertical tolerance (m) below the shallowest vertex. Defaults to 100.
- Return type:
list of int
- find_edge_patch_numbers(top=True, depth_tolerance=100)[source]
Return patch numbers for patches that lie on the outer edge of the fault.
- Parameters:
top (bool, optional) – If
True(default), include top-edge patches. IfFalse, exclude top-edge patches (returning only lateral and bottom edge patches).depth_tolerance (float or int, optional) – Vertical tolerance (m) for identifying top patches. Defaults to 100.
- Return type:
numpy.ndarray of int
- find_all_outside_edges()[source]
Return all boundary edges of the triangulated surface (edges belonging to only one triangle).
- Returns:
Each row is a sorted pair of vertex indices.
- Return type:
numpy.ndarray of shape (n_boundary_edges, 2)
- find_top_outside_vertex_indices(depth_tolerance=100)[source]
- Parameters:
depth_tolerance (float | int)
- find_bottom_outside_vertex_indices(depth_tolerance=100)[source]
- Parameters:
depth_tolerance (float | int)
- find_bottom_outside_vertices(depth_tolerance=100)[source]
- Parameters:
depth_tolerance (float | int)
- bottom_edge_point_cloud(depth_tolerance=100, num_interp=25)[source]
- Parameters:
depth_tolerance (float | int)
num_interp (int)
- grid_surface_rbf(resolution)[source]
Interpolate the fault surface depth onto a regular 2-D grid using RBF interpolation.
Projects the fault vertices into a regular (x, y) grid within the bounding box of the fault outline and interpolates the z (depth) coordinate using a radial basis function. Not suitable for near-vertical faults.
- Parameters:
resolution (float) – Grid spacing in NZTM metres.
- Returns:
xx (numpy.ndarray of shape (ny, nx)) – Easting grid.
yy (numpy.ndarray of shape (ny, nx)) – Northing grid.
z (numpy.ndarray of shape (ny, nx)) – Interpolated depth values.
- get_slip_vec_3d()[source]
Return the 3-D slip vector for each patch.
- Return type:
numpy.ndarray of shape (n_patches, 3)
- plot_2d(ax, cmap='viridis', max_slip=10.0, alpha=0.5, mm_yr=True)[source]
Plot the fault segment in map view coloured by slip rate.
- Parameters:
ax (matplotlib.axes.Axes) – Axes to draw onto.
cmap (str, optional) – Colormap name. Defaults to
"viridis".max_slip (float, optional) – Maximum slip rate for the colour scale. Defaults to 10.
alpha (float, optional) – Transparency (0–1). Defaults to 0.5.
mm_yr (bool, optional) – If
True(default), convert slip rate from m/s to mm/yr before plotting.
- Returns:
The tripcolor artist.
- Return type:
matplotlib.collections.PolyCollection
- to_mesh(write_slip=False, convert_to_mm_yr=False)[source]
Convert the segment to a
meshio.Meshobject.- Parameters:
write_slip (bool, optional) – If
True, attach slip and rake as cell data. Defaults toFalse.convert_to_mm_yr (bool, optional) – If
True, convert slip from m/s to mm/yr. Defaults toFalse.
- Return type:
meshio.Mesh
- to_stl(stl_name)[source]
Write the segment mesh to an STL file.
- Parameters:
stl_name (str) – Output STL file path.
- to_vtk(vtk_name, write_slip=False, convert_to_mm_yr=False)[source]
Write the segment mesh to a VTK file.
- Parameters:
vtk_name (str) – Output VTK file path.
write_slip (bool, optional) – Attach slip/rake as cell data. Defaults to
False.convert_to_mm_yr (bool, optional) – Convert slip to mm/yr. Defaults to
False.
- to_gpd(write_slip=False, crs=2193)[source]
Convert the segment to a GeoDataFrame of patch polygons.
- Parameters:
write_slip (bool, optional) – If
True, include"slip"and"rake"columns. Defaults toFalse.crs (int, optional) – EPSG code for the output CRS. Defaults to 2193 (NZTM).
- Return type:
geopandas.GeoDataFrame
- to_rsqsim_fault_file(flt_name, mm_yr=True)[source]
Write the segment to an RSQSim fault input file.
- Parameters:
flt_name (str) – Output file path.
mm_yr (bool, optional) – If
True(default), interpret slip rates as mm/yr (converts to m/s for the output file).
- to_rsqsim_fault_array(mm_yr=True)[source]
Build a DataFrame in RSQSim fault-file format for this segment.
Columns are the 9 vertex coordinates, rake, slip rate, fault number, and fault name.
- Parameters:
mm_yr (bool, optional) – If
True(default), treat stored slip rates as mm/yr and convert to m/s.- Return type:
pandas.DataFrame
- get_patch_centres()[source]
Return the centre coordinates of all patches.
- Return type:
numpy.ndarray of shape (n_patches, 3)
- get_average_dip(approx_spacing=5000.0)[source]
Estimate the average dip of the fault by fitting lines to along-dip cross-sections.
Samples cross-sections at approximately
approx_spacingintervals along the trace and returns the median dip.- Parameters:
approx_spacing (float, optional) – Approximate along-strike spacing (m) for cross-sections. Defaults to 5000.
- Returns:
Median dip angle in degrees.
- Return type:
float
- discretize_rectangular_tiles(tile_size=5000.0, interpolation_distance=1000.0)[source]
Discretize the fault into rectangular tiles of approximately uniform size.
Samples cross-sections along the trace, interpolates the down-dip profile, fits local best-fit planes, and constructs four-corner rectangular tile arrays centred on each interpolated point.
- Parameters:
tile_size (float, optional) – Target tile dimension along strike and down-dip (m). Defaults to 5000.
interpolation_distance (float, optional) – Spacing (m) used when interpolating the down-dip profile. Defaults to 1000.
- Returns:
Corner coordinates of each rectangular tile in NZTM (m). NaN-containing tiles are removed if any are produced.
- Return type:
numpy.ndarray of shape (n_tiles, 4, 3)
- class rsqsim_api.fault.segment.RsqSimFault(segments)[source]
Container for one or more
RsqSimSegmentobjects representing a single fault.- Parameters:
segments (RsqSimSegment or list of RsqSimSegment) – One or more segments composing this fault.