Module reference

Sift Functions

emd.sift.sift(X, sift_thresh=1e-08, max_imfs=None, verbose=None, imf_opts=None, envelope_opts=None, extrema_opts=None)[source]

Compute Intrinsic Mode Functions from an input data vector.

This function implements the original sift algorithm [1].

Parameters:
X : ndarray

1D input array containing the time-series data to be decomposed

sift_thresh : float

The threshold at which the overall sifting process will stop. (Default value = 1e-8)

max_imfs : int

The maximum number of IMFs to compute. (Default value = None)

Returns:
imf: ndarray

2D array [samples x nimfs] containing he Intrisic Mode Functions from the decomposition of X.

Other Parameters:
 
imf_opts : dict

Optional dictionary of keyword options to be passed to emd.get_next_imf

envelope_opts : dict

Optional dictionary of keyword options to be passed to emd.interp_envelope

extrema_opts : dict

Optional dictionary of keyword options to be passed to emd.get_padded_extrema

verbose : {None,’CRITICAL’,’WARNING’,’INFO’,’DEBUG’}

Option to override the EMD logger level for a call to this function.

References

[1]Huang, N. E., Shen, Z., Long, S. R., Wu, M. C., Shih, H. H., Zheng, Q., … Liu, H. H. (1998). The empirical mode decomposition and the Hilbert spectrum for nonlinear and non-stationary time series analysis. Proceedings of the Royal Society of London. Series A: Mathematical, Physical and Engineering Sciences, 454(1971), 903–995. https://doi.org/10.1098/rspa.1998.0193
emd.sift.ensemble_sift(X, nensembles=4, ensemble_noise=0.2, noise_mode='single', nprocesses=1, sift_thresh=1e-08, max_imfs=None, verbose=None, imf_opts=None, envelope_opts=None, extrema_opts=None)[source]

Compute Intrinsic Mode Functions with the ensemble EMD.

This function implements the ensemble empirical model decomposition algorithm defined in [1]. This approach sifts an ensemble of signals with white-noise added and treats the mean IMFs as the result. The resulting IMFs from the ensemble sift resembles a dyadic filter [2].

Parameters:
X : ndarray

1D input array containing the time-series data to be decomposed

nensembles : int

Integer number of different ensembles to compute the sift across.

ensemble_noise : float

Standard deviation of noise to add to each ensemble (Default value = .2)

noise_mode : {‘single’,’flip’}

Flag indicating whether to compute each ensemble with noise once or twice with the noise and sign-flipped noise (Default value = ‘single’)

nprocesses : int

Integer number of parallel processes to compute. Each process computes a single realisation of the total ensemble (Default value = 1)

sift_thresh : float

The threshold at which the overall sifting process will stop. (Default value = 1e-8)

max_imfs : int

The maximum number of IMFs to compute. (Default value = None)

Returns:
imf : ndarray

2D array [samples x nimfs] containing he Intrisic Mode Functions from the decomposition of X.

Other Parameters:
 
imf_opts : dict

Optional dictionary of keyword options to be passed to emd.get_next_imf.

envelope_opts : dict

Optional dictionary of keyword options to be passed to emd.interp_envelope

extrema_opts : dict

Optional dictionary of keyword options to be passed to emd.get_padded_extrema

verbose : {None,’CRITICAL’,’WARNING’,’INFO’,’DEBUG’}

Option to override the EMD logger level for a call to this function.

References

[1]Wu, Z., & Huang, N. E. (2009). Ensemble Empirical Mode Decomposition: A Noise-Assisted Data Analysis Method. Advances in Adaptive Data Analysis, 1(1), 1–41. https://doi.org/10.1142/s1793536909000047
[2]Wu, Z., & Huang, N. E. (2004). A study of the characteristics of white noise using the empirical mode decomposition method. Proceedings of the Royal Society of London. Series A: Mathematical, Physical and Engineering Sciences, 460(2046), 1597–1611. https://doi.org/10.1098/rspa.2003.1221
emd.sift.complete_ensemble_sift(X, nensembles=4, ensemble_noise=0.2, noise_mode='single', nprocesses=1, sift_thresh=1e-08, max_imfs=None, verbose=None, imf_opts=None, envelope_opts=None, extrema_opts=None)[source]

