Utility Functions (seislib.utils)

This module provides support for several useful functions for processing and analysing seismic data. These mainly build on ObsPy, and include:

  • Vectorized calculation of great-circle distances, azimuths, and back azimuths

  • Slicing / resampling of pairs of seismograms so as to adapt them to a common time window / sampling rate

  • Pre-processing operations such as zero-padding and bandpass filtering

  • Rotation of horizontal-component recordings to an arbitrary amount of degrees. SeisLib’s implementation also takes into account not uncommon ValueErrors raised by ObsPy due to the presence of slightly (sub-sample) differences in the time spanned by the two traces

  • Interpolation of scattered data on equal-area grids

  • Calculation of the Pearson correlation coefficient between physical parameters defined on equal-area grids

seislib.utils._utils.adapt_sampling_rate(st1, st2)

If the input streams (or traces) have different sampling rates, the one characterized by the largest sampling rate is downsampled to the sampling rate of the other stream (or trace).

The downsampling is carried out via the resample(), which modifies the input streams in place.

Parameters
st1, st2obspy.Stream, obspy.Trace
Returns
st1, st2obspy.Stream, obspy.Trace

Obspy stream or trace depending on the input. The input is permanently modified

seislib.utils._utils.adapt_timespan(st1, st2)

Slices all traces from the two input streams to the overlapping timerange. Then returns a copy of the sliced streams.

Parameters
st1, st2obspy.Stream, obspy.Trace
Returns
st1, st2obspy.Stream, obspy.Trace

Obspy stream or trace depending on the input. The original input is not permanently modified (a copy is returned)

Raises
TimeSpanException

If no overlap is available

Notes

The maximum precision achieved by this function is governed by the samling rate. If sub-sample precision is required, consider using adapt_timespan_interpolate()

seislib.utils._utils.adapt_timespan_interpolate(st1, st2, min_overlap=0)

Slices all traces from the two input streams to the overlapping timerange. Then returns a copy of the sliced streams. If the starttime of the sliced traces do not fit exactly (because of sub-sample time shifts), the traces are interpolated to remove the time shift.

Parameters
st1, st2obspy.Stream, obspy.Trace
Returns
st1_out, st2_outobspy.Stream, obspy.Trace

Obspy stream or trace depending on the input. The original input is not permanently modified (a copy is returned)

Raises
TimeSpanException

If no overlap is available

Notes

Interpolation can require a relatively long time depending on the size of the sliced stream. If speed is preferred to (sub-sample) precision, consider using adapt_timespan()

seislib.utils._utils.azimuth_backazimuth(lat1, lon1, lat2, lon2)

Calculates the azimuth and backazimuth (in degrees) between coordinate points (in degrees). This function calls directly the obspy gps2dist_azimuth, it only extends its functionality through the numpy.vectorize decorator

Parameters
lat1, lon1, lat2, lon2float or array-like of shape (n,)

Coordinates of the points on the Earth’s surface, in degrees.

Returns
tuple of shape (2,) or ndarray of shape (n, 2)

The columns consist of azimuth azimuth and backazimuth

seislib.utils._utils.bandpass_gaussian(data, dt, period, alpha)

Gaussian filter of real-valued data carried out in the frequency domain

The bandpass filter is carried out with a Gaussian filter centered at period, whose width is controlled by alpha:

exp(-alpha * ((f-f0)/f0)**2)

where f is frequency and f0 = 1 / period.

Parameters
datandarray of shape (n,)

Real-valued data to be filtered

dtfloat

Time sampling interval of the data

periodfloat

Central period, around which the (tight) bandapass filter is carried out

alphafloat

Parameter that controls the width of the Gaussian filter

Returns
numpy.ndarray of shape (n,)

Filtered data

seislib.utils._utils.gc_distance(lat1, lon1, lat2, lon2)

Calculates the great circle distance (in m) between coordinate points (in degrees). This function calls directly the obspy gps2dist_azimuth, it only extends its functionality through the numpy.vectorize decorator.

Parameters
lat1, lon1, lat2, lon2float or array-like of shape (n,)

