hyperspy._signals.signal1d module

class hyperspy._signals.signal1d.LazySignal1D(*args, **kwargs)

Bases: hyperspy._signals.lazy.LazySignal, hyperspy._signals.signal1d.Signal1D

class hyperspy._signals.signal1d.Signal1D(*args, **kwargs)

Bases: hyperspy.signal.BaseSignal, hyperspy._signals.common_signal1d.CommonSignal1D

align1D(start=None, end=None, reference_indices=None, max_shift=None, interpolate=True, number_of_interpolation_points=5, interpolation_method='linear', crop=True, expand=False, fill_value=nan, also_align=None, mask=None, show_progressbar=None)

Estimate the shifts in the signal axis using cross-correlation and use the estimation to align the data in place. This method can only estimate the shift by comparing unidimensional features that should not change the position.

To decrease memory usage, time of computation and improve accuracy it is convenient to select the feature of interest setting the start and end keywords. By default interpolation is used to obtain subpixel precision.

Parameters
  • end (start,) – The limits of the interval. If int they are taken as the axis index. If float they are taken as the axis value.

  • reference_indices (tuple of ints or None) – Defines the coordinates of the spectrum that will be used as eference. If None the spectrum at the current coordinates is used for this purpose.

  • max_shift (int) – “Saturation limit” for the shift.

  • interpolate (bool) – If True, interpolation is used to provide sub-pixel accuracy.

  • number_of_interpolation_points (int) – Number of interpolation points. Warning: making this number too big can saturate the memory

  • interpolation_method (str or int) – Specifies the kind of interpolation as a string (‘linear’, ‘nearest’, ‘zero’, ‘slinear’, ‘quadratic, ‘cubic’) or as an integer specifying the order of the spline interpolator to use.

  • crop (bool) – If True automatically crop the signal axis at both ends if needed.

  • expand (bool) – If True, the data will be expanded to fit all data after alignment. Overrides crop.

  • fill_value (float) – If crop is False fill the data outside of the original interval with the given value where needed.

  • also_align (list of signals, None) – A list of BaseSignal instances that has exactly the same dimensions as this one and that will be aligned using the shift map estimated using the this signal.

  • mask (BaseSignal or bool data type.) – It must have signal_dimension = 0 and navigation_shape equal to the current signal. Where mask is True the shift is not computed and set to nan.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

Returns

Return type

An array with the result of the estimation.

Raises

SignalDimensionError – If the signal dimension is not 1.

See also

None()

calibrate(display=True, toolkit=None)

Calibrate the spectral dimension using a gui. It displays a window where the new calibration can be set by:

  • setting the offset, units and scale directly

  • selecting a range by dragging the mouse on the spectrum figure and setting the new values for the given range limits

Parameters
  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

Notes

For this method to work the output_dimension must be 1.

Raises

SignalDimensionError – If the signal dimension is not 1.

create_model(dictionary=None)

Create a model for the current data.

Returns

model

Return type

Model1D instance.

crop_signal1D(*args, **kwargs)

Crop in place the spectral dimension.

Parameters

righ_value (left_value,) – If int the values are taken as indices. If float they are converted to indices using the spectral axis calibration. If left_value is None crops from the beginning of the axis. If right_value is None crops up to the end of the axis. If both are None the interactive cropping interface is activated enabling cropping the spectrum using a span selector in the signal plot.

Raises

SignalDimensionError – If the signal dimension is not 1.

estimate_peak_width(factor=0.5, window=None, return_interval=False, parallel=None, show_progressbar=None)

Estimate the width of the highest intensity of peak of the spectra at a given fraction of its maximum.

It can be used with asymmetric peaks. For accurate results any background must be previously substracted. The estimation is performed by interpolation using cubic splines.

Parameters
  • factor (0 < float < 1) – The default, 0.5, estimates the FWHM.

  • window (None or float) – The size of the window centred at the peak maximum used to perform the estimation. The window size must be chosen with care: if it is narrower than the width of the peak at some positions or if it is so wide that it includes other more intense peaks this method cannot compute the width and a NaN is stored instead.

  • return_interval (bool) – If True, returns 2 extra signals with the positions of the desired height fraction at the left and right of the peak.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

Returns

  • width or [width, left, right], depending on the value of

  • return_interval.