Compute Intrinsic Mode Functions with complete ensemble EMD.

This function implements the complete ensemble empirical model decomposition algorithm defined in [1]. This approach sifts an ensemble of signals with white-noise added taking a single IMF across all ensembles at before moving to the next IMF.

Parameters:
X : ndarray

1D input array containing the time-series data to be decomposed

nensembles : int

Integer number of different ensembles to compute the sift across.

ensemble_noise : float

Standard deviation of noise to add to each ensemble (Default value = .2)

noise_mode : {‘single’,’flip’}

Flag indicating whether to compute each ensemble with noise once or twice with the noise and sign-flipped noise (Default value = ‘single’)

nprocesses : int

Integer number of parallel processes to compute. Each process computes a single realisation of the total ensemble (Default value = 1)

sift_thresh : float

The threshold at which the overall sifting process will stop. (Default value = 1e-8)

max_imfs : int

The maximum number of IMFs to compute. (Default value = None)

Returns:
imf: ndarray

2D array [samples x nimfs] containing he Intrisic Mode Functions from the decomposition of X.

noise: array_like

The Intrisic Mode Functions from the decomposition of X.

Other Parameters:
 
imf_opts : dict

Optional dictionary of keyword options to be passed to emd.get_next_imf.

envelope_opts : dict

Optional dictionary of keyword options to be passed to emd.interp_envelope

extrema_opts : dict

Optional dictionary of keyword options to be passed to emd.get_padded_extrema

verbose : {None,’CRITICAL’,’WARNING’,’INFO’,’DEBUG’}

Option to override the EMD logger level for a call to this function.

References

[1]Torres, M. E., Colominas, M. A., Schlotthauer, G., & Flandrin, P. (2011). A complete ensemble empirical mode decomposition with adaptive noise. In 2011 IEEE International Conference on Acoustics, Speech and Signal Processing (ICASSP). IEEE. https://doi.org/10.1109/icassp.2011.5947265
emd.sift.mask_sift(X, mask_amp=1, mask_amp_mode='ratio_imf', mask_freqs='zc', mask_step_factor=2, ret_mask_freq=False, max_imfs=9, sift_thresh=1e-08, nphases=4, nprocesses=1, verbose=None, imf_opts=None, envelope_opts=None, extrema_opts=None)[source]

Compute Intrinsic Mode Functions using a mask sift.

This function implements a masked sift from a dataset using a set of masking signals to reduce mixing of components between modes [1], multiple masks of different phases can be applied when isolating each IMF [2].

This function can either compute the mask frequencies based on the fastest dynamics in the data (the properties of the first IMF from a standard sift) or apply a pre-specified set of masks.

Parameters:
X : ndarray

1D input array containing the time-series data to be decomposed

mask_amp : float or array_like

Amplitude of mask signals as specified by mask_amp_mode. If float the same value is applied to all IMFs, if an array is passed each value is applied to each IMF in turn (Default value = 1)

mask_amp_mode : {‘abs’,’ratio_imf’,’ratio_sig’}

Method for computing mask amplitude. Either in absolute units (‘abs’), or as a ratio of the standard deviation of the input signal (‘ratio_sig’) or previous imf (‘ratio_imf’) (Default value = ‘ratio_imf’)

mask_freqs : {‘zc’,’if’,float,,array_like}

Define the set of mask frequencies to use. If ‘zc’ or ‘if’ are passed, the frequency of the first mask is taken from either the zero-crossings or instantaneous frequnecy the first IMF of a standard sift on the data. If a float is passed this is taken as the first mask frequency. Subsequent masks are defined by the mask_step_factor. If an array_like vector is passed, the values in the vector will specify the mask frequencies.

mask_step_factor : float

Step in frequency between successive masks (Default value = 2)

mask_type : {‘all’,’sine’,’cosine’}

Which type of masking signal to use. ‘sine’ or ‘cosine’ options return the average of a +ve and -ve flipped wave. ‘all’ applies four masks: sine and cosine with +ve and -ve sign and returns the average of all four.

nphases : int > 0

The number of separate sinusoidal masks to apply for each IMF, the phase of masks are uniformly spread across a 0<=p<2pi range (Default=4).

ret_mask_freq : bool

Boolean flag indicating whether mask frequencies are returned (Default value = False)