Coordinates of the points on the Earth’s surface, in degrees.

Returns
Great-circle distance (in m)

If the input is an array (or list) of coordinates, an array of distances is returned

seislib.utils._utils.load_pickle(path)

Loads a .pickle file

Parameters
pathstr

Absolute path to the file

Returns
Object contained in the .pickle file
seislib.utils._utils.next_power_of_2(x)

Closest power of two larger or equal to x

Parameters
xint
Returns
int
seislib.utils._utils.pearson_corrcoef(v1, v2)

Pearson coerrelation coefficient between two vectors

Parameters
v1, v2lists or ndarrays (n,)
Returns
corrcoefffloat

Pearson correlation coefficient

pvaluefloat
Raises
ValueError

If v1 and v2 have different shapes

Notes

This function calls scipy.stats.pearsonr, bu handles the presence of nan values.

seislib.utils._utils.remove_file(file)

Removes a file from disk, handling eventual exceptions

Parameters
filestr

Absolute path to the file to be removed

seislib.utils._utils.resample(x, fs)

If the input streams (or traces) have different sampling rates, the one characterized by the largest sampling rate is downsampled to the sampling rate of the other stream (or trace).

Parameters
st1, st2obspy.Stream, obspy.Trace
Returns
st1, st2obspy.Stream, obspy.Trace

Obspy stream or trace depending on the input. The input is permanently modified

seislib.utils._utils.rotate(r, t, dphi)

Rotation of radial and transverse component of the seismogram by a specified angle, following the obspy signs convention.

Parameters
r, tnumpy.ndarray

Radial (r) and transverse (t) components

dphifloat

Angle in degrees

Returns
rnew, tnewnumpy.ndarray

Rotated components

seislib.utils._utils.rotate_stream(st, **kwargs)

The method calls the obspy.Stream.rotate method, which sometimes runs into errors if differences are present among the starttimes and/or endtimes of the traces constituting the stream. These are prevented by slicing the stream to a common time window and (if necessary) interpolating it so as to avoid sub-sample differences.

Parameters
stobspy.Stream
**kwargs

Optional arguments passed to obspy.Stream.rotate

Returns
stobspy.Stream

Rotated copy of the input Stream

seislib.utils._utils.running_mean(x, N)

Moving average

Parameters
xndarray of shape (m,)

Data vector

Nint

Controls the extent of the smoothing (larger values correspond to larger smoothing)

Returns
runmeanndarray of shape (m,)

Smoothed input

Notes

This is a simple implementation of a moving average. More sofisticated functions can be found, e.g., in scipy.signal.savgol_filter or scipy.ndimage.filters.uniform_filter1d

seislib.utils._utils.save_pickle(file, obj)

Saves an object to a .pickle file

Parameters
filestr

Absolute path to the resulting file

objpython object

Object to be saved (see documentation on the pickle module to know more on which Python objects can be stored into .pickle files)

seislib.utils._utils.scatter_to_mesh(lats, lons, c, mesh, method='linear')

Translates scattered data into a seislib mesh

Parameters
lats, lonsndarray (n,)

Coordinates of the scattered data

cndarray (n,)

Values of the scattered data

meshndarray (m, 4)

mesh in SeisLib’s format, where the four columns correspond to the boundaries of each pixel, i.e., lat1, lat2, lon1, lon2

method{‘linear’, ‘nearest’, ‘cubic’}

Interpolation method. Supported: ‘linear’ (default), ‘nearest’, and ‘cubic’. The three methods call LinearNDInterpolator, NearestNDInterpolator, and CloughTocher2DInterpolator of the scipy.interpolate module, respectively.

Returns
1-D ndarray containing the c values interpolated on mesh
Raises
NotImplementedError

If method is not supported

seislib.utils._utils.zeropad(tr, starttime, endtime)

Zeropads an obspy.Trace so as to cover the time window specified by starttime’and endtime

Parameters
trobspy.Trace
starttime, endtimeobspy.UTCDateTime
Returns
traceobspy.Trace

Zeropadded copy of the input trace.