clouddrift.wavelet#

This module provides functions for computing wavelet transforms and time-frequency analyses, notably using generalized Morse wavelets.

The Python code in this module was translated from the MATLAB implementation by J. M. Lilly in the jWavelet module of jLab (http://jmlilly.net/code.html).

Lilly, J. M. (2021), jLab: A data analysis package for Matlab, v.1.7.1, doi:10.5281/zenodo.4547006, http://www.jmlilly.net/software.

jLab is licensed under the Creative Commons Attribution-Noncommercial-ShareAlike License (https://creativecommons.org/licenses/by-nc-sa/4.0/). The code that is directly translated from jLab/jWavelet is licensed under the same license. Any other code that is added to this module and that is specific to Python and not the MATLAB implementation is licensed under CloudDrift’s MIT license.

Functions

morse_amplitude(gamma, beta[, order, ...])

Calculate the amplitude coefficient of the generalized Morse wavelets.

morse_freq(gamma, beta)

Frequency measures for generalized Morse wavelets.

morse_logspace_freq(gamma, beta, length[, ...])

Compute logarithmically-spaced frequencies for generalized Morse wavelets with parameters gamma and beta.

morse_properties(gamma, beta)

Calculate the properties of the demodulated generalized Morse wavelets.

morse_wavelet(length, gamma, beta, ...[, ...])

Compute the generalized Morse wavelets of Olhede and Walden (2002), doi: 10.1109/TSP.2002.804066.

morse_wavelet_transform(x, gamma, beta, ...)

Apply a continuous wavelet transform to an input signal using the generalized Morse wavelets of Olhede and Walden (2002).

wavelet_transform(x, wavelet[, boundary, ...])

Apply a continuous wavelet transform to an input signal using an input wavelet function.

clouddrift.wavelet.morse_amplitude(gamma: ndarray | float, beta: ndarray | float, order: int64 | None = 1, normalization: str | None = 'bandpass') float[source]#

Calculate the amplitude coefficient of the generalized Morse wavelets. By default, the amplitude is calculated such that the maximum of the frequency-domain wavelet is equal to 2, which is the bandpass normalization. Optionally, specify normalization="energy" in order to return the coefficient giving the wavelets unit energies. See Lilly and Olhede (2009), doi doi: 10.1109/TSP.2008.2007607.

Parameters#

gammanp.ndarray or float

Gamma parameter of the wavelets.

betanp.ndarray or float

Beta parameter of the wavelets.

orderint, optional

Order of wavelets, default is 1.

normalizationstr, optional

Normalization for the wavelets. By default it is assumed to be "bandpass" which uses a bandpass normalization, meaning that the FFT of the wavelets have peak value of 2 for all central frequencies radian_frequency. The other option is "energy" which uses the unit energy normalization. In this last case the time-domain wavelet energies np.sum(np.abs(wave)**2) are always unity.

Returns#

ampnp.ndarray or float

The amplitude coefficient of the wavelets.

Examples#

TODO

See Also#

morse_wavelet(), morse_freq(), morse_properties(), morse_logspace_freq().

clouddrift.wavelet.morse_freq(gamma: ndarray | float, beta: ndarray | float) tuple[ndarray] | tuple[float][source]#

Frequency measures for generalized Morse wavelets. This functions calculates three different measures fm, fe, and fi of the frequency of the lowest-order generalized Morse wavelet specified by parameters gamma and beta.

Note that all frequency quantities here are in radian as in cos(f t) and not cyclic as in np.cos(2 np.pi f t).

For beta=0, the corresponding wavelet becomes an analytic lowpass filter, and fm is not defined in the usual way but as the point at which the filter has decayed to one-half of its peak power.

For details see Lilly and Olhede (2009), doi: 10.1109/TSP.2008.2007607.

Parameters#

gammanp.ndarray or float

Gamma parameter of the wavelets.

betanp.ndarray or float

Beta parameter of the wavelets.

Returns#

fmnp.ndarray

The modal or peak frequency.

fenp.ndarray

The energy frequency.

finp.ndarray

The instantaneous frequency at the wavelets’ centers.

Examples#

>>> fm, fe, fi = morse_freq(3, 4)
>>> morse_freq(3, 4)
(array(1.10064242), 1.1025129235952809, 1.1077321674324723)
>>> morse_freq(3, np.array([10, 20, 30]))
(array([1.49380158, 1.88207206, 2.15443469]),
array([1.49421505, 1.88220264, 2.15450116]),
array([1.49543843, 1.88259299, 2.15470024]))
>>> morse_freq(np.array([3, 4, 5]), np.array([10, 20, 30]))
(array([1.49380158, 1.49534878, 1.43096908]),
array([1.49421505, 1.49080278, 1.4262489 ]),
array([1.49543843, 1.48652036, 1.42163583]))
>>> morse_freq(np.array([3, 4, 5]), 10)
(array([1.49380158, 1.25743343, 1.14869835]),
array([1.49421505, 1.25000964, 1.13759731]),
array([1.49543843, 1.24350315, 1.12739747]))

See Also#

morse_wavelet(), morse_amplitude()

clouddrift.wavelet.morse_logspace_freq(gamma: float, beta: float, length: int, highset: tuple[float] | None = (0.1, 3.141592653589793), lowset: tuple[float] | None = (5, 0), density: int | None = 4) ndarray[source]#

Compute logarithmically-spaced frequencies for generalized Morse wavelets with parameters gamma and beta. This is a useful function to obtain the frequencies needed for time-frequency analyses using wavelets. If radian_frequencies is the output, np.log(radian_frequencies) is uniformly spaced, following convention for wavelet analysis. See Lilly (2017), doi: 10.1098/rspa.2016.0776.

Default settings to compute the frequencies can be changed by passing optional arguments lowset, highset, and density. See below.

Parameters#

gammafloat

Gamma parameter of the Morse wavelets.

betafloat

Beta parameter of the Morse wavelets.

lengthint

Length of the Morse wavelets and input signals.

highsettuple of floats, optional.

Tuple of values (eta, high) used for high-frequency cutoff calculation. The highest frequency is set to be the minimum of a specified value and a cutoff frequency based on a Nyquist overlap condition: the highest frequency is the minimum of the specified value high, and the largest frequency for which the wavelet will satisfy the threshold level eta. Here eta be a number between zero and one specifying the ratio of a frequency-domain wavelet at the Nyquist frequency to its peak value. Default is (eta, high) = (0.1, np.pi).

lowsettuple of floats, optional.

Tupe of values (P, low) set used for low-frequency cutoff calculation based on an endpoint overlap condition. The lowest frequency is set such that the lowest-frequency wavelet will reach some number P, called the packing number, times its central window width at the ends of the time series. A choice of P=1 corresponds to roughly 95% of the time-domain wavelet energy being contained within the time series endpoints for a wavelet at the center of the domain. The second value of the tuple is the absolute lowest frequency. Default is (P, low) = (5, 0).

densityint, optional

This optional argument controls the number of points in the returned frequency array. Higher values of density mean more overlap in the frequency domain between transforms. When density=1, the peak of one wavelet is located at the half-power points of the adjacent wavelet. The default density=4 means that four other wavelets will occur between the peak of one wavelet and its half-power point.

Returns#

radian_frequencynp.ndarray

Logarithmically-spaced frequencies in radians cycles per unit time, sorted in descending order.

Examples#

Generate a frequency array for the generalized Morse wavelet with parameters gamma=3 and beta=5 for a time series of length n=1024:

>>> radian_frequency = morse_logspace_freq(3, 5, 1024)
>>> radian_frequency = morse_logspace_freq(3, 5, 1024, highset=(0.2, np.pi), lowset=(5, 0))
>>> radian_frequency = morse_logspace_freq(3, 5, 1024, highset=(0.2, np.pi), lowset=(5, 0), density=10)

See Also#

morse_wavelet(), morse_freq(), morse_properties()

clouddrift.wavelet.morse_properties(gamma: ndarray | float, beta: ndarray | float) tuple[ndarray] | tuple[float][source]#

Calculate the properties of the demodulated generalized Morse wavelets. See Lilly and Olhede (2009), doi: 10.1109/TSP.2008.2007607.

Parameters#

gammanp.ndarray or float

Gamma parameter of the wavelets.

betanp.ndarray or float

Beta parameter of the wavelets.

Returns#

widthnp.ndarray or float

Dimensionless time-domain window width of the wavelets.

skewnp.ndarray or float

Imaginary part of normalized third moment of the time-domain demodulate, or ‘demodulate skewness’.

kurtnp.ndarray or float

Normalized fourth moment of the time-domain demodulate, or ‘demodulate kurtosis’.

Examples#

TODO

See Also#

morse_wavelet(), morse_freq(), morse_amplitude(), morse_logspace_freq().

clouddrift.wavelet.morse_wavelet(length: int, gamma: float, beta: float, radian_frequency: ndarray, order: int | None = 1, normalization: str | None = 'bandpass') tuple[ndarray, ndarray][source]#

Compute the generalized Morse wavelets of Olhede and Walden (2002), doi: 10.1109/TSP.2002.804066.

Parameters#

lengthint

Length of the wavelets.

gammafloat

Gamma parameter of the wavelets.

betafloat

Beta parameter of the wavelets.

radian_frequencynp.ndarray

The radian frequencies at which the Fourier transform of the wavelets reach their maximum amplitudes. radian_frequency is between 0 and 2 * np.pi * 0.5, the normalized Nyquist radian frequency.

orderint, optional

Order of wavelets, default is 1.

normalizationstr, optional

Normalization for the wavelet output. By default it is assumed to be "bandpass" which uses a bandpass normalization, meaning that the FFT of the wavelets have peak value of 2 for all central frequencies radian_frequency. The other option is "energy"``which uses the unit energy normalization. In this last case, the time-domain wavelet energies ``np.sum(np.abs(wave)**2) are always unity.

Returns#

waveletnp.ndarray

Time-domain wavelets with shape (order, radian_frequency, length).

wavelet_fft: np.ndarray

Frequency-domain wavelets with shape (order, radian_frequency, length).

Examples#

Compute a Morse wavelet with gamma parameter 3, beta parameter 4, at radian frequency 0.2 cycles per unit time:

>>> wavelet, wavelet_fft = morse_wavelet(1024, 3, 4, np.array([2*np.pi*0.2]))
>>> np.shape(wavelet)
(1, 1, 1024)

Compute a suite of Morse wavelets with gamma parameter 3, beta parameter 4, up to order 3, at radian frequencies 0.2 and 0.3 cycles per unit time:

>>> wavelet, wavelet_fft = morse_wavelet(1024, 3, 4, np.array([2*np.pi*0.2, 2*np.pi*0.3]), order=3)
>>> np.shape(wavelet)
(3, 2, 1024)

Compute a Morse wavelet specifying an energy normalization : >>> wavelet, wavelet_fft = morse_wavelet(1024, 3, 4, np.array([2*np.pi*0.2]), normalization=”energy”)

Raises#

ValueError

If normalization optional argument is not in [“bandpass”, “energy”]``.

See Also#

wavelet_transform(), morse_wavelet_transform(), morse_freq(), morse_logspace_freq(), morse_amplitude(), morse_properties()

clouddrift.wavelet.morse_wavelet_transform(x: ndarray, gamma: float, beta: float, radian_frequency: ndarray, complex: bool | None = False, order: int | None = 1, normalization: str | None = 'bandpass', boundary: str | None = 'mirror', time_axis: int | None = -1) tuple[ndarray] | ndarray[source]#

Apply a continuous wavelet transform to an input signal using the generalized Morse wavelets of Olhede and Walden (2002). The wavelet transform is normalized differently for complex-valued input than for real-valued input, and this in turns depends on whether the optional argument normalization is set to "bandpass" or "energy" normalizations.

Parameters#

xnp.ndarray

Real- or complex-valued signals. The time axis is assumed to be the last. If not, specify optional argument time_axis.

gammafloat

Gamma parameter of the Morse wavelets.

betafloat

Beta parameter of the Morse wavelets.

radian_frequencynp.ndarray

An array of radian frequencies at which the Fourier transform of the wavelets reach their maximum amplitudes. radian_frequency is typically between 0 and 2 * np.pi * 0.5, the normalized Nyquist radian frequency.

complexboolean, optional

Specify explicitely if the input signal x is a complex signal. Default is False which means that the input is real but that is not explicitely tested by the function. This choice affects the normalization of the outputs and their interpretation. See examples below.

time_axisint, optional

Axis on which the time is defined for input x (default is last, or -1).

normalizationstr, optional

Normalization for the wavelet transforms. By default it is assumed to be "bandpass" which uses a bandpass normalization, meaning that the FFT of the wavelets have peak value of 2 for all central frequencies radian_frequency. However, if the optional argument complex=True is specified, the wavelets will be divided by 2 so that the total variance of the input complex signal is equal to the sum of the variances of the returned analytic (positive) and conjugate analytic (negative) parts. See examples below. The other option is "energy" which uses the unit energy normalization. In this last case, the time-domain wavelet energies np.sum(np.abs(wave)**2) are always unity.

boundarystr, optional

The boundary condition to be imposed at the edges of the input signal x. Allowed values are "mirror", "zeros", and "periodic". Default is "mirror".

orderint, optional

Order of Morse wavelets, default is 1.

Returns#

If the input signal is real as specificied by complex=False:

wtxnp.ndarray

Time-domain wavelet transform of input x with shape ((x shape without time_axis), orders, frequencies, time_axis) but with dimensions of length 1 removed (squeezed).

If the input signal is complex as specificied by complex=True, a tuple is returned:

wtx_pnp.array

Time-domain positive wavelet transform of input x with shape ((x shape without time_axis), frequencies, orders), but with dimensions of length 1 removed (squeezed).

wtx_nnp.array

Time-domain negative wavelet transform of input x with shape ((x shape without time_axis), frequencies, orders), but with dimensions of length 1 removed (squeezed).

Examples#

Apply a wavelet transform with a Morse wavelet with gamma parameter 3, beta parameter 4, at radian frequency 0.2 cycles per unit time:

>>> x = np.random.random(1024)
>>> wtx = morse_wavelet_transform(x, 3, 4, np.array([2*np.pi*0.2]))

Apply a wavelet transform with a Morse wavelet with gamma parameter 3, beta parameter 4, for a complex input signal at radian frequency 0.2 cycles per unit time. This case returns the analytic and conjugate analytic components:

>>> z = np.random.random(1024) + 1j*np.random.random(1024)
>>> wtz_p, wtz_n = morse_wavelet_transform(z, 3, 4, np.array([2*np.pi*0.2]), complex=True)

The same result as above can be otained by applying the Morse transform on the real and imaginary component of z and recombining the results as follows for the “bandpass” normalization: >>> wtz_real = morse_wavelet_transform(np.real(z)), 3, 4, np.array([2*np.pi*0.2])) >>> wtz_imag = morse_wavelet_transform(np.imag(z)), 3, 4, np.array([2*np.pi*0.2])) >>> wtz_p, wtz_n = (wtz_real + 1j*wtz_imag) / 2, (wtz_real - 1j*wtz_imag) / 2

