3.3. quakemigrate.io
The quakemigrate.io
module handles the various input/output operations performed
by QuakeMigrate. This includes:
Reading waveform data - The submodule data.py can handle any waveform data archive with a regular directory structure. It also provides functions for checking data quality and removing/simulating instrument reponse.
Reading station files, velocity model files, instrument response inventories and QuakeMigrate lookup tables.
The
Run
class encapsulates all i/o path information and logger configuration for a given QuakeMigrate run.The
Event
class encapsulates waveforms, coalescence information, picks and location information for a given event, and provides functionality to write “.event” files.Reading and writing results, including station availablity data and continuous coalescence output from detect; triggered event files from trigger, amplitude and local magnitude measurements and cut waveforms for located events.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
3.3.1. quakemigrate.io.amplitudes
Module to handle input/output of .amps files.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- quakemigrate.io.amplitudes.write_amplitudes(run, amplitudes, event)[source]
Write amplitude results to a new .amps file. This includes amplitude measurements, and the magnitude estimates derived from them (with station correction terms applied, if provided).
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.amplitudes (pandas.DataFrame object) –
P- and S-wave amplitude measurements for each component of each station in the station file, and individual local magnitude estimates derived from them. Columns = [“epi_dist”, “z_dist”, “P_amp”, “P_freq”, “P_time”,
”P_avg_amp”, “P_filter_gain”, “S_amp”, “S_freq”, “S_time”, “S_avg_amp”, “S_filter_gain”, “Noise_amp”, “is_picked”, “ML”, “ML_Err”]
Index = Trace ID (see obspy.Trace object property ‘id’)
event (
Event
object) – Light class encapsulating waveforms, coalescence information, picks and location information for a given event.
3.3.2. quakemigrate.io.availability
Module to handle input/output of StationAvailability.csv files.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- quakemigrate.io.availability.read_availability(run, starttime, endtime)[source]
Read in station availability data to a pandas.DataFrame from csv files split by Julian day.
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.starttime (obspy.UTCDateTime object) – Timestamp from which to read the station availability.
endtime (obspy.UTCDateTime object) – Timestamp up to which to read the station availability.
- Returns
availability – Details the availability of each station for each timestep of detect.
- Return type
pandas.DataFrame object
- quakemigrate.io.availability.write_availability(run, availability)[source]
Write out csv files (split by Julian day) containing station availability data.
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.availability (pandas.DataFrame object) – Details the availability of each station for each timestep of detect.
3.3.3. quakemigrate.io.core
Module to handle input/output for QuakeMigrate.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- class quakemigrate.io.core.Run(path, name, subname='', stage=None, loglevel='info')[source]
Bases:
object
Light class to encapsulate i/o path information for a given run.
- Parameters
stage (str) – Specifies run stage of QuakeMigrate (“detect”, “trigger”, or “locate”).
path (str) – Points to the top level directory containing all input files, under which the specific run directory will be created.
name (str) – Name of the current QuakeMigrate run.
subname (str, optional) – Optional name of a sub-run - useful when testing different trigger parameters, for example.
- path
Points to the top level directory containing all input files, under which the specific run directory will be created.
- Type
pathlib.Path object
- name
Name of the current QuakeMigrate run.
- Type
str
- run_path
Points to the run directory into which files will be written.
- Type
pathlib.Path object
- subname
Optional name of a sub-run - useful when testing different trigger parameters, for example.
- Type
str
- stage
Track which stage of QuakeMigrate is being run.
- Type
{“detect”, “trigger”, “locate”}, optional
- loglevel
Set the logging level. (Default “info”)
- Type
{“info”, “debug”}, optional
- logger(log)[source]
Configures the logging feature.
- Parameters
log (bool) – Toggle for logging. If True, will output to stdout and generate a log file.
- property name
Get the run name as a formatted string.
- quakemigrate.io.core.read_lut(lut_file)[source]
Read the contents of a pickle file and restore state of the lookup table object.
- Parameters
lut_file (str) – Path to pickle file to load.
- Returns
lut – Lookup table populated with grid specification and traveltimes.
- Return type
LUT
object
- quakemigrate.io.core.read_response_inv(response_file, sac_pz_format=False)[source]
Reads response information from file, returning it as a obspy.Inventory object.
- Parameters
response_file (str) – Path to response file. Please see the obspy.read_inventory() documentation for a full list of supported file formats. This includes a dataless.seed volume, a concatenated series of RESP files or a stationXML file.
sac_pz_format (bool, optional) – Toggle to indicate that response information is being provided in SAC Pole-Zero files. NOTE: not yet supported.
- Returns
response_inv – ObsPy response inventory.
- Return type
obspy.Inventory object
- Raises
NotImplementedError – If the user selects sac_pz_format=True.
TypeError – If the user provides a response file that is not readable by ObsPy.
- quakemigrate.io.core.read_stations(station_file, **kwargs)[source]
Reads station information from file.
- Parameters
station_file (str) –
Path to station file. File format (header line is REQUIRED, case sensitive, any order):
Latitude, Longitude, Elevation (units matching LUT grid projection; either metres or kilometres; positive upwards), Name
kwargs (dict) – Passthrough for pandas.read_csv kwargs.
- Returns
stn_data – Columns: “Latitude”, “Longitude”, “Elevation”, “Name”
- Return type
pandas.DataFrame object
- Raises
StationFileHeaderException – Raised if the input file is missing required entries in the header.
- quakemigrate.io.core.read_vmodel(vmodel_file, **kwargs)[source]
Reads velocity model information from file.
- Parameters
vmodel_file (str) –
Path to velocity model file. File format: (header line is REQUIRED, case sensitive, any order):
”Depth” of each layer in the model (units matching the LUT grid projection; positive-down) “V<phase>” velocity for each layer in the model, for each phase the user wishes to calculate traveltimes for (units matching the LUT grid projection). There are no required phases, and no maximum number of separate phases. E.g. “Vp”, “Vs”, “Vsh”.
kwargs (dict) – Passthrough for pandas.read_csv kwargs.
- Returns
vmodel_data –
- Columns:
”Depth” of each layer in model (positive down) “V<phase>” velocity for each layer in model (e.g. “Vp”)
- Return type
pandas.DataFrame object
- Raises
VelocityModelFileHeaderException – Raised if the input file is missing required entries in the header.
3.3.4. quakemigrate.io.cut_waveforms
Module to handle input/output of cut waveforms.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- quakemigrate.io.cut_waveforms.get_waveforms(st, event, waveform_type, units)[source]
Get real or simulated waveforms for a Stream.
- Parameters
st (obspy.Stream object) – Stream for which to get real or simulated waveforms.
event (
Event
object) – Light class encapsulating waveforms, coalescence information, picks and location information for a given event.waveform_type ({"real", "wa"}) – Whether to get real or Wood-Anderson simulated waveforms.
units ({"displacement", "velocity"}) – Units to return waveforms in.
- Returns
st_out – Stream of real or Wood-Anderson simulated waveforms in the requested units.
- Return type
obspy.Stream object
- quakemigrate.io.cut_waveforms.write_cut_waveforms(run, event, file_format, pre_cut=0.0, post_cut=0.0, waveform_type='raw', units='displacement')[source]
Output cut waveform data as a waveform file – defaults to miniSEED format.
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.event (
Event
object) – Light class encapsulating waveforms, coalescence information, picks and location information for a given event.file_format (str, optional) – File format to write waveform data to. Options are all file formats supported by ObsPy, including: “MSEED” (default), “SAC”, “SEGY”, “GSE2”
pre_cut (float or None, optional) – Specify how long before the event origin time to cut the waveform data from.
post_cut (float or None, optional) – Specify how long after the event origin time to cut the waveform data to.
waveform_type ({"raw", "real", "wa"}, optional) – Whether to output raw, real or Wood-Anderson simulated waveforms. Default: “raw”
units ({"displacement", "velocity"}, optional) – Whether to output displacement waveforms or velocity waveforms for real/ Wood-Anderson corrected traces. Default: displacement
- Raises
AttributeError – If real or wa waveforms are requested and no response inventory has been provided.
- quakemigrate.io.cut_waveforms.write_waveforms(st, fpath, fstem, file_format)[source]
Output waveform data as a waveform file – defaults to miniSEED format.
- Parameters
st (obspy.Stream object) – Waveforms to be written to file.
fpath (pathlib.Path object) – Path to output directory.
fstem (str) – File name (without suffix).
file_format (str) – File format to write waveform data to. Options are all file formats supported by ObsPy, including: “MSEED” (default), “SAC”, “SEGY”, “GSE2”
3.3.5. quakemigrate.io.data
Module for processing waveform files stored in a data archive.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- class quakemigrate.io.data.Archive(archive_path, stations, archive_format=None, **kwargs)[source]
Bases:
object
The Archive class handles the reading of archived waveform data.
It is capable of handling any regular archive structure. Requests to read waveform data are served up as a
WaveformData
object.If provided, a response inventory for the archive will be stored with the waveform data for response removal, if needed (e.g. for local magnitude calculation, or to output real cut waveforms).
By default, data with mismatched sampling rates will only be decimated. If necessary, and if the user specifies resample = True and an upfactor to upsample by upfactor = int for the waveform archive, data can also be upsampled and then, if necessary, subsequently decimated to achieve the desired sampling rate.
For example, for raw input data sampled at a mix of 40, 50 and 100 Hz, to achieve a unified sampling rate of 50 Hz, the user would have to specify an upfactor of 5; 40 Hz x 5 = 200 Hz, which can then be decimated to 50 Hz - see
resample()
.- Parameters
archive_path (str) – Location of seismic data archive: e.g.: “./DATA_ARCHIVE”.
stations (pandas.DataFrame object) – Station information. Columns [“Latitude”, “Longitude”, “Elevation”, “Name”]. See
read_stations()
archive_format (str, optional) – Sets directory structure and file naming format for different archive formats. See
path_structure()
kwargs (**dict) – See Archive Attributes for details.
- archive_path
Location of seismic data archive: e.g.: ./DATA_ARCHIVE.
- Type
pathlib.Path object
- stations
Series object containing station names.
- Type
pandas.Series object
- format
Directory structure and file naming format of data archive.
- Type
str
- read_all_stations
If True, read all stations in archive for that time period. Else, only read specified stations.
- Type
bool, optional
- resample
If true, perform resampling of data which cannot be decimated directly to the desired sampling rate. See
resample()
- Type
bool, optional
- response_inv
ObsPy response inventory for this waveform archive, containing response information for each channel of each station of each network.
- Type
obspy.Inventory object, optional
- pre_filt
Pre-filter to apply during the instrument response removal. E.g. (0.03, 0.05, 30., 35.) - all in Hz. (Default None)
- Type
tuple of floats
- water_level
Water level to use in instrument response removal. (Default 60.)
- Type
float
- remove_full_response
Whether to remove the full response (including the effect of digital FIR filters) or just the instrument transform function (as defined by the PolesZeros Response Stage). Significantly slower. (Default False)
- Type
bool
- upfactor
Factor by which to upsample the data to enable it to be decimated to the desired sampling rate, e.g. 40Hz -> 50Hz requires upfactor = 5. See
resample()
- Type
int, optional
- interpolate
If data is timestamped “off-sample” (i.e. a non-integer number of samples after midnight), whether to interpolate the data to apply the necessary correction. Default behaviour is to just alter the metadata, resulting in a sub-sample timing offset. See
shift_to_sample()
.- Type
bool, optional
- path_structure(archive_format, channels='*')[source]
Set the directory structure and file naming format of the data archive.
- path_structure(archive_format='YEAR/JD/STATION', channels='*')[source]
Define the directory structure and file naming format of the data archive.
- Parameters
archive_format (str, optional) – Directory structure and file naming format of the data archive. This may be the name of a generic archive format (e.g. SeisComp3), or one of a selection of additional formats built into QuakeMigrate.
channels (str, optional) – Channel codes to include. E.g. channels=”[B,H]H*”. (Default “*”)
- Raises
ArchivePathStructureError – If the archive_format specified by the user is not a valid option.
- read_waveform_data(starttime, endtime, pre_pad=0.0, post_pad=0.0)[source]
Read in waveform data from the archive between two times.
Supports all formats currently supported by ObsPy, including: “MSEED”, “SAC”, “SEGY”, “GSE2”.
Optionally, read data with some pre- and post-pad, and for all stations in the archive - this will be stored in data.raw_waveforms, while data.waveforms will contain only data for selected stations between starttime and endtime.
- Parameters
starttime (obspy.UTCDateTime object) – Timestamp from which to read waveform data.
endtime (obspy.UTCDateTime object) – Timestamp up to which to read waveform data.
pre_pad (float, optional) – Additional pre pad of data to read. Defaults to 0.
post_pad (float, optional) – Additional post pad of data to read. Defaults to 0.
- Returns
data – Object containing the waveform data read from the archive that satisfies the query.
- Return type
WaveformData
object- Raises
ArchiveEmptyException – If no data files are found in the archive for this day(s).
DataAvailabilityException – If no data is found in the archive for the specified stations within the specified time window.
- class quakemigrate.io.data.WaveformData(starttime, endtime, stations=None, response_inv=None, water_level=60.0, pre_filt=None, remove_full_response=False, read_all_stations=False, resample=False, upfactor=None, pre_pad=0.0, post_pad=0.0)[source]
Bases:
object
The WaveformData class encapsulates the waveform data returned by an Archive query.
It also provides a number of utility functions. These include removing instrument response and checking data availability against a flexible set of data quality criteria.
- Parameters
starttime (obspy.UTCDateTime object) – Timestamp of first sample of waveform data requested from the archive.
endtime (obspy.UTCDateTime object) – Timestamp of last sample of waveform data requested from the archive.
stations (pandas.Series object, optional) – Series object containing station names.
read_all_stations (bool, optional) – If True, raw_waveforms contain all stations in archive for that time period. Else, only selected stations will be included.
resample (bool, optional) – If true, allow resampling of data which cannot be decimated directly to the desired sampling rate. See
resample()
Default: Falseupfactor (int, optional) – Factor by which to upsample the data to enable it to be decimated to the desired sampling rate, e.g. 40Hz -> 50Hz requires upfactor = 5. See
resample()
response_inv (obspy.Inventory object, optional) – ObsPy response inventory for this waveform data, containing response information for each channel of each station of each network.
pre_filt (tuple of floats) – Pre-filter to apply during the instrument response removal. E.g. (0.03, 0.05, 30., 35.) - all in Hz. (Default None)
water_level (float) – Water level to use in instrument response removal. (Default 60.)
remove_full_response (bool) – Whether to remove the full response (including the effect of digital FIR filters) or just the instrument transform function (as defined by the PolesZeros Response Stage). Significantly slower. (Default False)
pre_pad (float, optional) – Additional pre pad of data included in raw_waveforms.
post_pad (float, optional) – Additional post pad of data included in raw_waveforms.
- starttime
Timestamp of first sample of waveform data requested from the archive.
- Type
obspy.UTCDateTime object
- endtime
Timestamp of last sample of waveform data requested from the archive.
- Type
obspy.UTCDateTime object
- stations
Series object containing station names.
- Type
pandas.Series object
- read_all_stations
If True, raw_waveforms contain all stations in archive for that time period. Else, only selected stations will be included.
- Type
bool
- raw_waveforms
Raw seismic data read in from the archive. This may be for all stations in the archive, or only those specified by the user. See read_all_stations. It may also cover the time period between starttime and endtime, or feature an additional pre- and post-pad. See pre_pad and post_pad.
- Type
obspy.Stream object
- waveforms
Seismic data read in from the archive for the specified list of stations, between starttime and endtime.
- Type
obspy.Stream object
- pre_pad
Additional pre pad of data included in raw_waveforms.
- Type
float
- post_pad
Additional post pad of data included in raw_waveforms.
- Type
float
- check_availability(stream, \*\*data_quality_params)[source]
Check data availability against a set of data quality criteria.
- get_wa_waveform(trace, \*\*response_removal_params)[source]
Calculate the Wood-Anderson corrected waveform for a obspy.Trace object.
- Raises
NotImplementedError – If the user attempts to use the get_real_waveform() method.
- check_availability(st, all_channels=False, n_channels=None, allow_gaps=False, full_timespan=True, check_sampling_rate=False, sampling_rate=None, check_start_end_times=False)[source]
Check waveform availability against data quality criteria.
There are a number of hard-coded checks: for whether any data is present; for whether the data is a flatline (all samples have the same value); and for whether the data contains overlaps. There are a selection of additional optional checks which can be specified according to the onset function / user preference.
- Parameters
st (obspy.Stream object) – Stream containing the waveform data to check against the availability criteria.
all_channels (bool, optional) – Whether all supplied channels (distinguished by SEED id) need to meet the availability criteria to mark the data as ‘available’.
n_channels (int, optional) – If all_channels=True, this argument is required (in order to specify the number of channels expected to be present).
allow_gaps (bool, optional) – Whether to allow gaps.
full_timespan (bool, optional) – Whether to ensure the data covers the entire timespan requested; note that this implicitly requires that there be no gaps. Checks the number of samples in the trace, not the start and end times; for that see check_start_end_times.
check_sampling_rate (bool, optional) – Check that all channels are at the desired sampling rate.
sampling_rate (float, optional) – If check_sampling_rate=True, this argument is required to specify the sampling rate that the data should be at.
check_start_end_times (bool, optional) – A stricter alternative to full_timespan; checks that the first and last sample of the trace have exactly the requested timestamps.
- Returns
available (int) – 0 if data doesn’t meet the availability requirements; 1 if it does.
availability (dict) – Dict of {tr_id : available} for each unique SEED ID in the input stream (available is again 0 or 1).
- Raises
TypeError – If the user specifies all_channels=True but does not specify n_channels.
- get_real_waveform(tr, velocity=True)[source]
Calculate the real waveform for a Trace by removing the instrument response.
- Parameters
tr (obspy.Trace object) – Trace containing the waveform for which to remove the instrument response.
velocity (bool, optional) – Output velocity waveform (as opposed to displacement). Default: True.
- Returns
tr – Trace with instrument response removed.
- Return type
obspy.Trace object
- Raises
AttributeError – If no response inventory has been supplied.
ResponseNotFoundError – If the response information for a trace can’t be found in the supplied response inventory.
ResponseRemovalError – If the deconvolution of the instrument response is unsuccessful.
- get_wa_waveform(tr, velocity=False)[source]
Calculate simulated Wood Anderson displacement waveform for a Trace.
- Parameters
tr (obspy.Trace object) – Trace containing the waveform to be corrected to a Wood-Anderson response.
velocity (bool, optional) – Output velocity waveform, instead of displacement. Default: False. NOTE: all attenuation functions provided within the QM local_mags module are calculated for displacement seismograms.
- Returns
tr – Trace corrected to Wood-Anderson response.
- Return type
obspy.Trace object
3.3.6. quakemigrate.io.event
Module containing the Event class, which stores information related to an individual event.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- class quakemigrate.io.event.Event(marginal_window, triggered_event=None)[source]
Bases:
object
Light class to encapsulate information about an event, including waveform data, coalescence information, origin time, locations, picks, magnitudes.
- Parameters
marginal_window (float) – Estimate of the uncertainty in the event origin time; time window over which the 4-D coalescence image is marginalised around the peak coalescence time (event origin time) to produce the 3-D coalescence map.
triggered_event (pandas.Series object, optional) – Contains information on the candidate event identified by
trigger()
- coa_data
Event coalescence data computed during locate.
- DTnumpy.ndarray of obspy.UTCDateTime objects, shape(nsamples)
Timestamps for the coalescence data.
- COAnumpy.ndarray of floats, shape(nsamples)
Max coalescence value in the grid at each timestep.
- COA_NORMnumpy.ndarray of floats, shape(nsamples)
Normalised max coalescence value in the grid at each timestep.
- Xnumpy.ndarray of floats, shape(nsamples)
X coordinate of maximum coalescence value in the grid at each timestep, in input (geographic) projection coordinates.
- Ynumpy.ndarray of floats, shape(nsamples)
Y coordinate of maximum coalescence value in the grid at each timestep, in input (geographic) projection coordinates.
- Znumpy.ndarray of floats, shape(nsamples)
Z coordinate of maximum coalescence value in the grid at each timestep, in input (geographic) projection coordinates.
- Type
pandas.DataFrame object
- data
Light class encapsulating waveform data returned from an archive query.
- Type
WaveformData
object
- hypocentre
[X, Y, Z]; Geographical coordinates of the event hypocentre (default is interpolated peak of a spline function fitted to the marginalised 3-D coalescence map).
- Type
numpy.ndarray of floats
- locations
Information on the various locations and reported uncertainties.
- splinedict
The location of the peak coalescence value in the marginalised 3-D coalescence map, interpolated using a 3-D spline. If no spline fit was able to be made, it is just the gridded peak location.
- gaussiandict
The location and uncertainty as determined by fitting a 3-D Gaussian to the marginalised 3-D coalescence map in a small region around the (gridded) peak coalescence location.
- covariancedict
The location and uncertainty as determined by calculating the covariance of the coalescence values in X, Y, and Z above some percentile of the max coalescence value in the marginalised 3-D coalescence map.
- Type
dict
- map4d
4-D coalescence map generated in
locate()
.- Type
numpy.ndarray, shape(nx, ny, nz, nsamp), optional
- max_coalescence
Dictionary containing the raw and normalised maximum coalescence values in the 3-D grid at the timestamp corresponding to the instantaneous (non-marginalised) maximum coalescence value in the 4-D grid (i.e. the event origin time).
- Type
dict
- onset_data
Light class encapsulating data generated during onset calculation.
- Type
OnsetData
object
- otime
Timestamp of the instantaneous peak in the 4-D coalescence function generated in
locate()
- best estimate of the event origin time.- Type
obspy.UTCDateTime object
- trigger_info
Useful information about the triggered event to be fed forward.
- TRIG_COAfloat
The peak value of the coalescence stream used to trigger the event.
- DEC_COAfloat
The coalescence value of the “raw” maximum coalsecence stream at the trigger_time.
- DEC_COA_NORMfloat
The coalescence value of the normalised maximum coalsecence stream at the trigger_time.
- Type
dict
- trigger_time
The time of the peak in the continuous coalescence stream (output by detect) corresponding to the triggered event.
- Type
obspy.UTCDateTime object
- uid
A unique identifier for the event based on the event trigger time.
- Type
str
- add_compute_output(times, max_coa, max_coa_n, coord, map4d, onset_data)[source]
Add values returned by
_compute()
to the event.
- add_covariance_location(xyz, xyz_unc)[source]
Add the covariance location and uncertainty to the event.
- add_gaussian_location(xyz, xyz_unc)[source]
Add the gaussian location and uncertainty to the event.
- add_waveform_data(data)[source]
Add waveform data read from the archive to the event (as a
WaveformData
object).
- in_marginal_window(marginal_window)[source]
Simple test to see if event is within the marginal window around the event origin time (time of max instantaneous coalescence value).
- mw_times(marginal_window, sampling_rate)[source]
Generates timestamps for data in the window around the event trigger scanned by
_compute()
; trigger_time +/- 2*`marginal_window`.
- trim2window(marginal_window)[source]
Trim the coalescence data and map4d to the marginal window about the event origin time.
- get_hypocentre(method)[source]
Get the event hypocentre estimate calculated by a specific method; {“gaussian”, “covariance”, “spline”}.
- add_compute_output(times, max_coa, max_coa_n, coord, map4d, onset_data)[source]
Append outputs of compute to the Event object. This includes time series of the maximum coalescence values in the 3-D grid at each timestep, and their locations, the full 4-D coalescence map, and the onset data generated for migration.
- Parameters
times (numpy.ndarray of obspy.UTCDateTime objects, shape(nsamples)) – Timestamps for the coalescence data.
max_coa (numpy.ndarray of floats, shape(nsamples)) – Max coalescence value in the grid at each timestep.
max_coa_n (numpy.ndarray of floats, shape(nsamples)) – Normalised max coalescence value in the grid at each timestep.
coord (numpy.ndarray of floats, shape(nsamples, 3)) – [x, y, z] Location of maximum coalescence in the grid at each timestep, in input (geographic) projection coordinates
map4d (numpy.ndarry, shape(nx, ny, nz, nsamp)) – 4-D coalescence map.
onset_data (
OnsetData
object) – Light class encapsulating data generated during onset calculation.
- add_covariance_location(xyz, xyz_unc)[source]
Add the location determined by calculating the 3-D covariance of the marginalised coalescence map filtered above a percentile threshold.
- Parameters
xyz (numpy.ndarray of floats, shape(3)) – Geographical coordinates (lon/lat/depth) of covariance location.
xyz_unc (`numpy.ndarray’ of floats, shape(3)) – One sigma uncertainties on the covariance location (units determined by the LUT projection units).
- add_gaussian_location(xyz, xyz_unc)[source]
Add the location determined by fitting a 3-D Gaussian to a small window around the Gaussian smoothed maximum coalescence location.
- Parameters
xyz (numpy.ndarray of floats, shape(3)) – Geographical coordinates (lon/lat/depth) of Gaussian location.
xyz_unc (`numpy.ndarray’ of floats, shape(3)) – One sigma uncertainties on the Gaussian location (units determined by the LUT projection units).
- add_local_magnitude(mag, mag_err, mag_r2)[source]
Add outputs from local magnitude calculation to the Event object.
- Parameters
mag (float) – Network-averaged local magnitude estimate for the event.
mag_err (float) – (Weighted) standard deviation of the magnitude estimates from amplitude measurements on individual stations/channels.
mag_r2 (float) – r-squared statistic describing the fit of the amplitude vs. distance curve predicted by the calculated mean_mag and chosen attenuation model to the measured amplitude observations. This is intended to be used to help discriminate between ‘real’ events, for which the predicted amplitude vs. distance curve should provide a good fit to the observations, from artefacts, which in general will not.
- add_picks(pick_df, **kwargs)[source]
Add phase picks, and a selection of picker outputs and parameters.
- Parameters
pick_df (pandas.DataFrame object) – DataFrame that contains the measured picks with columns: [“Name”, “Phase”, “ModelledTime”, “PickTime”, “PickError”, “SNR”] Each row contains the phase pick from one station/phase.
**kwargs –
For
GaussianPicker
:- gaussfitsdict of dicts
Keys “station”[“phase”], each containing:
”popt” : popt “xdata” : x_data “xdata_dt” : x_data_dt “PickValue” : max_onset “PickThreshold” : threshold
- pick_windowsdict
{station : phase{window}}
window: [min_time, modelled_arrival, max_time] - all ints, referring to indices of the onset function.
- add_spline_location(xyz)[source]
Add the location determined by fitting a 3-D spline to a small window around the maximum coalescence location and interpolating.
- Parameters
xyz (numpy.ndarray of floats, shape(3)) – Geographical coordinates (lon/lat/depth) of best-fitting location.
- add_waveform_data(data)[source]
Add waveform data in the form of a
WaveformData
object.- Parameters
data (
WaveformData
object) – Contains cut waveforms - raw_waveforms may be for all stations in the archive, and include an additional pre- and post-pad; waveforms contains data only for the stations and time period required for migration.
- get_hypocentre(method='spline')[source]
Get an estimate of the event hypocentre location.
- Parameters
method ({"spline", "gaussian", "covariance"}, optional) – Which location result to return. (Default “spline”).
- Returns
ev_loc – [x_coordinate, y_coordinate, z_coordinate] of event hypocentre, in the global (geographic) coordinate system.
- Return type
numpy.ndarray of floats
- get_loc_uncertainty(method='gaussian')[source]
Get an estimate of the hypocentre location uncertainty.
- Parameters
method ({"gaussian", "covariance"}, optional) – Which location result to return. (Default “gaussian”).
- Returns
ev_loc_unc – [x_uncertainty, y_uncertainty, z_uncertainty] of event hypocentre; units are determined by the LUT projection units.
- Return type
numpy.ndarray of floats
- property hypocentre
Get an estimate of the event hypocentre location.
- Parameters
method ({"spline", "gaussian", "covariance"}, optional) – Which location result to return. (Default “spline”).
- Returns
ev_loc – [x_coordinate, y_coordinate, z_coordinate] of event hypocentre, in the global (geographic) coordinate system.
- Return type
numpy.ndarray of floats
- in_marginal_window()[source]
Test if triggered event time is within marginal window around the maximum coalescence time (origin time).
- Returns
cond – Result of test.
- Return type
bool
- property loc_uncertainty
Get an estimate of the hypocentre location uncertainty.
- Parameters
method ({"gaussian", "covariance"}, optional) – Which location result to return. (Default “gaussian”).
- Returns
ev_loc_unc – [x_uncertainty, y_uncertainty, z_uncertainty] of event hypocentre; units are determined by the LUT projection units.
- Return type
numpy.ndarray of floats
- property local_magnitude
Get the local magnitude, if it exists.
- property max_coalescence
Get information related to the maximum coalescence.
- mw_times(sampling_rate)[source]
Utility function to generate timestamps for the time period around the trigger time for which the 4-D coalescence function is calculated in
_compute()
.- Returns
times – Timestamps for time range trigger_time +/- 2 * marginal_window.
- Return type
numpy.ndarray of obspy.UTCDateTime, shape(nsamples)
3.3.7. quakemigrate.io.scanmseed
Module to handle input/output of .scanmseed files.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- class quakemigrate.io.scanmseed.ScanmSEED(run, continuous_write, sampling_rate)[source]
Bases:
object
Light class to encapsulate the data output by the detect stage of QuakeMigrate. This data is stored in an obspy.Stream object with the channels: [“COA”, “COA_N”, “X”, “Y”, “Z”].
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.continuous_write (bool) – Option to continuously write the .scanmseed file output by
detect()
at the end of every time step. Default behaviour is to write in day chunks where possible.sampling_rate (int) – Desired sampling rate of input data; sampling rate at which to compute the coalescence function. Default: 50 Hz.
- stream
Output of
detect()
stored in obspy.Stream object. The values have been multiplied by a factor to make use of more efficient compression. Channels: [“COA”, “COA_N”, “X”, “Y”, “Z”]- Type
obspy.Stream object
- written
Tracker for whether the data appended has been written recently.
- Type
bool
- append(times, max_coa, max_coa_n, coord, map4d=None)[source]
Append the output of
_compute()
to the coalescence stream.
- empty(starttime, timestep, i, msg)[source]
Create an set of empty arrays for a given timestep and append to the coalescence stream.
- append(starttime, max_coa, max_coa_n, coord, ucf)[source]
Append latest timestep of
detect()
output to obspy.Stream object.Multiply channels [“COA”, “COA_N”, “X”, “Y”, “Z”] by factors of [“1e5”, “1e5”, “1e6”, “1e6”, “1e3”] respectively, round and convert to int32 as this dramatically reduces memory usage, and allows the coastream data to be saved in mSEED format with STEIM2 compression. The multiplication factor is removed when the data is read back in.
- Parameters
starttime (obspy.UTCDateTime object) – Timestamp of first sample of coalescence data.
max_coa (numpy.ndarray of floats, shape(nsamples)) – Coalescence value through time.
max_coa_n (numpy.ndarray of floats, shape(nsamples)) – Normalised coalescence value through time.
coord (numpy.ndarray of floats, shape(nsamples)) – Location of maximum coalescence through time in input projection space.
ucf (float) – A conversion factor based on the lookup table grid projection. Used to ensure the same level of precision (millimetre) is retained during compression, irrespective of the units of the grid projection.
- empty(starttime, timestep, i, msg, ucf)[source]
Create an empty set of arrays to write to .scanmseed; used where there is no data available to run
_compute()
.- Parameters
starttime (obspy.UTCDateTime object) – Timestamp of first sample in the given timestep.
timestep (float) – Length (in seconds) of timestep used in detect().
i (int) – The ith timestep of the continuous compute.
msg (str) – Message to output to log giving details as to why this timestep is empty.
ucf (float) – A conversion factor based on the lookup table grid projection. Used to ensure the same level of precision (millimetre) is retained during compression, irrespective of the units of the grid projection.
- write(write_start=None, write_end=None)[source]
Write a new .scanmseed file from an obspy.Stream object containing the data output from detect(). Note: values have been multiplied by a power of ten, rounded and converted to an int32 array so the data can be saved as mSEED with STEIM2 compression. This multiplication factor is removed when the data is read back in with read_scanmseed().
- Parameters
write_start (obspy.UTCDateTime object, optional) – Timestamp from which to write the coalescence stream to file.
write_end (obspy.UTCDateTime object, optional) – Timestamp up to which to write the coalescence stream to file.
- quakemigrate.io.scanmseed.read_scanmseed(run, starttime, endtime, pad, ucf)[source]
Read .scanmseed files between two time stamps. Files are labelled by year and Julian day.
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.starttime (obspy.UTCDateTime object) – Timestamp from which to read the coalescence stream.
endtime (obspy.UTCDateTime object) – Timestamp up to which to read the coalescence stream.
pad (float) – Read in “pad” seconds of additional data on either end.
ucf (float) – A conversion factor based on the lookup table grid projection. Used to ensure the same level of precision (millimetre) is retained during compression, irrespective of the units of the grid projection.
- Returns
data (pandas.DataFrame object) – Data output by detect() – decimated scan. Columns: [“DT”, “COA”, “COA_N”, “X”, “Y”, “Z”] - X/Y/Z as lon/lat/units where units is the user-selected units of the lookup table grid projection (either metres or kilometres).
stats (obspy.trace.Stats object) – Container for additional header information for coalescence trace. Contains keys: network, station, channel, starttime, endtime,
sampling_rate, delta, npts, calib, _format, mseed
3.3.8. quakemigrate.io.triggered_events
Module to handle input/output of TriggeredEvents.csv files.
- copyright
2020–2023, QuakeMigrate developers.
- license
GNU General Public License, Version 3 (https://www.gnu.org/licenses/gpl-3.0.html)
- quakemigrate.io.triggered_events.read_triggered_events(run, **kwargs)[source]
Read triggered events from .csv file.
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.starttime (obspy.UTCDateTime object, optional) – Timestamp from which to include events in the locate scan.
endtime (obspy.UTCDateTime object, optional) – Timestamp up to which to include events in the locate scan.
trigger_file (str, optional) – File containing triggered events to be located.
- Returns
events – Triggered events information. Columns: [“EventID”, “CoaTime”, “TRIG_COA”, “COA_X”, “COA_Y”, “COA_Z”, “COA”, “COA_NORM”].
- Return type
pandas.DataFrame object
- quakemigrate.io.triggered_events.write_triggered_events(run, events, starttime)[source]
Write triggered events to a .csv file.
- Parameters
run (
Run
object) – Light class encapsulating i/o path information for a given run.events (pandas.DataFrame object) – Triggered events information. Columns: [“EventID”, “CoaTime”, “TRIG_COA”, “COA_X”, “COA_Y”, “COA_Z”, “COA”, “COA_NORM”].
starttime (obspy.UTCDateTime object) – Timestamp from which events have been triggered.