twistpy.polarization.TimeFrequencyAnalysis3C#

class TimeFrequencyAnalysis3C(N: Trace, E: Trace, Z: Trace, window: dict, dsfacf: int = 1, dsfact: int = 1, frange: Optional[Tuple[float, float]] = None, trange: Optional[Tuple[UTCDateTime, UTCDateTime]] = None, k: float = 1, verbose: bool = True, timeaxis: str = 'utc')[source]#

Time-frequency domain three-component polarization analysis.

Single-station three-component polarization analysis on time-frequency decomposed signals. The time-frequency representation is obtained via the S-transform (twistpy.utils.s_transform).

Parameters
NTrace

North component seismogram

ETrace

East component seismogram

ZTrace

Vertical component seismogram

windowdict

2-D time-frequency window within which the covariance matrices are averaged:

window = {‘number_of_periods’: float, ‘frequency_extent’: float}
The width of the window in time is frequency-dependent and stretches over ‘number_of_periods’ periods, where the period is defined as 1/f. In frequency, the window stretches over the specified ‘frequency_extent’ in Hz.
dsfactint, default=1

Down-sampling factor of the time-axis prior to the computation of polarization attributes.

dsfacfint, default=1

Down-sampling factor of the frequency-axis prior to the computation of polarization attributes.

Warning

Downsampling of the frequency axis (dsfacf > 1) prevents the accurate computation of the inverse S-transform! Down-sampling should be avoided for filtering applications.

frangetuple, (f_min, f_max)

Limit the analysis to the specified frequency band in Hz

trangetuple, (t_min, t_max)

Limit the analysis to the specified time window. t_min and t_max should either be specified as a tuple of two UTCDateTime objects (if timeaxis=’utc’), or as a tuple of float (if timeaxis=’rel’).

kfloat, default=1.

k-value used in the S-transform, corresponding to a scaling factor that controls the number of oscillations in the Gaussian window. When k increases, the frequency resolution increases, with a corresponding loss of time resolution.

See also

twistpy.utils.s_transform

verbosebool, default=True

Run in verbose mode.

timeaxis:obj:’str’, default=’utc’

Specify whether the time axis of plots is shown in UTC (timeaxis=’utc’) or in seconds relative to the first sample (timeaxis=’rel’).

Attributes
t_pollist of UTCDateTime

Time axis of the computed polarization attributes.

f_polndarray of float

Frequency axis of the computed polarization attributes.

Cndarray of complex128

Complex covariance matrices at each window position

timelist of UTCDateTime

Time axis of the input traces

deltafloat

Sampling interval of the input data in seconds

dop2D numpy.ndarray of float

Degree of polarization estimated at each window position in time and frequency.

Hint

The definition of the degree of polarization follows the one from Samson & Olson (1980): Some comments on the descriptions of the polarization states of waves, Geophysical Journal of the Royal Astronomical Society, https://doi.org/10.1111/j.1365-246X.1980.tb04308.x. It is defined as:

\[P^2=\sum_{j,k=1}^{n}(\lambda_j-\lambda_k)^2/[2(n-1)(\sum_{j=1}^{n}(\lambda_j)^2)]\]

where \(P^2\) is the degree of polarization and \(\lambda_j\) are the eigenvalues of the covariance matrix.

elli2D numpy.ndarray of float

Ellipticity estimated at each window position

inc12D numpy.ndarray of float

Inclination of the major semi-axis of the polarization ellipse estimated at each window position

inc22D numpy.ndarray of float

Inclination of the minor semi-axis of the polarization ellipse estimated at each window position

azi12D numpy.ndarray of float

Azimuth of the major semi-axis of the polarization ellipse estimated at each window position

azi22D numpy.ndarray of float

Azimuth of the minor semi-axis of the polarization ellipse estimated at each window position

Methods

filter([plot_filtered_attributes, suppress, ...])

Filter data based on polarization attributes.

plot([clip, major_semi_axis, show, alpha, ...])

Plot polarization analysis.

polarization_analysis()

Perform polarization analysis.

save(name)

Save the current TimeFrequencyAnalaysis3C object to a file on the disk in the current working directory.

Examples using twistpy.polarization.TimeFrequencyAnalysis3C#

3-C Polarization analysis and filtering in the time-frequency domain

3-C Polarization analysis and filtering in the time-frequency domain

3-C Polarization analysis and filtering: Real data example using a marsquake recorded by InSight

3-C Polarization analysis and filtering: Real data example using a marsquake recorded by InSight