max_imfs : int

The maximum number of IMFs to compute. (Default value = None)

sift_thresh : float

The threshold at which the overall sifting process will stop. (Default value = 1e-8)

Returns:
imf : ndarray

2D array [samples x nimfs] containing he Intrisic Mode Functions from the decomposition of X.

mask_freqs : ndarray

1D array of mask frequencies, if ret_mask_freq is set to True.

Other Parameters:
 
imf_opts : dict

Optional dictionary of keyword arguments to be passed to emd.get_next_imf

envelope_opts : dict

Optional dictionary of keyword options to be passed to emd.interp_envelope

extrema_opts : dict

Optional dictionary of keyword options to be passed to emd.get_padded_extrema

verbose : {None,’CRITICAL’,’WARNING’,’INFO’,’DEBUG’}

Option to override the EMD logger level for a call to this function.

Notes

Here are some example mask_sift variants you can run:

A mask sift in which the mask frequencies are determined with zero-crossings and mask amplitudes by a ratio with the amplitude of the previous IMF (note - this is also the default):

>>> imf = emd.sift.mask_sift(X, mask_amp_mode='ratio_imf', mask_freqs='zc')

A mask sift in which the first mask is set at .4 of the sampling rate and subsequent masks found by successive division of this mask_freq by 3:

>>> imf = emd.sift.mask_sift(X, mask_freqs=.4, mask_step_factor=3)

A mask sift using user specified frequencies and amplitudes:

>>> mask_freqs = np.array([.4,.2,.1,.05,.025,0])
>>> mask_amps = np.array([2,2,1,1,.5,.5])
>>> imf = emd.sift.mask_sift(X, mask_freqs=mask_freqs, mask_amp=mask_amps, mask_amp_mode='abs')

References

[1]Ryan Deering, & James F. Kaiser. (2005). The Use of a Masking Signal to Improve Empirical Mode Decomposition. In Proceedings. (ICASSP ’05). IEEE International Conference on Acoustics, Speech, and Signal Processing, 2005. IEEE. https://doi.org/10.1109/icassp.2005.1416051
[2]Tsai, F.-F., Fan, S.-Z., Lin, Y.-S., Huang, N. E., & Yeh, J.-R. (2016). Investigating Power Density and the Degree of Nonlinearity in Intrinsic Components of Anesthesia EEG by the Hilbert-Huang Transform: An Example Using Ketamine and Alfentanil. PLOS ONE, 11(12), e0168108. https://doi.org/10.1371/journal.pone.0168108

Sift Utilities

emd.sift.get_config(siftname='sift')[source]

Return a SiftConfig with default options for a specified sift variant.

Helper function for specifying config objects specifying parameters to be used in a sift. The functions used during the sift areinspected automatically and default values are populated into a nested dictionary which can be modified and used as input to one of the sift functions.

Parameters:
siftname : str

Name of the sift function to find configuration from

Returns:
SiftConfig

A modified dictionary containing the sift specification

Notes

The sift config acts as a nested dictionary which can be modified to specify parameters for different parts of the sift. This is initialised using this function

config = emd.sift.get_config()

The first level of the dictionary contains six sub-dicts configuring different parts of the algorithm:

config[‘sift’] - top level sift options, mostly specific to the particular sift algorithm config[‘imf’] - options for detecting IMFs config[‘envelope’] - options for upper and lower envelope interpolation config[‘extrema’] - options for extrema detection config[‘mag_pad’] - options for y-values of padded extrema at edges config[‘loc_pad’] - options for x-values of padded extrema at edges

Specific values can be modified in the dictionary

config[‘extrema’][‘parabolic_extrema’] = True

or using this shorthand

config[‘imf/env_step_factor’] = 1/3

Finally, the SiftConfig dictionary should be nested before being passed as keyword arguments to a sift function.

imfs = emd.sift.sift(X, **config)

emd.sift.get_next_imf(X, env_step_size=1, max_iters=1000, energy_thresh=None, stop_method='sd', sd_thresh=0.1, rilling_thresh=(0.05, 0.5, 0.05), envelope_opts=None, extrema_opts=None)[source]

Compute the next IMF from a data set.

This is a helper function used within the more general sifting functions.

Parameters:
X : ndarray [nsamples x 1]