For the “energy” normalization, the analytic and conjugate analytic components are obtained as follows with this alternative method: >>> wtz_real = morse_wavelet_transform(np.real(z)), 3, 4, np.array([2*np.pi*0.2])) >>> wtz_imag = morse_wavelet_transform(np.imag(z)), 3, 4, np.array([2*np.pi*0.2])) >>> wtz_p, wtz_n = (wtz_real + 1j*wtz_imag) / np.sqrt(2), (wtz_real - 1j*wtz_imag) / np.sqrt(2)

The input signal can have an arbitrary number of dimensions but its time_axis must be specified if it is not the last:

>>> x = np.random.random((1024,10,15))
>>> wtx = morse_wavelet_transform(x, 3, 4, np.array([2*np.pi*0.2]), time_axis=0)

The default way to handle the boundary conditions is to mirror the ends points but this can be changed by specifying the chosen boundary method:

>>> x = np.random.random((10,15,1024))
>>> wtx = morse_wavelet_transform(x, 3, 4, np.array([2*np.pi*0.2]), boundary="periodic")

This function can be used to conduct a time-frequency analysis of the input signal by specifying a range of randian frequencies using the morse_logspace_freq function as an example:

>>> x = np.random.random(1024)
>>> gamma = 3
>>> beta = 4
>>> radian_frequency = morse_logspace_freq(gamma, beta, np.shape(x)[0])
>>> wtx = morse_wavelet_transform(x, gamma, beta, radian_frequency)

