rsqsim_api.catalogue.event

RSQSim earthquake event representation.

Provides RsqSimEvent, which stores the catalogue parameters (time, magnitude, location) and the per-patch slip distribution for a single simulated rupture, together with methods for computing derived quantities and exporting slip distributions in a variety of formats.

Also provides OpenQuakeMultiSquareRupture for packaging events as multi-planar ruptures compatible with the OpenQuake engine.

Attributes

transformer_nztm2wgs

Classes

RsqSimEvent

A single earthquake event from an RSQSim catalogue.

OpenQuakeMultiSquareRupture

Multi-planar rupture representation for the OpenQuake engine.

Module Contents

rsqsim_api.catalogue.event.transformer_nztm2wgs[source]
class rsqsim_api.catalogue.event.RsqSimEvent[source]

A single earthquake event from an RSQSim catalogue.

Stores basic seismic catalogue parameters and, when populated, the patch-level slip distribution including patch objects, slip magnitudes, and rupture timing.

event_id[source]

Unique integer identifier for the event in the catalogue.

Type:

int or None

t0[source]

Origin time in seconds from the start of the simulation.

Type:

float or None

m0[source]

Scalar seismic moment in N·m.

Type:

float or None

mw[source]

Moment magnitude.

Type:

float or None

x, y, z

Hypocentre coordinates in NZTM (metres, EPSG:2193). z is negative downward.

Type:

float or None

area[source]

Total rupture area in m².

Type:

float or None

dt[source]

Rupture duration in seconds.

Type:

float or None

patches[source]

List of RsqSimTriangularPatch objects that slipped in this event.

Type:

list or None

patch_slip[source]

Slip magnitude (m) for each entry in patches.

Type:

numpy.ndarray or None

faults[source]

Unique RsqSimSegment objects involved in this event.

Type:

list or None

patch_time[source]

Absolute rupture time (s) for each patch.

Type:

numpy.ndarray or None

patch_numbers[source]

Global patch IDs corresponding to patches.

Type:

numpy.ndarray or None

length[source]

Estimated rupture length (m) along fault traces.

Type:

float or None

event_id = None[source]
t0 = None[source]
m0 = None[source]
mw = None[source]
area = None[source]
dt = None[source]
patches = None[source]
patch_slip = None[source]
faults = None[source]
patch_time = None[source]
patch_numbers = None[source]
length = None[source]
property num_faults[source]

Number of unique fault segments involved in this event.

property bounds[source]

Axis-aligned bounding box of the rupture in NZTM coordinates.

Returns:

[x_min, y_min, x_max, y_max] in metres (NZTM).

Return type:

list of float

property exterior[source]

Shapely geometry representing the union of all ruptured patch polygons.

Returns:

Union of all patch outlines as a Shapely geometry object.

Return type:

shapely.geometry.base.BaseGeometry

property mean_slip[source]

Mean slip (m) averaged over all ruptured patches.

property mean_strike[source]

Mean strike (degrees, 0–360) averaged over all ruptured patches.

property mean_strike_180[source]

Mean strike (degrees, 0–180) averaged over all ruptured patches.

property mean_dip[source]

Mean dip (degrees) averaged over all ruptured patches.

property mean_rake[source]

Mean rake (degrees) averaged over all ruptured patches.

find_first_fault(fault_model, name=True)[source]

Return the fault that first ruptured in this event.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model used to look up patch-to-fault associations.

  • name (bool, optional) – If True (default), return the fault name string. If False, return the fault segment object.

Returns:

The name (or object) of the fault that hosted the first ruptured patch.

Return type:

str or RsqSimSegment

classmethod from_catalogue_array(t0, m0, mw, x, y, z, area, dt, event_id=None)[source]

Construct an event from basic catalogue parameters.

Creates an event with catalogue-level metadata only (no patch slip distribution).

Parameters:
  • t0 (float) – Origin time in seconds from the simulation start.

  • m0 (float) – Scalar seismic moment in N·m.

  • mw (float) – Moment magnitude.

  • x (float) – Hypocentre easting in NZTM (m).

  • y (float) – Hypocentre northing in NZTM (m).

  • z (float) – Hypocentre depth in NZTM (m, negative downward).

  • area (float) – Total rupture area in m².

  • dt (float) – Rupture duration in seconds.

  • event_id (int, optional) – Unique catalogue event ID.