1D input array containing the time-series data to be decomposed

env_step_size : float

Scaling of envelope prior to removal at each iteration of sift. The average of the upper and lower envelope is muliplied by this value before being subtracted from the data. Values should be between 0 > x >= 1 (Default value = 1)

max_iters : int > 0

Maximum number of iterations to compute before throwing an error

energy_thresh : float > 0

Threshold for energy difference (in decibels) between IMF and residual to suggest stopping overall sift. (Default is None, recommended value is 50)

stop_method : {‘sd’,’rilling’,’fixed’}

Flag indicating which metric to use to stop sifting and return an IMF.

sd_thresh : float

Used if ‘stop_method’ is ‘sd’. The threshold at which the sift of each IMF will be stopped. (Default value = .1)

rilling_thresh : tuple

Used if ‘stop_method’ is ‘rilling’, needs to contain three values (sd1, sd2, alpha). An evaluation function (E) is defined by dividing the residual by the mode amplitude. The sift continues until E < sd1 for the fraction (1-alpha) of the data, and E < sd2 for the remainder. See section 3.2 of http://perso.ens-lyon.fr/patrick.flandrin/NSIP03.pdf

Returns:
proto_imf : ndarray

1D vector containing the next IMF extracted from X

continue_flag : bool

Boolean indicating whether the sift can be continued beyond this IMF

Other Parameters:
 
envelope_opts : dict

Optional dictionary of keyword arguments to be passed to emd.interp_envelope

extrema_opts : dict

Optional dictionary of keyword options to be passed to emd.get_padded_extrema

emd.sift.get_next_imf_mask(X, z, amp, nphases=4, nprocesses=1, imf_opts=None, envelope_opts=None, extrema_opts=None)[source]

Compute the next IMF from a data set a mask sift.

This is a helper function used within the more general sifting functions.

Parameters:
X : ndarray

1D input array containing the time-series data to be decomposed

z : float

Mask frequency as a proportion of the sampling rate, values between 0->z->.5

amp : float

Mask amplitude

nphases : int > 0

The number of separate sinusoidal masks to apply for each IMF, the phase of masks are uniformly spread across a 0<=p<2pi range (Default=4).

nprocesses : int

Integer number of parallel processes to compute. Each process computes an IMF from the signal plus a mask. nprocesses should be less than or equal to nphases, no additional benefit from setting nprocesses > nphases (Default value = 1)

Returns:
proto_imf : ndarray

1D vector containing the next IMF extracted from X

continue_sift : bool

Boolean indicating whether the sift can be continued beyond this IMF

Other Parameters:
 
imf_opts : dict

Optional dictionary of keyword arguments to be passed to emd.get_next_imf

envelope_opts : dict

Optional dictionary of keyword options to be passed to emd.interp_envelope

extrema_opts : dict

Optional dictionary of keyword options to be passed to emd.get_padded_extrema

emd.sift.interp_envelope(X, mode='upper', interp_method='splrep', extrema_opts=None, ret_extrema=False)[source]

Interpolate the amplitude envelope of a signal.

Parameters:
X : ndarray

Input signal

mode : {‘upper’,’lower’,’combined’}

Flag to set which envelope should be computed (Default value = ‘upper’)

interp_method : {‘splrep’,’pchip’,’mono_pchip’}

Flag to indicate which interpolation method should be used (Default value = ‘splrep’)

Returns:
ndarray

Interpolated amplitude envelope

emd.sift.get_padded_extrema(X, pad_width=2, mode='peaks', parabolic_extrema=False, loc_pad_opts=None, mag_pad_opts=None)[source]

Identify and pad the extrema in a signal.

This function returns a set of extrema from a signal including padded extrema at the edges of the signal. Padding is carried out using numpy.pad.

Parameters:
X : ndarray

Input signal

pad_width : int >= 0

Number of additional extrema to add to the start and end

mode : {‘peaks’, ‘troughs’, ‘abs_peaks’}

Switch between detecting peaks, troughs or peaks in the abs signal

parabolic_extrema : bool

Flag indicating whether extrema positions should be refined by parabolic interpolation

loc_pad_opts : dict

Optional dictionary of options to be passed to np.pad when padding extrema locations

mag_pad_opts : dict

Optional dictionary of options to be passed to np.pad when padding extrema magnitudes