Raises#

ValueError

If the time axis is outside of the valid range ([-1, np.ndim(x)-1]). If boundary optional argument is not in [“mirror”, “zeros”, “periodic”]``. If normalization optional argument is not in [“bandpass”, “energy”]``.

See Also#

morse_wavelet(), wavelet_transform(), morse_logspace_freq()

clouddrift.wavelet.wavelet_transform(x: ndarray, wavelet: ndarray, boundary: str | None = 'mirror', time_axis: int | None = -1, freq_axis: int | None = -2, order_axis: int | None = -3) ndarray[source]#

Apply a continuous wavelet transform to an input signal using an input wavelet function. Such wavelet can be provided by the function morse_wavelet.

Parameters#

xnp.ndarray

Real- or complex-valued signals.

waveletnp.ndarray

A suite of time-domain wavelets, typically returned by the function morse_wavelet. The length of the time axis of the wavelets must be the last one and matches the length of the time axis of x. The other dimensions (axes) of the wavelets (such as orders and frequencies) are typically organized as orders, frequencies, and time, unless specified by optional arguments freq_axis and order_axis. The normalization of the wavelets is assumed to be “bandpass”, if not, use kwarg normalization=”energy”, see morse_wavelet.

boundarystr, optional

The boundary condition to be imposed at the edges of the input signal x. Allowed values are "mirror", "zeros", and "periodic". Default is "mirror".