estimate_shift1D(start=None, end=None, reference_indices=None, max_shift=None, interpolate=True, number_of_interpolation_points=5, mask=None, show_progressbar=None, parallel=None)

Estimate the shifts in the current signal axis using cross-correlation. This method can only estimate the shift by comparing unidimensional features that should not change the position in the signal axis. To decrease the memory usage, the time of computation and the accuracy of the results it is convenient to select the feature of interest providing sensible values for start and end. By default interpolation is used to obtain subpixel precision.

Parameters
  • end (start,) – The limits of the interval. If int they are taken as the axis index. If float they are taken as the axis value.

  • reference_indices (tuple of ints or None) – Defines the coordinates of the spectrum that will be used as eference. If None the spectrum at the current coordinates is used for this purpose.

  • max_shift (int) – “Saturation limit” for the shift.

  • interpolate (bool) – If True, interpolation is used to provide sub-pixel accuracy.

  • number_of_interpolation_points (int) – Number of interpolation points. Warning: making this number too big can saturate the memory

  • mask (BaseSignal of bool.) – It must have signal_dimension = 0 and navigation_shape equal to the current signal. Where mask is True the shift is not computed and set to nan.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

Returns

Return type

An array with the result of the estimation in the axis units. Although the computation is performed in batches if the signal is lazy, the result is computed in memory because it depends on the current state of the axes that could change later on in the workflow.

Raises

SignalDimensionError – If the signal dimension is not 1.

filter_butterworth(cutoff_frequency_ratio=None, type='low', order=2, display=True, toolkit=None)

Butterworth filter in place.

Parameters
  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

Raises

SignalDimensionError – If the signal dimension is not 1.

find_peaks1D_ohaver(xdim=None, slope_thresh=0, amp_thresh=None, subchannel=True, medfilt_radius=5, maxpeakn=30000, peakgroup=10, parallel=None)

Find positive peaks along a 1D Signal. It detects peaks by looking for downward zero-crossings in the first derivative that exceed ‘slope_thresh’.

‘slope_thresh’ and ‘amp_thresh’, control sensitivity: higher values will neglect broad peaks (slope) and smaller features (amp), respectively.

peakgroup is the number of points around the top of the peak that are taken to estimate the peak height. For spikes or very narrow peaks, set peakgroup to 1 or 2; for broad or noisy peaks, make peakgroup larger to reduce the effect of noise.

Parameters
  • slope_thresh (float, optional) – 1st derivative threshold to count the peak; higher values will neglect broader features; default is set to 0.

  • amp_thresh (float, optional) – intensity threshold below which peaks are ignored; higher values will neglect smaller features; default is set to 10% of max(y).

  • medfilt_radius (int, optional) – median filter window to apply to smooth the data (see scipy.signal.medfilt); if 0, no filter will be applied; default is set to 5.

  • peakgroup (int, optional) – number of points around the “top part” of the peak that are taken to estimate the peak height; default is set to 10

  • maxpeakn (int, optional) – number of maximum detectable peaks; default is set to 5000.

  • subchannel (bool, optional) – default is set to True.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

Returns

  • structured array of shape (npeaks) containing fields (‘position’,)

  • ’width’, and ‘height’ for each peak.

Raises

SignalDimensionError – If the signal dimension is not 1.

gaussian_filter(FWHM)

Applies a Gaussian filter in the spectral dimension in place.

Parameters

FWHM (float) – The Full Width at Half Maximum of the gaussian in the spectral axis units

Raises
hanning_taper(side='both', channels=None, offset=0)

Apply a hanning taper to the data in place.

Parameters
  • side ('left', 'right' or 'both') – Specify which side to use.

  • channels (None or int) – The number of channels to taper. If None 5% of the total number of channels are tapered.

  • offset (int) –

Returns

Return type

channels

Raises

SignalDimensionError – If the signal dimension is not 1.

integrate_in_range(signal_range='interactive', display=True, toolkit=None)

Sums the spectrum over an energy range, giving the integrated area. The energy range can either be selected through a GUI or the command line.

Parameters

signal_range (a tuple of this form (l, r) or "interactive") – l and r are the left and right limits of the range. They can be numbers or None, where None indicates the extremes of the interval. If l and r are floats the signal_range will be in axis units (for example eV). If l and r are integers the signal_range will be in index units. When signal_range is “interactive” (default) the range is selected using a GUI.