Returns:
locs : ndarray

location of extrema in samples

mags : ndarray

Magnitude of each extrema

Frequency Functions

emd.spectra.frequency_transform(imf, sample_rate, method, smooth_phase=5)[source]

Compute instantaneous phase, frequency and amplitude from a set of IMFs.

Several approaches are implemented from [1] and [2].

Parameters:
imf : ndarray

Input array of IMFs.

sample_rate : float

Sampling frequency of the signal in Hz

method : {‘hilbert’,’quad’,’direct_quad’,’nht’}

The method for computing the frequency stats

smooth_phase : int

Length of window when smoothing the unwrapped phase (Default value = 31)

Returns:
IP : ndarray

Array of instantaneous phase estimates

IF : ndarray

Array of instantaneous frequency estimates

IA : ndarray

Array of instantaneous amplitude estimates

References

[1]Huang, N. E., Shen, Z., Long, S. R., Wu, M. C., Shih, H. H., Zheng, Q., … Liu, H. H. (1998). The empirical mode decomposition and the Hilbert spectrum for nonlinear and non-stationary time series analysis. Proceedings of the Royal Society of London. Series A: Mathematical, Physical and Engineering Sciences, 454(1971), 903–995. https://doi.org/10.1098/rspa.1998.0193
[2]Huang, N. E., Wu, Z., Long, S. R., Arnold, K. C., Chen, X., & Blank, K. (2009). On Instantaneous Frequency. Advances in Adaptive Data Analysis, 1(2), 177–229. https://doi.org/10.1142/s1793536909000096
emd.spectra.quadrature_transform(X)[source]

Compute the quadrature transform on a set of time-series.

This algorithm is defined in equation 34 of [1]. The return is a complex array with the input data as the real part and the quadrature transform as the imaginary part.

Parameters:
X : ndarray

Array containing time-series to transform

Returns:
quad_signal : ndarray

Complex valued array containing the quadrature transformed signal

References

[1]Huang, N. E., Wu, Z., Long, S. R., Arnold, K. C., Chen, X., & Blank, K. (2009). On Instantaneous Frequency. Advances in Adaptive Data Analysis, 1(2), 177–229. https://doi.org/10.1142/s1793536909000096
emd.spectra.phase_from_complex_signal(complex_signal, smoothing=None, ret_phase='wrapped', phase_jump='ascending')[source]

Compute the instantaneous phase from a complex signal.

The complex input may be obtained from either the Hilbert Transform or by Direct Quadrature.

Parameters:
complex_signal : complex ndarray

Complex valued input array

smoothing : int

Integer window length used in phase smoothing (Default value = None)

ret_phase : {‘wrapped’,’unwrapped’}

Flag indicating whether to return the wrapped or unwrapped phase (Default value = ‘wrapped’)

phase_jump : {‘ascending’,’peak’,’descending’,’trough’}

Flag indicating where in the cycle the phase jump should be (Default value = ‘ascending’)

Returns:
IP : ndarray

Array of instantaneous phase values

emd.spectra.freq_from_phase(iphase, sample_rate)[source]

Compute the instantaneous frequency from the instantaneous phase.

Parameters:
iphase : ndarray

Input array containing the unwrapped instantaneous phase time-course

sample_rate : float

The sampling frequency of the data

Returns:
IF : ndarray

Array containing the instantaneous frequencies

emd.spectra.phase_from_freq(ifrequency, sample_rate, phase_start=-3.141592653589793)[source]

Compute the instantaneous phase of a signal from its instantaneous frequency.

Parameters:
ifrequency : ndarray

Input array containing the instantaneous frequencies of a signal

sample_rate : float

The sampling frequency of the data

phase_start : float

Start value of the phase output (Default value = -np.pi)

Returns:
IP : ndarray

The instantaneous phase of the signal

emd.spectra.direct_quadrature(fm)[source]

Compute the quadrature transform on a set of time-series.

This algorithm is defined in equation 35 of [1].

Section 3.2 of ‘on instantaneous frequency’

THIS IS IN DEVELOPMENT

Parameters:
fm : ndarray

Input signal containing a frequency-modulated signal

References