Returns:

Event populated with catalogue parameters; patch attributes remain None.

Return type:

RsqSimEvent

property patch_outline_gs[source]

2193).

Returns:

One Shapely polygon per ruptured patch, with CRS set to EPSG:2193 (NZTM).

Return type:

geopandas.GeoSeries

Type:

GeoSeries of individual patch outlines in NZTM (EPSG

classmethod from_earthquake_list(t0, m0, mw, x, y, z, area, dt, patch_numbers, patch_slip, patch_time, fault_model, filter_single_patches=True, min_patches=10, min_slip=1, event_id=None)[source]

Construct an event from catalogue parameters and patch slip data.

Filters out fault segments that contribute fewer than min_patches patches, then resolves patch IDs to patch objects and fault objects using fault_model.

Parameters:
  • t0 (float) – Origin time in seconds.

  • m0 (float) – Scalar seismic moment in N·m.

  • mw (float) – Moment magnitude.

  • x (float) – Hypocentre easting in NZTM (m).

  • y (float) – Hypocentre northing in NZTM (m).

  • z (float) – Hypocentre depth in NZTM (m, negative downward).

  • area (float) – Total rupture area in m².

  • dt (float) – Rupture duration in seconds.

  • patch_numbers (array-like) – Global patch IDs that slipped in this event.

  • patch_slip (array-like) – Slip magnitude (m) for each entry in patch_numbers.

  • patch_time (array-like) – Absolute rupture time (s) for each entry in patch_numbers.

  • fault_model (RsqSimMultiFault) – Fault model providing the patch and fault dictionaries.

  • filter_single_patches (bool, optional) – Unused; retained for API compatibility.

  • min_patches (int, optional) – Minimum number of patches a fault must contribute to be retained. Defaults to 10.

  • min_slip (float or int, optional) – Unused; retained for API compatibility.

  • event_id (int, optional) – Unique catalogue event ID.

Returns:

Fully populated event including patch slip distribution and fault objects.

Return type:

RsqSimEvent

classmethod from_multiprocessing(t0, m0, mw, x, y, z, area, dt, patch_numbers, patch_slip, patch_time, fault_model, mask, event_id=None)[source]

Construct an event from pre-computed multiprocessing results.

Similar to from_earthquake_list(), but accepts a boolean mask that has already been applied to filter patches, rather than recomputing it here.

Parameters:
  • t0 (float) – Origin time in seconds.

  • m0 (float) – Scalar seismic moment in N·m.

  • mw (float) – Moment magnitude.

  • x (float) – Hypocentre easting in NZTM (m).

  • y (float) – Hypocentre northing in NZTM (m).

  • z (float) – Hypocentre depth in NZTM (m, negative downward).

  • area (float) – Total rupture area in m².

  • dt (float) – Rupture duration in seconds.

  • patch_numbers (array-like) – Global patch IDs that slipped in this event.

  • patch_slip (array-like) – Slip magnitude (m) for each entry in patch_numbers.

  • patch_time (array-like) – Absolute rupture time (s) for each entry in patch_numbers.

  • fault_model (RsqSimMultiFault) – Fault model providing patch and fault dictionaries.

  • mask (list of bool) – Boolean mask selecting the patches to retain from patch_numbers.

  • event_id (int, optional) – Unique catalogue event ID.

Returns:

Fully populated event.

Return type:

RsqSimEvent

sub_events_by_fault(fault_model, min_slip=0.1)[source]

Split the event into per-fault sub-events.

For each fault involved in the event, a new RsqSimEvent is created containing only the patches on that fault. Faults where no patch has slip >= min_slip are skipped.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch and fault dictionaries.

  • min_slip (float, optional) – Minimum slip threshold (m) for a patch to be included in a sub-event. Defaults to 0.1 m.

Returns:

One event per fault that contributed qualifying slip.

Return type:

list of RsqSimEvent

find_mean_slip()[source]

Compute and cache the mean slip across all ruptured patches.

Sets _mean_slip (accessed via mean_slip) to the arithmetic mean of patch_slip. Does nothing if no patches are loaded.

find_mean_strike()[source]

Compute and cache the arithmetic mean strike (0–360°).

Sets _mean_strike (accessed via mean_strike). Does nothing if no patches are loaded.

find_mean_strike_180()[source]

Compute and cache the mean strike folded into the range 0–180°.

Strikes in 180–360° are reduced by 180° before averaging, so that conjugate orientations are treated as equivalent. Sets _mean_strike_180 (accessed via mean_strike_180). Does nothing if no patches are loaded.

find_mean_dip()[source]

Compute and cache the arithmetic mean dip (degrees).

Sets _mean_dip (accessed via mean_dip). Does nothing if no patches are loaded.

find_mean_rake()[source]

Compute and cache the arithmetic mean rake (degrees).

Sets _mean_rake (accessed via mean_rake). Does nothing if no patches are loaded.

find_length(min_slip_percentile=None)[source]

Estimate and cache the rupture length along fault surface traces.

For each fault, the range of distances of patch centroids projected onto the fault trace is summed to give a cumulative rupture length. Sets length.

Parameters:

min_slip_percentile (float or None, optional) – Unused; reserved for future filtering by slip percentile.

make_fault_moment_dict(fault_model, mu=30000000000.0, by_cfm_names=True, min_m0=0.0)[source]

Build a dictionary of seismic moment released on each fault.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch and fault dictionaries.

  • mu (float, optional) – Shear modulus in Pa. Defaults to 3×10¹⁰ Pa (30 GPa).

  • by_cfm_names (bool, optional) – If True (default), merge sub-segments by stripping trailing digits and hyphens to recover CFM fault names. If False, use the raw fault segment names.

  • min_m0 (float, optional) – Minimum moment (N·m) for a fault to be included. Defaults to 0 (include all faults).

Returns:

Mapping of fault name (str) to scalar moment (float, N·m).

Return type:

dict

make_moment_prop_dict(fault_model, mu=30000000000.0, by_cfm_names=True)[source]

Build a dictionary of each fault’s fractional contribution to total moment.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch and fault dictionaries.

  • mu (float, optional) – Shear modulus in Pa. Defaults to 3×10¹⁰ Pa.

  • by_cfm_names (bool, optional) – Passed to make_fault_moment_dict(); if True (default), merge sub-segments to CFM fault names.

Returns:

Mapping of fault name (str) to proportion of total event moment (float, 0–1), sorted in descending order.

Return type:

dict

plot_slip_2d(subduction_cmap='plasma', crustal_cmap='viridis', show=True, extra_sub_list=None, write=None, subplots=None, global_max_sub_slip=0, global_max_slip=0, figsize=(6.4, 4.8), hillshading_intensity=0.0, bounds=None, plot_rivers=True, plot_lakes=True, plot_highways=True, plot_boundaries=False, create_background=False, coast_only=True, hillshade_cmap=cm.terrain, plot_log_scale=False, log_cmap='magma', log_min=1.0, log_max=100.0, plot_traces=True, trace_colour='pink', land_color='antiquewhite', min_slip_percentile=None, min_slip_value=None, plot_zeros=True, wgs=False, title=None, plot_edge_label=True, plot_cbars=True, alpha=1.0, coast_on_top=False)[source]

Plot a 2-D map of the slip distribution for this event.

Subduction-interface faults and crustal faults are coloured with separate colourmaps. Optionally overlays a coastal/background map and saves to file.

Parameters:
  • subduction_cmap (str, optional) – Matplotlib colourmap name for subduction-interface patches. Defaults to "plasma".

  • crustal_cmap (str, optional) – Matplotlib colourmap name for crustal patches. Defaults to "viridis".

  • show (bool, optional) – If True (default), call plt.show() after plotting.

  • extra_sub_list (list of str, optional) – Additional fault names to treat as subduction interface, supplementing the built-in bruce_subduction list.

  • write (str or None, optional) – File path to save the figure (e.g. "event_001.png"). If None, the figure is not saved.

  • subplots (tuple or str or None, optional) – Existing (fig, ax) tuple or path to a pickled figure to plot on top of. If None, a new figure is created.

  • global_max_sub_slip (float, optional) – Fixed colourbar maximum for subduction slip. If 0, the per-event maximum is used.

  • global_max_slip (float, optional) – Fixed colourbar maximum for crustal slip. If 0, the per-event maximum is used.

  • figsize (tuple of float, optional) – Figure size in inches (width, height).

  • hillshading_intensity (float, optional) – Hillshading intensity passed to plot_background.

  • bounds (tuple or None, optional) – Map extent (x_min, y_min, x_max, y_max) in NZTM. Defaults to the event bounding box.

  • plot_rivers (bool, optional) – Toggle background map layers.

  • plot_lakes (bool, optional) – Toggle background map layers.

  • plot_highways (bool, optional) – Toggle background map layers.

  • plot_boundaries (bool, optional) – Toggle background map layers.

  • create_background (bool, optional) – If True, render a full background map.

  • coast_only (bool, optional) – If True (default), only render the coastline as background.

  • hillshade_cmap (LinearSegmentedColormap, optional) – Colourmap for hillshading.

  • plot_log_scale (bool, optional) – If True, use logarithmic colour scaling.

  • log_cmap (str, optional) – Colourmap for log-scale plots. Defaults to "magma".

  • log_min (float, optional) – Colour scale limits for log-scale plots.

  • log_max (float, optional) – Colour scale limits for log-scale plots.

  • plot_traces (bool, optional) – If True, plot fault surface traces.

  • trace_colour (str, optional) – Colour for fault traces.

  • land_color (str, optional) – Background land colour.

  • min_slip_percentile (float or None, optional) – If set, only show patches above this slip percentile.

  • min_slip_value (float or None, optional) – If set, only show patches with slip >= this value (m).

  • plot_zeros (bool, optional) – If True (default), plot zero-slip patches.

  • wgs (bool, optional) – If True, use WGS84 coordinates rather than NZTM.

  • title (str or None, optional) – Figure super-title.

  • plot_edge_label (bool, optional) – If True, show axis edge labels on the background map.

  • plot_cbars (bool, optional) – If True (default), add colourbars to the figure.

  • alpha (float, optional) – Patch transparency (0–1). Defaults to 1.

  • coast_on_top (bool, optional) – If True, redraw the coastline on top of the slip patches.

Returns:

Matplotlib PolyCollection objects for each fault plotted.

Return type:

list

plot_slip_evolution(subduction_cmap='plasma', crustal_cmap='viridis', show=True, step_size=1, write=None, fps=20, file_format='gif', figsize=(6.4, 4.8), extra_sub_list=None)[source]

Animate the temporal evolution of slip propagation for this event.

Produces a FuncAnimation that steps through rupture time, progressively revealing each patch as it slips. Can be saved as a GIF or video file.

Parameters:
  • subduction_cmap (str, optional) – Matplotlib colourmap for subduction-interface patches. Defaults to "plasma".

  • crustal_cmap (str, optional) – Matplotlib colourmap for crustal patches. Defaults to "viridis".

  • show (bool, optional) – If True (default), display the animation interactively.

  • step_size (int, optional) – Time step (s) between animation frames. Defaults to 1.

  • write (str or None, optional) – Output file path prefix (without extension). If None the animation is not saved.

  • fps (int, optional) – Frames per second for the saved animation. Defaults to 20.

  • file_format (str, optional) – Output format: "gif", "mov", "avi", or "mp4". Defaults to "gif".

  • figsize (tuple of float, optional) – Figure size in inches (width, height).

  • extra_sub_list (list of str, optional) – Additional fault names to treat as subduction interface.

find_surface_faults(fault_model, min_slip=0.1, method='vertex', n_patches=1, max_depth=-1000.0, faults2ignore='hikurangi')[source]

Identify fault segments with surface-rupturing patches.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing the patch dictionary.

  • min_slip (float, optional) – Minimum slip (m) for a patch to count as ruptured. Defaults to 0.1 m.

  • method (str, optional) – Depth criterion: "vertex" uses the shallowest vertex of a patch; "centroid" uses the patch centroid depth. Defaults to "vertex".

  • n_patches (int, optional) – Minimum number of qualifying surface patches required for a fault to be included in the output. Defaults to 1.

  • max_depth (float, optional) – Depth threshold (m, negative downward). Patches shallower than this are considered surface-rupturing. Defaults to -1000 m.

  • faults2ignore (list of str or str, optional) – Fault name(s) to exclude from the search (e.g. subduction interface). Defaults to "hikurangi".

Returns:

Names of fault segments that have surface rupture.

Return type:

list of str

split_by_fault(fault_model, min_slip=0.1, min_patches=1)[source]

Create per-fault sub-event objects for this event.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model used to resolve patch IDs and fault names.

  • min_slip (float, optional) – Minimum moment (N·m) for a fault to be included (passed as min_m0 to make_fault_moment_dict()). Defaults to 0.1.

  • min_patches (int, optional) – Unused; retained for API compatibility.

Returns:

Mapping of fault name (str) to its corresponding RsqSimEvent sub-event.

Return type:

dict

slip_dist_array(include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False)[source]

Build a dense array of the slip distribution across all patches.

Each row corresponds to one triangular patch and contains the three vertex coordinates followed by slip, rake, and rupture time.

Parameters:
  • include_zeros (bool, optional) – If True (default), include zero-slip patches (patches on involved faults that did not slip in this event).

  • min_slip_percentile (float or None, optional) – If provided (and min_slip_value is None), patches with slip below this percentile are excluded or zeroed.

  • min_slip_value (float or None, optional) – If provided, patches with slip below this value (m) are excluded or zeroed.

  • nztm_to_lonlat (bool, optional) – If True, transform vertex coordinates from NZTM to WGS84 longitude/latitude before including them.

Returns:

Columns: [x1,y1,z1, x2,y2,z2, x3,y3,z3, slip_m, rake_deg, time_s].

Return type:

numpy.ndarray of shape (n_patches, 12)

slip_dist_bounds(include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False)[source]

Compute the bounding box of the slip distribution array.

Parameters:
Returns:

(x_min, y_min, x_max, y_max) extents of the slip distribution patch vertices.

Return type:

tuple of float

slip_dist_to_mesh(include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False)[source]

Convert the slip distribution to a meshio Mesh with cell data.

Parameters:
Returns:

Triangulated mesh with cell_data containing "slip", "rake", and "time" arrays.

Return type:

meshio.Mesh

slip_dist_to_vtk(vtk_file, include_zeros=True, min_slip_percentile=None, min_slip_value=None)[source]

Write the slip distribution to a VTK file.

Parameters:
slip_dist_to_obj(obj_file, include_zeros=True, min_slip_percentile=None, min_slip_value=None)[source]

Write the slip distribution to a Wavefront OBJ file.

Parameters:
slip_dist_to_txt(txt_file, include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False)[source]

Write the slip distribution to a space-delimited text file.

Each row contains the three vertex coordinates of a patch followed by slip, rake, and rupture time.

Parameters:
  • txt_file (str or path-like) – Output text file path.

  • include_zeros (bool, optional) – Passed to slip_dist_array().

  • min_slip_percentile (float or None, optional) – Passed to slip_dist_array().

  • min_slip_value (float or None, optional) – Passed to slip_dist_array().

  • nztm_to_lonlat (bool, optional) – If True, vertex coordinates are in WGS84 lon/lat and the header is adjusted accordingly.

slip_dist_to_gdf(gdf_file, include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False, crs='2193')[source]

Build a GeoDataFrame of the slip distribution.

Creates a GeoDataFrame where each row is a triangular patch polygon with slip, rake, and time attributes.

Parameters:
  • gdf_file (str or None) – Unused; retained for API compatibility.

  • include_zeros (bool, optional) – Passed to slip_dist_array().

  • min_slip_percentile (float or None, optional) – Passed to slip_dist_array().

  • min_slip_value (float or None, optional) – Passed to slip_dist_array().

  • nztm_to_lonlat (bool, optional) – If True, reproject the result from NZTM to WGS84.

  • crs (str, optional) – CRS for the output GeoDataFrame. Defaults to "2193" (NZTM / EPSG:2193).

Returns:

GeoDataFrame with columns ["slip", "rake", "time"] and a geometry column of triangular patch polygons.

Return type:

geopandas.GeoDataFrame

slip_dist_to_geojson(geojson_file, include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False)[source]

Write the slip distribution to a GeoJSON file.

Parameters:
  • geojson_file (str) – Output GeoJSON file path.

  • include_zeros (bool, optional) – Passed to slip_dist_to_gdf().

  • min_slip_percentile (float or None, optional) – Passed to slip_dist_to_gdf().

  • min_slip_value (float or None, optional) – Passed to slip_dist_to_gdf().

  • nztm_to_lonlat (bool, optional) – If True, reproject patches to WGS84 before writing.

slip_dist_to_shapefile(shapefile_file, include_zeros=True, min_slip_percentile=None, min_slip_value=None, nztm_to_lonlat=False)[source]

Write the slip distribution to an ESRI Shapefile.

Parameters:
  • shapefile_file (str) – Output shapefile path (without extension).

  • include_zeros (bool, optional) – Passed to slip_dist_to_gdf().

  • min_slip_percentile (float or None, optional) – Passed to slip_dist_to_gdf().

  • min_slip_value (float or None, optional) – Passed to slip_dist_to_gdf().

  • nztm_to_lonlat (bool, optional) – If True, reproject patches to WGS84 before writing.

slip_dist_to_gmt(fault_model, gmt_prefix, min_slip_value=None, nztm_to_lonlat=False, subduction_names=('hikkerm', 'puysegur'))[source]

Write the slip distribution to GMT multi-segment text files.

Writes two files: <gmt_prefix>_crustal.gmt for crustal faults and <gmt_prefix>_subduction.gmt for subduction-interface faults. Each patch is written as a >-Z<slip> segment header followed by vertex coordinates.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing the patch dictionary.

  • gmt_prefix (str) – Prefix for the output GMT file names.

  • min_slip_value (float or None, optional) – If set, only patches with slip >= this value (m) are written.

  • nztm_to_lonlat (bool, optional) – Unused; reserved for future coordinate transformation.

  • subduction_names (iterable of str, optional) – Fault names to classify as subduction interface. Defaults to ("hikkerm", "puysegur").

discretize_tiles(tile_list, probability, rake)[source]

Select tiles that overlap the rupture exterior by at least 50 %.

Parameters:
  • tile_list (list of Polygon) – Candidate rectangular tiles in NZTM coordinates.

  • probability (float) – Unused; retained for API compatibility with discretize_openquake().

  • rake (float) – Unused; retained for API compatibility.

Returns:

GeoSeries (EPSG:2193) of tiles where the intersection with the rupture exterior covers at least half the tile area.

Return type:

geopandas.GeoSeries

discretize_openquake(tile_list, probability, rake)[source]

Build an OpenQuake multi-planar rupture from overlapping tiles.

Selects tiles from tile_list that overlap the rupture exterior by at least 50 %, then packages them as an OpenQuakeMultiSquareRupture.

Parameters:
Returns:

Rupture object if any tiles overlap; None otherwise.

Return type:

OpenQuakeMultiSquareRupture or None

discretize_openquake_ktree(fault_model, quads_dict, probability, subduction_names=('hikkerm', 'puysegur'), min_moment=1e+18, min_slip=0.1, tile_size=5000.0, write_mesh=False, write_geojson=False, xml_dir=None, threshold=0.5)[source]

Build OpenQuake XML ruptures using a KD-tree tile assignment.

Uses slip_dist_quads_ktree() to find ruptured quadrilateral tiles per fault, then writes separate XML files for crustal and subduction components.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch and fault dictionaries.

  • quads_dict (dict) – Mapping of fault segment name to a (n, 4, 3) array of pre-computed quadrilateral tile corner coordinates.

  • probability (float) – Probability of occurrence for the OpenQuake rupture.

  • subduction_names (iterable of str, optional) – Fault names classified as subduction interface.

  • min_moment (float, optional) – Minimum moment (N·m) for a fault to be included.

  • min_slip (float, optional) – Minimum slip (m) for a patch to count as ruptured.

  • tile_size (float, optional) – Tile size hint (m); unused here but passed through for context.

  • write_mesh (bool, optional) – If True, also write VTK mesh files for each component.

  • write_geojson (bool, optional) – If True, also write GeoJSON files for each component.

  • xml_dir (str or None, optional) – Directory in which to write output files. Defaults to the current working directory.

  • threshold (float, optional) – Fraction of patches in a tile that must be ruptured for the tile to be included. Passed to slip_dist_quads_ktree().

event_to_json(fault_model, path2cfm, catalogue_version='v1', xml_dir='OQ-events', wgs84=False, subd_tile_size=15000.0, tile_size=5000.0, tectonic_region='NZ')[source]

Export event slip distribution to a ShakeMap-compatible GeoJSON file.

Discretises the rupturing faults into rectangular tiles, selects those near slipping patches, converts coordinates to WGS84, and writes a GeoJSON file with metadata to <xml_dir>/event_<id>/<id>.json.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model used to merge and discretise fault segments.

  • path2cfm (str) – Path to the Community Fault Model directory (required for catalogue_version="v2" name lookups).

  • catalogue_version (str, optional) – "v1" or "v2". Affects which fault name dictionary is used. Defaults to "v1".

  • xml_dir (str, optional) – Output directory. Defaults to "OQ-events".

  • wgs84 (bool, optional) – If True, hypocentre coordinates are already in WGS84 (lon/lat). If False (default), NZTM coordinates are transformed to WGS84.

  • subd_tile_size (float, optional) – Tile size (m) used when discretising subduction faults. Defaults to 15000 m.

  • tile_size (float, optional) – Tile size (m) used when discretising crustal faults. Defaults to 5000 m.

  • tectonic_region (str, optional) – Tectonic region tag written into the metadata. Defaults to "NZ".

event_to_OQ_xml(fault_model, path2cfm, catalogue_version='v2', xml_dir='OQ_events', subd_tile_size=15000.0, tile_size=5000.0, probability=0.9, tectonic_region='NZ', min_mag=6.0, hypocentre=None, nztm2wgs=True)[source]

Export the event as an OpenQuake multi-planar rupture XML file.

Discretises faults exceeding min_mag into rectangular tiles, selects tiles near ruptured patches, and writes an OpenQuake NRML XML file to <xml_dir>/<event_id>/event_<event_id>.xml.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model for segment merging and discretisation.

  • path2cfm (str) – Path to the Community Fault Model directory.

  • catalogue_version (str, optional) – "v1" or "v2". Defaults to "v2".

  • xml_dir (str, optional) – Output directory. Defaults to "OQ_events".

  • subd_tile_size (float, optional) – Tile size (m) for subduction faults. Defaults to 15000 m.

  • tile_size (float, optional) – Tile size (m) for crustal faults. Defaults to 5000 m.

  • probability (float, optional) – Probability of occurrence for the rupture. Defaults to 0.9.

  • tectonic_region (str, optional) – Tectonic region tag. Defaults to "NZ".

  • min_mag (float, optional) – Only include faults whose moment corresponds to at least this magnitude. Defaults to 6.0.

  • hypocentre (list or None, optional) – [x, y, z] hypocentre override. If None, uses the event hypocentre.

  • nztm2wgs (bool, optional) – If True (default), convert NZTM coordinates to WGS84 before writing.

slip_dist_quads_ktree(fault_model, quads_dict, min_moment=1e+18, min_slip=0.0, threshold_for_inclusion=0.5, slip_per_quad=False)[source]

Find ruptured quadrilateral tiles for each fault using a KD-tree.

For each fault segment that exceeds min_moment, uses a KD-tree on pre-computed quad centroids to assign triangular patches to quads, then returns quads where the fraction of ruptured patches exceeds threshold_for_inclusion.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch dictionaries and patch centres.

  • quads_dict (dict) – Mapping of fault segment name to a (n, 4, 3) array of quad corner coordinates (NZTM, m).

  • min_moment (float, optional) – Minimum moment (N·m) for a fault to be considered. Defaults to 1×10¹⁸ N·m.

  • min_slip (float, optional) – Minimum patch slip (m) for a patch to count as ruptured. Defaults to 0.

  • threshold_for_inclusion (float, optional) – Fraction of a quad’s assigned patches that must be ruptured for the quad to be included. Defaults to 0.5.

  • slip_per_quad (bool, optional) – If True, compute and store average slip per quad (not yet returned in output).

Returns:

Mapping of fault segment name to a (n_ruptured, 4, 3) array of ruptured quad corner coordinates.

Return type:

dict

to_oq_points(fault_model)[source]

Convert event to OpenQuake point-source representation.

Not yet implemented.

Parameters:

fault_model (RsqSimMultiFault) – Fault model (reserved for future use).

get_crustal_component(fault_model, crustal_names, min_moment=1e+18, min_slip=0.0)[source]

Compute moment, mean rake, and hypocentre for crustal fault components.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch dictionaries and patch moments.

  • crustal_names (list of str) – Names of fault segments to treat as crustal.

  • min_moment (float, optional) – Minimum moment (N·m) for a fault to be considered. Defaults to 1×10¹⁸ N·m.

  • min_slip (float, optional) – Minimum patch slip (m) to count as ruptured. Defaults to 0.

Returns:

(mw, mean_rake, hypocentre_xyz) for the crustal component, or None if no crustal faults qualify.

Return type:

tuple of (float, float, numpy.ndarray) or None

get_subduction_component(fault_model, subduction_names, min_moment=1e+18, min_slip=0.0)[source]

Compute moment, mean rake, and hypocentre for subduction fault components.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model providing patch dictionaries and patch moments.

  • subduction_names (list of str) – Names of fault segments to treat as subduction interface.

  • min_moment (float, optional) – Minimum moment (N·m) for a fault to be considered. Defaults to 1×10¹⁸ N·m.

  • min_slip (float, optional) – Minimum patch slip (m) to count as ruptured. Defaults to 0.

Returns:

(mw, mean_rake, hypocentre_xyz) for the subduction component, or None if no subduction faults qualify.

Return type:

tuple of (float, float, numpy.ndarray) or None

slip_dist_to_quads(fault_model, path2cfm, catalogue_version='v2', vtk_dir='fault_vtks', subd_tile_size=15000.0, tile_size=5000.0)[source]

Write per-fault slip distributions as quadrilateral VTK meshes.

Discretises each involved fault into rectangular tiles, assigns the average slip of the nearest triangular patches to each tile, and writes a VTK file per fault to <vtk_dir>/event_<id>_quad/<fault_name>_<event_id>.vtk.

Parameters:
  • fault_model (RsqSimMultiFault) – Fault model for segment merging and discretisation.

  • path2cfm (str) – Path to the Community Fault Model directory.

  • catalogue_version (str, optional) – "v1" or "v2". Defaults to "v2".

  • vtk_dir (str, optional) – Output directory for VTK files. Defaults to "fault_vtks".

  • subd_tile_size (float, optional) – Tile size (m) for subduction faults. Defaults to 15000 m.

  • tile_size (float, optional) – Tile size (m) for crustal faults. Defaults to 5000 m.

class rsqsim_api.catalogue.event.OpenQuakeMultiSquareRupture(tile_list, probability, magnitude, rake, hypocentre, event_id, name='Subduction earthquake', tectonic_region='subduction', nztm2wgs=True)[source]

Multi-planar rupture representation for the OpenQuake engine.

Packages a list of rectangular tile polygons (each representing a fault-surface patch) together with event metadata into the structure needed to write OpenQuake NRML multiPlanesRupture XML.

Parameters:
  • tile_list (list[shapely.geometry.Polygon])

  • probability (float)

  • magnitude (float)

  • rake (float)

  • hypocentre (numpy.ndarray)

  • event_id (int)

  • name (str)

  • tectonic_region (str)

  • nztm2wgs (bool)

patches[source]

Individual rectangular patches converted from input polygons.

Type:

list of OpenQuakeRectangularPatch

prob[source]

Probability of occurrence.

Type:

float

magnitude[source]

Moment magnitude.

Type:

float

rake[source]

Mean rake (degrees).

Type:

float

hypocentre[source]

Hypocentre coordinates [x, y, z] in the input CRS.

Type:

numpy.ndarray

inv_prob[source]

Complementary probability (1 - prob).

Type:

float

hyp_depth[source]

Hypocentre depth in km (positive downward).

Type:

float

hyp_lon, hyp_lat

Hypocentre longitude and latitude in WGS84 degrees.

Type:

float

event_id[source]

Catalogue event identifier.

Type:

int

name[source]

Descriptive rupture name.

Type:

str

tectonic_region[source]

Tectonic region tag for OpenQuake.

Type:

str

patches[source]
prob[source]
magnitude[source]
rake[source]
hypocentre[source]
inv_prob[source]
hyp_depth[source]
event_id[source]
name = 'Subduction earthquake'[source]
tectonic_region = 'subduction'[source]
to_oq_xml(write=None)[source]

Build (and optionally write) an OpenQuake NRML rupture XML element.

Constructs a <multiPlanesRupture> element containing magnitude, rake, hypocenter, and one <planarSurface> per patch.

Parameters:

write (str or None, optional) – Output file path. If the path does not end with ".xml" the extension is appended automatically. If None, the XML is not written to disk.

Returns:

The root <nrml> element (regardless of whether the file was written).

Return type:

xml.etree.ElementTree.Element