openquake.hmtk.seismicity package#
Subpackages#
- openquake.hmtk.seismicity.completeness package
- openquake.hmtk.seismicity.declusterer package
- openquake.hmtk.seismicity.max_magnitude package
- Submodules
- openquake.hmtk.seismicity.max_magnitude.base module
- openquake.hmtk.seismicity.max_magnitude.cumulative_moment_release module
- openquake.hmtk.seismicity.max_magnitude.kijko_nonparametric_gaussian module
- openquake.hmtk.seismicity.max_magnitude.kijko_sellevol_bayes module
- openquake.hmtk.seismicity.max_magnitude.kijko_sellevol_fixed_b module
- Module contents
- openquake.hmtk.seismicity.occurrence package
- Submodules
- openquake.hmtk.seismicity.occurrence.aki_maximum_likelihood module
- openquake.hmtk.seismicity.occurrence.b_maximum_likelihood module
- openquake.hmtk.seismicity.occurrence.base module
- openquake.hmtk.seismicity.occurrence.kijko_smit module
- openquake.hmtk.seismicity.occurrence.penalized_mle module
- openquake.hmtk.seismicity.occurrence.utils module
- openquake.hmtk.seismicity.occurrence.weichert module
- Module contents
- openquake.hmtk.seismicity.smoothing package
Submodules#
openquake.hmtk.seismicity.catalogue module#
Prototype of a ‘Catalogue’ class
- class openquake.hmtk.seismicity.catalogue.Catalogue[source]#
Bases:
object
General Catalogue Class
- FLOAT_ATTRIBUTE_LIST = ['second', 'timeError', 'longitude', 'latitude', 'SemiMajor90', 'SemiMinor90', 'ErrorStrike', 'depth', 'depthError', 'magnitude', 'sigmaMagnitude']#
- INT_ATTRIBUTE_LIST = ['year', 'month', 'day', 'hour', 'minute', 'flag']#
- SORTED_ATTRIBUTE_LIST = ['eventID', 'Agency', 'year', 'month', 'day', 'hour', 'minute', 'second', 'timeError', 'longitude', 'latitude', 'SemiMajor90', 'SemiMinor90', 'ErrorStrike', 'depth', 'depthError', 'magnitude', 'sigmaMagnitude', 'magnitudeType']#
- STRING_ATTRIBUTE_LIST = ['eventID', 'Agency', 'magnitudeType', 'comment']#
- TOTAL_ATTRIBUTE_LIST = ['magnitude', 'year', 'ErrorStrike', 'sigmaMagnitude', 'comment', 'day', 'second', 'depthError', 'eventID', 'Agency', 'month', 'flag', 'SemiMinor90', 'longitude', 'latitude', 'magnitudeType', 'timeError', 'SemiMajor90', 'hour', 'depth', 'minute']#
- catalogue_mt_filter(mt_table, flag=None)[source]#
Filter the catalogue using a magnitude-time table. The table has two columns and n-rows.
- Parameters:
mt_table (nump.ndarray) – Magnitude time table with n-rows where column 1 is year and column 2 is magnitude
- concatenate(catalogue)[source]#
This method attaches one catalogue to the current one
- Parameters:
catalogue – An instance of
htmk.seismicity.catalogue.Catalogue
- get_bounding_box()[source]#
Returns the bounding box of the catalogue
- Returns:
(West, East, South, North)
- get_depth_distribution(depth_bins, normalisation=False, bootstrap=None)[source]#
Gets the depth distribution of the earthquake catalogue to return a single histogram. Depths may be normalised. If uncertainties are found in the catalogue the distrbution may be bootstrap sampled
- Parameters:
depth_bins (numpy.ndarray) – getBin edges for the depths
normalisation (bool) – Choose to normalise the results such that the total contributions sum to 1.0 (True) or not (False)
bootstrap (int) – Number of bootstrap samples
- Returns:
Histogram of depth values
- get_depth_pmf(depth_bins, default_depth=5.0, bootstrap=None)[source]#
Returns the depth distribution of the catalogue as a probability mass function
- get_magnitude_depth_distribution(magnitude_bins, depth_bins, normalisation=False, bootstrap=None)[source]#
Returns a 2-D magnitude-depth histogram for the catalogue
- Parameters:
magnitude_bins (numpy.ndarray) – Bin edges for the magnitudes
depth_bins (numpy.ndarray) – Bin edges for the depths
normalisation (bool) – Choose to normalise the results such that the total contributions sum to 1.0 (True) or not (False)
bootstrap (int) – Number of bootstrap samples
- Returns:
2D histogram of events in magnitude-depth bins
- get_magnitude_time_distribution(magnitude_bins, time_bins, normalisation=False, bootstrap=None)[source]#
Returns a 2-D histogram indicating the number of earthquakes in a set of time-magnitude bins. Time is in decimal years!
- Parameters:
magnitude_bins (numpy.ndarray) – Bin edges for the magnitudes
time_bins (numpy.ndarray) – Bin edges for the times
normalisation (bool) – Choose to normalise the results such that the total contributions sum to 1.0 (True) or not (False)
bootstrap (int) – Number of bootstrap samples
- Returns:
2D histogram of events in magnitude-year bins
- get_observed_mmax_sigma(default=None)[source]#
- Returns:
the sigma for the maximum observed magnitude
- load_from_array(keys, data_array)[source]#
This loads the data contained in an array into the catalogue object
- Parameters:
keys – A list of keys explaining the content of the columns in the array
- load_to_array(keys)[source]#
This loads the data contained in the catalogue into a numpy array. The method works only for float data
- Parameters:
keys – A list of keys to be uploaded into the array
- purge_catalogue(flag_vector)[source]#
Purges present catalogue with invalid events defined by flag_vector
- Parameters:
flag_vector (numpy.ndarray) – Boolean vector showing if events are selected (True) or not (False)
- select_catalogue_events(id0)[source]#
Orders the events in the catalogue according to an indexing vector.
- Parameters:
id0 (np.ndarray) – Pointer array indicating the locations of selected events
- update_end_year()[source]#
NOTE: To be called only when the catalogue is loaded (not when it is modified by declustering or completeness-based filtering)
- update_start_year()[source]#
NOTE: To be called only when the catalogue is loaded (not when it is modified by declustering or completeness-based filtering)
- write_catalogue(output_file, key_list=['eventID', 'Agency', 'year', 'month', 'day', 'hour', 'minute', 'second', 'timeError', 'longitude', 'latitude', 'SemiMajor90', 'SemiMinor90', 'ErrorStrike', 'depth', 'depthError', 'magnitude', 'sigmaMagnitude', 'magnitudeType'])[source]#
Writes the catalogue to file using HTMK format (CSV).
- Parameters:
output_file – Name of the output file
key_list – Optional list of attribute keys to be exported
openquake.hmtk.seismicity.gcmt_catalogue module#
Implements sets of classes for mapping components of the focal mechanism
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTCatalogue(start_year=None, end_year=None)[source]#
Bases:
Catalogue
Class to hold a catalogue of moment tensors
- FLOAT_ATTRIBUTE_LIST = ['second', 'timeError', 'longitude', 'latitude', 'SemiMajor90', 'SemiMinor90', 'ErrorStrike', 'depth', 'depthError', 'magnitude', 'sigmaMagnitude', 'moment', 'strike1', 'rake1', 'dip1', 'strike2', 'rake2', 'dip2', 'eigenvalue_b', 'azimuth_b', 'plunge_b', 'eigenvalue_p', 'azimuth_p', 'plunge_p', 'eigenvalue_t', 'azimuth_t', 'plunge_t', 'f_clvd', 'e_rel']#
- INT_ATTRIBUTE_LIST = ['eventID', 'year', 'month', 'day', 'hour', 'minute', 'flag']#
- STRING_ATTRIBUTE_LIST = ['Agency', 'magnitudeType', 'comment', 'centroidID']#
- TOTAL_ATTRIBUTE_LIST = ['azimuth_b', 'depthError', 'eventID', 'centroidID', 'Agency', 'plunge_p', 'month', 'flag', 'longitude', 'e_rel', 'magnitudeType', 'eigenvalue_b', 'timeError', 'dip2', 'azimuth_p', 'SemiMajor90', 'strike1', 'minute', 'magnitude', 'year', 'ErrorStrike', 'sigmaMagnitude', 'comment', 'day', 'second', 'f_clvd', 'rake1', 'plunge_t', 'eigenvalue_p', 'SemiMinor90', 'moment', 'latitude', 'plunge_b', 'rake2', 'azimuth_t', 'hour', 'eigenvalue_t', 'dip1', 'depth', 'strike2']#
- gcmt_to_simple_array(centroid_location=True)[source]#
Converts the GCMT catalogue to a simple array of [ID, year, month, day, hour, minute, second, long., lat., depth, Mw, strike1, dip1, rake1, strike2, dip2, rake2, b-plunge, b-azimuth, b-eigenvalue, p-plunge, p-azimuth, p-eigenvalue, t-plunge, t-azimuth, t-eigenvalue, moment, f_clvd, erel]
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTCentroid(reference_date, reference_time)[source]#
Bases:
object
Representation of a GCMT centroid
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTEvent[source]#
Bases:
object
Class to represent full GCMT moment tensor in ndk format
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTHypocentre[source]#
Bases:
object
Simple representation of the hypocentre
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTMomentTensor(reference_frame=None)[source]#
Bases:
object
Class to represent a moment tensor and its associated methods
- eigendecompose(normalise=False)[source]#
Performs and eigendecomposition of the tensor and orders into descending eigenvalues
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTNodalPlanes[source]#
Bases:
object
Class to represent the nodal plane distribution of the tensor Each nodal plane is represented as a dictionary of the form: {‘strike’:, ‘dip’:, ‘rake’:}
- class openquake.hmtk.seismicity.gcmt_catalogue.GCMTPrincipalAxes[source]#
Bases:
object
Class to represent the eigensystem of the tensor in terms of T-, B- and P- plunge and azimuth
openquake.hmtk.seismicity.gcmt_utils module#
Set of moment tensor utility functions
- openquake.hmtk.seismicity.gcmt_utils.eigendecompose(tensor, normalise=False)[source]#
Performs and eigendecomposition of the tensor and orders into descending eigenvalues
- openquake.hmtk.seismicity.gcmt_utils.get_azimuth_plunge(vect, degrees=True)[source]#
For a given vector in USE format, retrieve the azimuth and plunge
- openquake.hmtk.seismicity.gcmt_utils.moment_magnitude_scalar(moment)[source]#
Uses Hanks & Kanamori formula for calculating moment magnitude from a scalar moment (Nm)
- openquake.hmtk.seismicity.gcmt_utils.ned_to_use(tensor)[source]#
Converts a tensor in NED coordinate sytem to USE
- openquake.hmtk.seismicity.gcmt_utils.normalise_tensor(tensor)[source]#
Normalise the tensor by dividing it by its norm, defined such that np.sqrt(X:X)
- openquake.hmtk.seismicity.gcmt_utils.tensor_components_to_ned(mrr, mtt, mpp, mrt, mrp, mtp)[source]#
Converts components to North, East, Down definition:
NED = [[mtt, -mtp, mrt], [-mtp, mpp, -mrp], [mrt, -mtp, mrr]]
- openquake.hmtk.seismicity.gcmt_utils.tensor_components_to_use(mrr, mtt, mpp, mrt, mrp, mtp)[source]#
Converts components to Up, South, East definition:
USE = [[mrr, mrt, mrp], [mtt, mtt, mtp], [mrp, mtp, mpp]]
- openquake.hmtk.seismicity.gcmt_utils.tensor_to_6component(tensor, frame='USE')[source]#
Returns a tensor to six component vector [Mrr, Mtt, Mpp, Mrt, Mrp, Mtp]
- openquake.hmtk.seismicity.gcmt_utils.unique_euler(alpha, beta, gamma)[source]#
Uniquify euler angle triplet. Put euler angles into ranges compatible with (dip,strike,-rake) in seismology: alpha (dip) : [0, pi/2] beta (strike) : [0, 2*pi) gamma (-rake) : [-pi, pi) If alpha is near to zero, beta is replaced by beta+gamma and gamma is set to zero, to prevent that additional ambiguity.
If alpha is near to pi/2, beta is put into the range [0,pi).
openquake.hmtk.seismicity.selector module#
Class to implement set of functionalities for selecting events from and earthquake catalogue
- class openquake.hmtk.seismicity.selector.CatalogueSelector(master_catalogue, create_copy=True)[source]#
Bases:
object
Class to implement methods for selecting subsets of the catalogue according to various attribute criteria.
- Attr catalogue:
The catalogue to which the selection is applied as instance of openquake.hmtk.seismicity.catalogue.Catalogue
- Attr create_copy:
Boolean to indicate whether to create copy of the original catalogue before selecting {default = True}
- cartesian_square_centred_on_point(point, distance, **kwargs)[source]#
Select earthquakes from within a square centered on a point
- Parameters:
point – Centre point as instance of nhlib.geo.point.Point class
distance – Distance (km)
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
class containing only selected events
- circular_distance_from_point(point, distance, **kwargs)[source]#
Select earthquakes within a distance from a Point
- Parameters:
point – Centre point as instance of nhlib.geo.point.Point class
distance (float) – Distance (km)
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
containing only selected events
- create_cluster_set(vcl)[source]#
For a given catalogue and list of cluster IDs this function splits the catalogue into a dictionary containing an individual catalogue of events within each cluster
- Parameters:
vcl (numpy.ndarray) – Cluster ID list
- Returns:
Dictionary of instances of the :class: openquake.hmtk.seismicity.catalogue.Catalogue, where each instance if the catalogue of each cluster
- select_catalogue(valid_id)[source]#
Method to post-process the catalogue based on the selection options
- Parameters:
valid_id (numpy.ndarray) – Boolean vector indicating whether each event is selected (True) or not (False)
- Returns:
Catalogue of selected events as instance of openquake.hmtk.seismicity.catalogue.Catalogue class
- within_bounding_box(limits)[source]#
Selects the earthquakes within a bounding box.
- Parameters:
limits –
- A list or a numpy array with four elements in the following order:
min x (longitude)
min y (latitude)
max x (longitude)
max y (latitude)
- Returns:
Returns a :class:htmk.seismicity.catalogue.Catalogue` instance
- within_depth_range(lower_depth=None, upper_depth=None)[source]#
Selects events within a specified depth range
- Parameters:
lower_depth (float) – Lower depth for consideration
upper_depth (float) – Upper depth for consideration
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
containing only selected events
- within_joyner_boore_distance(surface, distance, **kwargs)[source]#
Select events within a Joyner-Boore distance of a fault
- Parameters:
surface – Fault surface as instance of nhlib.geo.surface.base.SimpleFaultSurface or as instance of nhlib.geo.surface.ComplexFaultSurface
distance (float) – Rupture distance (km)
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
containing only selected events
- within_magnitude_range(lower_mag=None, upper_mag=None)[source]#
- Parameters:
lower_mag (float) – Lower magnitude for consideration
upper_mag (float) – Upper magnitude for consideration
- Returns:
Instance of openquake.hmtk.seismicity.catalogue.Catalogue class containing only selected events
- within_polygon(polygon, distance=None, **kwargs)[source]#
Select earthquakes within polygon
- Parameters:
polygon – Centre point as instance of nhlib.geo.polygon.Polygon class
distance (float) – Buffer distance (km) (can take negative values)
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
containing only selected events
- within_rupture_distance(surface, distance, **kwargs)[source]#
Select events within a rupture distance from a fault surface
- Parameters:
surface – Fault surface as instance of nhlib.geo.surface.base.BaseSurface
distance (float) – Rupture distance (km)
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
containing only selected events
- within_time_period(start_time=None, end_time=None)[source]#
Select earthquakes occurring within a given time period
- Parameters:
start_time – Earliest time (as datetime.datetime object)
end_time – Latest time (as datetime.datetime object)
- Returns:
Instance of
openquake.hmtk.seismicity.catalogue.Catalogue
containing only selected events
openquake.hmtk.seismicity.utils module#
Utility functions for seismicity calculations
- openquake.hmtk.seismicity.utils.area_of_polygon(polygon)[source]#
Returns the area of an OpenQuake polygon in square kilometres
- openquake.hmtk.seismicity.utils.bootstrap_histogram_1D(values, intervals, uncertainties=None, normalisation=False, number_bootstraps=None, boundaries=None, random_seed=42)[source]#
Bootstrap samples a set of vectors
- Parameters:
values (numpy.ndarray) – The data values
intervals (numpy.ndarray) – The bin edges
uncertainties (numpy.ndarray) – The standard deviations of each observation
normalisation (bool) – If True then returns the histogram as a density function
number_bootstraps (int) – Number of bootstraps
boundaries (tuple) – (Lower, Upper) bounds on the data
random_seed – Seed used in the random number generator
returns – 1-D histogram of data
- openquake.hmtk.seismicity.utils.bootstrap_histogram_2D(xvalues, yvalues, xbins, ybins, boundaries=[None, None], xsigma=None, ysigma=None, normalisation=False, number_bootstraps=None, random_seed=42)[source]#
Calculates a 2D histogram of data, allowing for normalisation and bootstrap sampling
- Parameters:
xvalues (numpy.ndarray) – Data values of the first variable
yvalues (numpy.ndarray) – Data values of the second variable
xbins (numpy.ndarray) – Bin edges for the first variable
ybins (numpy.ndarray) – Bin edges for the second variable
boundaries (list) – List of (Lower, Upper) tuples corresponding to the bounds of the two data sets
xsigma (numpy.ndarray) – Error values (standard deviatons) on first variable
ysigma (numpy.ndarray) – Error values (standard deviatons) on second variable
normalisation (bool) – If True then returns the histogram as a density function
number_bootstraps (int) – Number of bootstraps
random_seed – Seed used in the random number generator
returns – 2-D histogram of data
- openquake.hmtk.seismicity.utils.decimal_time(year, month, day, hour, minute, second)[source]#
Returns the full time as a decimal value
- Parameters:
year – Year of events (integer numpy.ndarray)
month – Month of events (integer numpy.ndarray)
day – Days of event (integer numpy.ndarray)
hour – Hour of event (integer numpy.ndarray)
minute – Minute of event (integer numpy.ndarray)
second – Second of event (float numpy.ndarray)
- Returns decimal_time:
Decimal representation of the time (as numpy.ndarray)
- openquake.hmtk.seismicity.utils.decimal_year(year, month, day)[source]#
Allows to calculate the decimal year for a vector of dates (TODO this is legacy code kept to maintain comparability with previous declustering algorithms!)
- Parameters:
year (numpy.ndarray) – year column from catalogue matrix
month (numpy.ndarray) – month column from catalogue matrix
day (numpy.ndarray) – day column from catalogue matrix
- Returns:
decimal year column
- Return type:
numpy.ndarray
- openquake.hmtk.seismicity.utils.greg2julian(year, month, day, hour, minute, second)[source]#
Function to convert a date from Gregorian to Julian format
- Parameters:
year – Year of events (integer numpy.ndarray)
month – Month of events (integer numpy.ndarray)
day – Days of event (integer numpy.ndarray)
hour – Hour of event (integer numpy.ndarray)
minute – Minute of event (integer numpy.ndarray)
second – Second of event (float numpy.ndarray)
- Returns julian_time:
Julian representation of the time (as float numpy.ndarray)
- openquake.hmtk.seismicity.utils.haversine(lon1, lat1, lon2, lat2, radians=False, earth_rad=6371.227)[source]#
Allows to calculate geographical distance using the haversine formula.
- Parameters:
lon1 (numpy.ndarray) – longitude of the first set of locations
lat1 (numpy.ndarray) – latitude of the frist set of locations
lon2 (numpy.float64) – longitude of the second set of locations
lat2 (numpy.float64) – latitude of the second set of locations
radians (bool) – states if locations are given in terms of radians
earth_rad (float) – radius of the earth in km
- Returns:
geographical distance in km
- Return type:
numpy.ndarray
- openquake.hmtk.seismicity.utils.hmtk_histogram_1D(values, intervals, offset=1e-10)[source]#
So, here’s the problem. We tend to refer to certain data (like magnitudes) rounded to the nearest 0.1 (or similar, i.e. 4.1, 5.7, 8.3 etc.). We also like our tables to fall on on the same interval, i.e. 3.1, 3.2, 3.3 etc. We usually assume that the counter should correspond to the low edge, i.e. 3.1 is in the group 3.1 to 3.2 (i.e. L <= M < U). Floating point precision can be a bitch! Because when we read in magnitudes from files 3.1 might be represented as 3.0999999999 or as 3.1000000000001 and this is seemingly random. Similarly, if np.arange() is used to generate the bin intervals then we see similar floating point problems emerging. As we are frequently encountering density plots with empty rows or columns where data should be but isn’t because it has been assigned to the wrong group.
Instead of using numpy’s own historgram function we use a slower numpy version that allows us to offset the intervals by a smaller amount and ensure that 3.0999999999, 3.0, and 3.10000000001 would fall in the group 3.1 - 3.2!
- Parameters:
values (numpy.ndarray) – Values of data
intervals (numpy.ndarray) – Data bins
offset (float) – Small amount to offset the bins for floating point precision
- Returns:
Count in each bin (as float)
- openquake.hmtk.seismicity.utils.hmtk_histogram_2D(xvalues, yvalues, bins, x_offset=1e-10, y_offset=1e-10)[source]#
See the explanation for the 1D case - now applied to 2D.
- Parameters:
xvalues (numpy.ndarray) – Values of x-data
yvalues (numpy.ndarray) – Values of y-data
bins (tuple) – Tuple containing bin intervals for x-data and y-data (as numpy arrays)
x_offset (float) – Small amount to offset the x-bins for floating point precision
y_offset (float) – Small amount to offset the y-bins for floating point precision
- Returns:
Count in each bin (as float)
- openquake.hmtk.seismicity.utils.leap_check(year)[source]#
Returns logical array indicating if year is a leap year
- openquake.hmtk.seismicity.utils.lonlat_to_laea(lon, lat, lon0, lat0, f_e=0.0, f_n=0.0)[source]#
Converts vectors of longitude and latitude into Lambert Azimuthal Equal Area projection (km), with respect to an origin point
- Parameters:
lon (numpy.ndarray) – Longitudes
lat (numpy.ndarray) – Latitude
lon0 (float) – Central longitude
lat0 (float) – Central latitude
f_e (float) – False easting (km)
f_e – False northing (km)
- Returns:
easting (km)
northing (km)
- openquake.hmtk.seismicity.utils.piecewise_linear_scalar(params, xval)[source]#
Piecewise linear function for a scalar variable xval (float).
- Parameters:
params – Piecewise linear parameters (numpy.ndarray) in the following form: [slope_i,… slope_n, turning_point_i, …, turning_point_n, intercept] Length params === 2 * number_segments, e.g. [slope_1, slope_2, slope_3, turning_point1, turning_point_2, intercept]
xval – Value for evaluation of function (float)
- Returns:
Piecewise linear function evaluated at point xval (float)
- openquake.hmtk.seismicity.utils.sample_truncated_gaussian_vector(data, uncertainties, bounds=None)[source]#
Samples a Gaussian distribution subject to boundaries on the data
- Parameters:
data (numpy.ndarray) – Vector of N data values
uncertainties (numpy.ndarray) – Vector of N data uncertainties
number_bootstraps (int) – Number of bootstrap samples
bounds (tuple) – (Lower, Upper) bound of data space