[1]Huang, N. E., Wu, Z., Long, S. R., Arnold, K. C., Chen, X., & Blank, K. (2009). On Instantaneous Frequency. Advances in Adaptive Data Analysis, 1(2), 177–229. https://doi.org/10.1142/s1793536909000096
emd.spectra.phase_angle(fm)[source]

Compute the phase angle of a set of time-series.

This algorithm is defined in equation 35 of [1].

THIS IS IN DEVELOPMENT

Parameters:
X : ndarray

Array containing time-series to transform

Returns:
quad_signal : ndarray

Complex valued array containing the quadrature transformed signal

References

[1]Huang, N. E., Wu, Z., Long, S. R., Arnold, K. C., Chen, X., & Blank, K. (2009). On Instantaneous Frequency. Advances in Adaptive Data Analysis, 1(2), 177–229. https://doi.org/10.1142/s1793536909000096

Spectrum Functions

emd.spectra.hilberthuang_1d(infr, inam, freq_edges, mode='energy')[source]

Compute the Hilbert-Huang transform with no time-dimension.

The Hilbert-Huang transform instataneous frequency and instantaneous amplitude time-series and represents the energy of a signal across frequency [1].

Parameters:
infr : ndarray

2D first level instantaneous frequencies

inam : ndarray

2D first level instantaneous amplitudes

freq_edges : ndarray

Vector of frequency bins for carrier frequencies

mode : {‘energy’,’amplitude’}

Flag indicating whether to sum the energy or amplitudes (Default value = ‘energy’)

Returns:
specs : ndarray

2D array containing Hilbert-Huang Spectrum [ frequencies x imfs ]

References

[1]Huang, N. E., Shen, Z., Long, S. R., Wu, M. C., Shih, H. H., Zheng, Q., … Liu, H. H. (1998). The empirical mode decomposition and the Hilbert spectrum for nonlinear and non-stationary time series analysis. Proceedings of the Royal Society of London. Series A: Mathematical, Physical and Engineering Sciences, 454(1971), 903–995. https://doi.org/10.1098/rspa.1998.0193
emd.spectra.hilberthuang(infr, inam, freq_edges, mode='energy', return_sparse=False)[source]

Compute a Hilbert-Huang transform.

The Hilbert-Huang transform instataneous frequency and instantaneous amplitude time-series and represents the energy of a signal across time and frequency [1].

Parameters:
infr : ndarray

2D first level instantaneous frequencies

inam : ndarray

2D first level instantaneous amplitudes

freq_edges : ndarray

Vector of frequency bins for carrier frequencies

mode : {‘energy’,’amplitude’}

Flag indicating whether to sum the energy or amplitudes (Default value = ‘energy’)

return_sparse : bool

Flag indicating whether to return the full or sparse form(Default value = True)

Returns:
hht : ndarray

2D array containing the Hilbert-Huang Transform

Notes

If return_sparse is set to True the returned array is a sparse matrix in COOrdinate form (scipy.sparse.coo_matrix), also known as ‘ijv’ or ‘triplet’ form. This is much more memory efficient than the full form but may not behave as expected in functions expecting full arrays.

References

[1]Huang, N. E., Shen, Z., Long, S. R., Wu, M. C., Shih, H. H., Zheng, Q., … Liu, H. H. (1998). The empirical mode decomposition and the Hilbert spectrum for nonlinear and non-stationary time series analysis. Proceedings of the Royal Society of London. Series A: Mathematical, Physical and Engineering Sciences, 454(1971), 903–995. https://doi.org/10.1098/rspa.1998.0193
emd.spectra.holospectrum(infr, infr2, inam2, freq_edges, freq_edges2, mode='energy', squash_time='sum')[source]

Compute a Holospectrum.

Holospectra are computed from the first and second layer frequecy statistics of a dataset. The Holospectrum represents the energy of a signal across time, carrier frequency and amplitude-modulation frequency [1].

Parameters:
infr : ndarray

2D first level instantaneous frequencies

infr2 : ndarray

3D second level instantaneous frequencies

inam2 : ndarray

3D second level instantaneous amplitudes

freq_edges : ndarray

Vector of frequency bins for carrier frequencies

freq_edges2 :

Vector of frequency bins for amplitude-modulation frequencies

mode : {‘energy’,’amplitude’}

Flag indicating whether to sum the energy or amplitudes (Default value = ‘energy’)