Returns

integrated_spectrum

Return type

BaseSignal subclass

See also

None()

Examples

Using the GUI

>>> s = hs.signals.Signal1D(range(1000))
>>> s.integrate_in_range() #doctest: +SKIP

Using the CLI

>>> s_int = s.integrate_in_range(signal_range=(560,None))

Selecting a range in the axis units, by specifying the signal range with floats.

>>> s_int = s.integrate_in_range(signal_range=(560.,590.))

Selecting a range using the index, by specifying the signal range with integers.

>>> s_int = s.integrate_in_range(signal_range=(100,120))
interpolate_in_between(start, end, delta=3, parallel=None, show_progressbar=None, **kwargs)

Replace the data in a given range by interpolation. The operation is performed in place.

Parameters
  • end (start,) – The limits of the interval. If int they are taken as the axis index. If float they are taken as the axis value.

  • delta (int or float) – The windows around the (start, end) to use for interpolation

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

  • extra keyword arguments are passed to (All) –

  • See the function documentation (scipy.interpolate.interp1d.) –

  • details. (for) –

Raises

SignalDimensionError – If the signal dimension is not 1.

remove_background(signal_range='interactive', background_type='Power Law', polynomial_order=2, fast=True, zero_fill=False, plot_remainder=True, show_progressbar=None, display=True, toolkit=None)

Remove the background, either in place using a gui or returned as a new spectrum using the command line.

Parameters
  • signal_range ("interactive", tuple of ints or floats, optional) – If this argument is not specified, the signal range has to be selected using a GUI. And the original spectrum will be replaced. If tuple is given, the a spectrum will be returned.

  • background_type (str) – The type of component which should be used to fit the background. Possible components: PowerLaw, Gaussian, Offset, Polynomial If Polynomial is used, the polynomial order can be specified

  • polynomial_order (int, default 2) – Specify the polynomial order if a Polynomial background is used.

  • fast (bool) – If True, perform an approximative estimation of the parameters. If False, the signal is fitted using non-linear least squares afterwards.This is slower compared to the estimation but possibly more accurate.

  • zero_fill (bool) – If True, all spectral channels lower than the lower bound of the fitting range will be set to zero (this is the default behavior of Gatan’s DigitalMicrograph). Setting this value to False allows for inspection of the quality of background fit throughout the pre-fitting region.

  • plot_remainder (bool) – If True, add a (green) line previewing the remainder signal after background removal. This preview is obtained from a Fast calculation so the result may be different if a NLLS calculation is finally performed.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

Examples

Using gui, replaces spectrum s

>>> s = hs.signals.Signal1D(range(1000))
>>> s.remove_background() #doctest: +SKIP

Using command line, returns a spectrum

>>> s1 = s.remove_background(signal_range=(400,450), background_type='PowerLaw')

Using a full model to fit the background

>>> s1 = s.remove_background(signal_range=(400,450), fast=False)
Raises

SignalDimensionError – If the signal dimension is not 1.

shift1D(shift_array, interpolation_method='linear', crop=True, expand=False, fill_value=nan, parallel=None, show_progressbar=None)

Shift the data in place over the signal axis by the amount specified by an array.

Parameters
  • shift_array (numpy array) – An array containing the shifting amount. It must have axes_manager._navigation_shape_in_array shape.

  • interpolation_method (str or int) – Specifies the kind of interpolation as a string (‘linear’, ‘nearest’, ‘zero’, ‘slinear’, ‘quadratic, ‘cubic’) or as an integer specifying the order of the spline interpolator to use.

  • crop (bool) – If True automatically crop the signal axis at both ends if needed.

  • expand (bool) – If True, the data will be expanded to fit all data after alignment. Overrides crop.

  • fill_value (float) – If crop is False fill the data outside of the original interval with the given value where needed.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

Raises

SignalDimensionError – If the signal dimension is not 1.

smooth_lowess(smoothing_parameter=None, number_of_iterations=None, show_progressbar=None, parallel=None, display=True, toolkit=None)

Lowess data smoothing in place. If smoothing_parameter or number_of_iterations are None the method is run in interactive mode.

