4.7.1.7. eqcorrscan.utils.mag_calc.SVD_moments

eqcorrscan.utils.mag_calc.SVD_moments(U, s, V, stachans, event_list, n_SVs=4)[source]

Calculate relative moments/amplitudes using singular-value decomposition.

Convert basis vectors calculated by singular value decomposition (see the SVD functions in clustering) into relative moments.

For more information see the paper by Rubinstein & Ellsworth (2010).

Parameters:
  • U (list) List of the numpy.ndarray input basis vectors from the SVD, one array for each channel used.
  • s (list) List of the numpy.ndarray of singular values, one array for each channel.
  • V (list) List of numpy.ndarray of output basis vectors from SVD, one array per channel.
  • stachans (list) List of station.channel input
  • event_list (list) List of events for which you have data, such that event_list[i] corresponds to stachans[i], U[i] etc. and event_list[i][j] corresponds to event j in U[i]. These are a series of indexes that map the basis vectors to their relative events and channels - if you have every channel for every event generating these is trivial (see example).
  • n_SVs (int) Number of singular values to use, defaults to 4.
Returns:

M, array of relative moments

Return type:

numpy.ndarray

Returns:

events_out, list of events that relate to M (in order), does not include the magnitude information in the events, see note.

Return type:

obspy.core.event.event.Event

Note

M is an array of relative moments (or amplitudes), these cannot be directly compared to true moments without calibration.

Example

>>> from eqcorrscan.utils.mag_calc import SVD_moments
>>> from obspy import read
>>> import glob
>>> import os
>>> from eqcorrscan.utils.clustering import SVD
>>> import numpy as np
>>> # Do the set-up
>>> testing_path = 'eqcorrscan/tests/test_data/similar_events'
>>> stream_files = glob.glob(os.path.join(testing_path, '*'))
>>> stream_list = [read(stream_file) for stream_file in stream_files]
>>> event_list = []
>>> for i, stream in enumerate(stream_list):
...     st_list = []
...     for tr in stream:
...         if (tr.stats.station, tr.stats.channel) not in
...           [('WHAT2', 'SH1'), ('WV04', 'SHZ'), ('GCSZ', 'EHZ')]:
...             stream.remove(tr)
...             continue
...         tr.detrend('simple')
...         tr.filter('bandpass', freqmin=5.0, freqmax=15.0)
...         tr.trim(tr.stats.starttime + 40, tr.stats.endtime - 45)
...         st_list.append(i)
...     event_list.append(st_list) 
>>> event_list = np.asarray(event_list).T.tolist()
>>> SVec, SVal, U, stachans = SVD(stream_list=stream_list) 
['GCSZ.EHZ', 'WV04.SHZ', 'WHAT2.SH1']
>>> M, events_out = SVD_moments(U=U, s=SVal, V=SVec, stachans=stachans,
...                             event_list=event_list)