squash_time : {‘sum’,’mean’,False}

Flag indicating whether to marginalise over the time dimension (Default value = ‘sum’)

Returns:
holo : ndarray

Holospectrum of input data.

Notes

Output will be a 3D [samples x am_freq x carrier_freq] array if squash_time is False and a 2D [ am_freq x carrier_freq ] array if squash_time is true.

References

[1]Huang, N. E., Hu, K., Yang, A. C. C., Chang, H.-C., Jia, D., Liang, W.-K., … Wu, Z. (2016). On Holo-Hilbert spectral analysis: a full informational spectral representation for nonlinear and non-stationary data. Philosophical Transactions of the Royal Society A: Mathematical, Physical and Engineering Sciences, 374(2065), 20150206. https://doi.org/10.1098/rsta.2015.0206

Spectrum Utilities

emd.spectra.define_hist_bins(data_min, data_max, nbins, scale='linear')[source]

Define the bin edges and centre values for use in a histogram.

Parameters:
data_min : float

Value for minimum edge

data_max : float

Value for maximum edge

nbins : int

Number of bins to create

scale : {‘linear’,’log’}

Flag indicating whether to use a linear or log spacing between bins (Default value = ‘linear’)

Returns:
edges : ndarray

1D array of bin edges

centres : ndarray

1D array of bin centres

Notes

>> edges,centres = emd.spectra.define_hist_bins( 1, 5, 3 ) >> print(edges) [1. 2. 3. 4. 5.] >> print(centres) [1.5 2.5 3.5 4.5]

emd.spectra.define_hist_bins_from_data(X, nbins=None, mode='sqrt', scale='linear')[source]

Find the bin edges and centre frequencies for use in a histogram.

If nbins is defined, mode is ignored

Parameters:
X : ndarray

Dataset whose summary stats will define the histogram

nbins : int

Number of bins to create, if undefined this is derived from the data (Default value = None)

mode : {‘sqrt’}

Method for deriving number of bins if nbins is undefined (Default value = ‘sqrt’)

scale : {‘linear’,’log’}

(Default value = ‘linear’)

Returns:
edges : ndarray

1D array of bin edges

centres : ndarray

1D array of bin centres

Cycle Analysis

emd.cycles.get_cycle_vector(phase, return_good=True, mask=None, imf=None, phase_step=4.71238898038469, phase_edge=0.2617993877991494)[source]

Identify cycles within a instantaneous phase time-course.

Cycles are identified by large phase jumps and can optionally be tested to remove ‘bad’ cycles by criteria in Notes.

Parameters:
phase : ndarray

Input vector of Instantaneous Phase values

return_good : bool

Boolean indicating whether ‘bad’ cycles should be removed (Default value = True)

mask : ndarray

Vector of mask values that should be ignored (Default value = None)

imf : ndarray

Optional array of IMFs to used for control point identification when identifying good/bad cycles (Default value = None)

phase_step : scalar

Minimum value in the differential of the wrapped phase to identify a cycle transition (Default value = 1.5*np.pi)

phase_edge : scalar

Maximum distance from 0 or 2pi for the first and last phase value in a good cycle. Only used when return_good is True (Default value = np.pi/12)

Returns:
ndarray

Vector of integers indexing the location of each cycle

Notes

Good cycles are those with 1 : A strictly positively increasing phase 2 : A phase starting within phase_step of zero (ie 0 < x < phase_edge) 3 : A phase ending within phase_step of 2pi (is 2pi-phase_edge < x < 2pi) 4 : A set of 4 unique control points (asc-zero, peak, desc-zero & trough)

Good cycles can be idenfied with: >> good_cycles = emd.utils.get_cycle_vector( phase )

The total number of cycles is then >> good_cycles.max()

Indices where good cycles is zero do not contain a valid cycle bad_segments = good_cycles>0

A single cycle can be isolated by matching its index, eg for the 5th cycle cycle_5_inds = good_cycles==5

emd.cycles.get_cycle_stat(cycles, values, mode='cycle', out=None, func=<function mean>)[source]

Compute the average of a set of observations for each cycle.

Parameters:
cycles : ndarray

array whose content index cycle locations

values : ndarray

array of observations to average within each cycle

mode : {‘compressed’,’full’}