time_axisint, optional

Axis on which the time is defined for input x (default is last, or -1). Note that the time axis of the wavelets must be last.

freq_axisint, optional

Axis of wavelet for the frequencies (default is second or 1)

order_axisint, optional

Axis of wavelet for the orders (default is first or 0)

Returns#

wtxnp.ndarray

Time-domain wavelet transform of x with shape ((x shape without time_axis), orders, frequencies, time_axis) but with dimensions of length 1 removed (squeezed).

Examples#

Apply a wavelet transform with a Morse wavelet with gamma parameter 3, beta parameter 4, at radian frequency 0.2 cycles per unit time:

>>> x = np.random.random(1024)
>>> wavelet, _ = morse_wavelet(1024, 3, 4, np.array([2*np.pi*0.2]))
>>> wtx = wavelet_transform(x, wavelet)

The input signal can have an arbitrary number of dimensions but its time_axis must be specified if it is not the last:

>>> x = np.random.random((1024,10,15))
>>> wavelet, _ = morse_wavelet(1024, 3, 4, np.array([2*np.pi*0.2]))
>>> wtx = wavelet_transform(x, wavelet,time_axis=0)

Raises#

ValueError

If the time axis is outside of the valid range ([-1, N-1]). If the shape of time axis is different for input signal and wavelet. If boundary optional argument is not in [“mirror”, “zeros”, “periodic”]``.

See Also#

morse_wavelet(), morse_wavelet_transform(), morse_freq()