Parameters
  • smoothing_parameter (float or None) – Between 0 and 1. The fraction of the data used when estimating each y-value.

  • number_of_iterations (int or None) – The number of residual-based reweightings to perform.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

Raises

Notes

This method uses the lowess algorithm from the statsmodels library, which needs to be installed to use this method.

smooth_savitzky_golay(polynomial_order=None, window_length=None, differential_order=0, parallel=None, display=True, toolkit=None)

Apply a Savitzky-Golay filter to the data in place. If polynomial_order or window_length or differential_order are None the method is run in interactive mode.

Parameters
  • polynomial_order (int, optional) – The order of the polynomial used to fit the samples. polyorder must be less than window_length.

  • window_length (int, optional) – The length of the filter window (i.e. the number of coefficients). window_length must be a positive odd integer.

  • differential_order (int, optional) – The order of the derivative to compute. This must be a nonnegative integer. The default is 0, which means to filter the data without differentiating.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

Notes

More information about the filter in scipy.signal.savgol_filter.

smooth_tv(smoothing_parameter=None, show_progressbar=None, parallel=None, display=True, toolkit=None)

Total variation data smoothing in place.

Parameters
  • smoothing_parameter (float or None) – Denoising weight relative to L2 minimization. If None the method is run in interactive mode.

  • show_progressbar (None or bool) – If True, display a progress bar. If None, the default from the preferences settings is used.

  • parallel (None or bool) – If True, perform computation in parallel using multiple cores. If None, the default from the preferences settings is used.

  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

Raises

SignalDimensionError – If the signal dimension is not 1.

spikes_removal_tool(signal_mask=None, navigation_mask=None, display=True, toolkit=None)

Graphical interface to remove spikes from EELS spectra.

Parameters
  • signal_mask (boolean array) – Restricts the operation to the signal locations not marked as True (masked)

  • navigation_mask (boolean array) – Restricts the operation to the navigation locations not marked as True (masked)

  • display (bool) – If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

  • toolkit (str, iterable of strings or None) – If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

See also

None()

hyperspy._signals.signal1d.find_peaks_ohaver(y, x=None, slope_thresh=0.0, amp_thresh=None, medfilt_radius=5, maxpeakn=30000, peakgroup=10, subchannel=True)

Find peaks along a 1D line.

Function to locate the positive peaks in a noisy x-y data set. Detects peaks by looking for downward zero-crossings in the first derivative that exceed ‘slope_thresh’. Returns an array containing position, height, and width of each peak. Sorted by position. ‘slope_thresh’ and ‘amp_thresh’, control sensitivity: higher values will neglect wider peaks (slope) and smaller features (amp), respectively.

Parameters
  • y (array) – 1D input array, e.g. a spectrum

  • x (array (optional)) – 1D array describing the calibration of y (must have same shape as y)

  • slope_thresh (float (optional)) – 1st derivative threshold to count the peak; higher values will neglect broader features; default is set to 0.

  • amp_thresh (float (optional)) – intensity threshold below which peaks are ignored; higher values will neglect smaller features; default is set to 10% of max(y).

  • medfilt_radius (int (optional)) – median filter window to apply to smooth the data (see scipy.signal.medfilt); if 0, no filter will be applied; default is set to 5.

  • peakgroup (int (optional)) – number of points around the “top part” of the peak that are taken to estimate the peak height; for spikes or very narrow peaks, keep PeakGroup=1 or 2; for broad or noisy peaks, make PeakGroup larger to reduce the effect of noise; default is set to 10.

  • maxpeakn (int (optional)) – number of maximum detectable peaks; default is set to 30000.

  • subchannel (bool (optional)) – default is set to True.

Returns

P – contains fields: ‘position’, ‘width’, and ‘height’ for each peak.

Return type

structured array of shape (npeaks)

Examples

>>> x = np.arange(0,50,0.01)
>>> y = np.cos(x)
>>> peaks = find_peaks_ohaver(y, x, 0, 0)

Notes

Original code from T. C. O’Haver, 1995. Version 2 Last revised Oct 27, 2006 Converted to Python by Michael Sarahan, Feb 2011. Revised to handle edges better. MCS, Mar 2011

hyperspy._signals.signal1d.interpolate1D(number_of_interpolation_points, data)