Flag to indicate whether to return a single value per cycle or the average values filled within a vector of the same size as values (Default value = ‘compressed’)

func : function

Function to call on the data in values for each cycle (Default np.mean). This can be any function, built-in or user defined, that processes a single vector of data returning a single value.

Returns:
ndarray

Array containing the cycle-averaged values

emd.cycles.get_control_points(x, cycles, interp=False, mode='cycle')[source]

Identify sets of control points from identified cycles. The control points are the ascending zero, peak, descending zero & trough.

Parameters:
x : ndarray

Input array of oscillatory data

good_cycles : ndarray

array whose content index cycle locations

Returns:
ndarray

The control points for each cycle in x

emd.cycles.phase_align(ip, x, cycles=None, npoints=48, interp_kind='linear', ii=None, mode='cycle')[source]

Compute phase alignment of a vector of observed values across a set of cycles.

Parameters:
ip : ndarray

Input array of Instantaneous Phase values to base alignment on

x : ndarray

Input array of observed values to phase align

cycles : ndarray (optional)

Optional set of cycles within IP to use (Default value = None)

npoints : int

Number of points in the phase cycle to align to (Default = 48)

interp_kind : {‘linear’,’nearest’,’zero’,’slinear’, ‘quadratic’,’cubic’,’previous’, ‘next’}

Type of interpolation to perform. Argument is passed onto scipy.interpolate.interp1d. (Default = ‘linear’)

Returns:
ndarray :

array containing the phase aligned observations

emd.cycles.bin_by_phase(ip, x, nbins=24, weights=None, variance_metric='variance', bin_edges=None)[source]

Compute distribution of x by phase-bins in the Instantaneous Frequency.

Parameters:
ip : ndarray

Input vector of instataneous phase values

x : ndarray

Input array of values to be binned, first dimension much match length of IP

nbins : integer

number of phase bins to define (Default value = 24)

weights : ndarray (optional)

Optional set of linear weights to apply before averaging (Default value = None)

variance_metric : {‘variance’,’std’,’sem’}

Flag to select whether the variance, standard deviation or standard error of the mean in computed across cycles (Default value = ‘variance’)

bin_edges : ndarray (optional)

Optional set of bin edges to override automatic bin specification (Default value = None)

Returns:
avg : ndarray

Vector containing the average across cycles as a function of phase

var : ndarray

Vector containing the selected variance metric across cycles as a function of phase

bin_centres : ndarray

Vector of bin centres

emd.cycles.mean_vector(IP, X, mask=None)[source]

Compute the mean vector of a set of values wrapped around the unit circle.

Parameters:
IP : ndarray

Instantaneous Phase values

X : ndarray

Observations corresponding to IP values

mask :

(Default value = None)

Returns:
mv : ndarray

Set of mean vectors

emd.cycles.kdt_match(x, y, K=15, distance_upper_bound=inf)[source]

Find unique nearest-neighbours between two n-dimensional feature sets. Useful for matching two sets of cycles on one or more features (ie amplitude and average frequency).

Rows in x are matched to rows in y. As such - it is good to have (many) more rows in y than x if possible.

This uses a k-dimensional tree to query for the K nearest neighbours and returns the closest unique neighbour. If no unique match is found - the row is not returned. Increasing K will find more matches but allow matches between more distant observations.

Not advisable for use with more than a handful of features.

Parameters:
x : ndarray

[ num observations x num features ] array to match to

y : ndarray

[ num observations x num features ] array of potential matches

K : int

number of potential nearest-neigbours to query

Returns:
ndarray

indices of matched observations in x

ndarray

indices of matched observations in y

Package Utilities

emd.support.get_install_dir()[source]

Get directory path of currently installed & imported emd.

emd.support.get_installed_version()[source]

Read version of currently installed & imported emd.

Version is determined according to local setup.py. If a user has made local changes this version may not be exactly the same as the online package.

emd.logger.set_up(prefix='', log_file='', level=None, console_format=None)[source]

Initialise the EMD module logger.

Parameters:
prefix : str

Optional prefix to attach to logger output

log_file : str

Optional path to a log file to record logger output

level : {‘CRITICAL’, ‘WARNING’, ‘INFO’, ‘DEBUG’}

String indicating initial logging level

console_format : str

Formatting string for console logging.