clouddrift.wavelet.wavelet_transform

clouddrift.wavelet.wavelet_transform#

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()