hyperspy._signals.signal2d module

class hyperspy._signals.signal2d.LazySignal2D(*args, **kwargs)

Bases: hyperspy._signals.lazy.LazySignal, hyperspy._signals.signal2d.Signal2D

Create a Signal from a numpy array.

Parameters
  • data (numpy.ndarray) – The signal data. It can be an array of any dimensions.

  • axes (dict, optional) – Dictionary to define the axes (see the documentation of the AxesManager class for more details).

  • attributes (dict, optional) – A dictionary whose items are stored as attributes.

  • metadata (dict, optional) – A dictionary containing a set of parameters that will to stores in the metadata attribute. Some parameters might be mandatory in some cases.

  • original_metadata (dict, optional) – A dictionary containing a set of parameters that will to stores in the original_metadata attribute. It typically contains all the parameters that has been imported from the original data file.

class hyperspy._signals.signal2d.Signal2D(*args, **kw)

Bases: hyperspy.signal.BaseSignal, hyperspy._signals.common_signal2d.CommonSignal2D

Create a Signal from a numpy array.

Parameters
  • data (numpy.ndarray) – The signal data. It can be an array of any dimensions.

  • axes (dict, optional) – Dictionary to define the axes (see the documentation of the AxesManager class for more details).

  • attributes (dict, optional) – A dictionary whose items are stored as attributes.

  • metadata (dict, optional) – A dictionary containing a set of parameters that will to stores in the metadata attribute. Some parameters might be mandatory in some cases.

  • original_metadata (dict, optional) – A dictionary containing a set of parameters that will to stores in the original_metadata attribute. It typically contains all the parameters that has been imported from the original data file.

add_ramp(ramp_x, ramp_y, offset=0)

Add a linear ramp to the signal.

Parameters
  • ramp_x (float) – Slope of the ramp in x-direction.

  • ramp_y (float) – Slope of the ramp in y-direction.

  • offset (float, optional) – Offset of the ramp at the signal fulcrum.

Notes

The fulcrum of the linear ramp is at the origin and the slopes are given in units of the axis with the according scale taken into account. Both are available via the axes_manager of the signal.

align2D(crop=True, fill_value=nan, shifts=None, expand=False, interpolation_order=1, show_progressbar=None, parallel=None, max_workers=None, **kwargs)

Align the images in-place using scipy.ndimage.shift().

The images can be aligned using either user-provided shifts or by first estimating the shifts.

See estimate_shift2D() for more details on estimating image shifts.

Parameters
  • crop (bool) – If True, the data will be cropped not to include regions with missing data

  • fill_value (int, float, nan) – The areas with missing data are filled with the given value. Default is nan.

  • shifts (None or list of tuples) – If None the shifts are estimated using estimate_shift2D().

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

  • interpolation_order (int, default 1.) – The order of the spline interpolation. Default is 1, linear 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 multithreading. If None, the default from the preferences settings is used. The number of threads is controlled by the max_workers argument.

  • max_workers (None or int) – Maximum number of threads used when parallel=True. If None, defaults to min(32, os.cpu_count()).

  • **kwargs – Keyword arguments passed to estimate_shift2D()

Returns

shifts – The estimated shifts are returned only if shifts is None

Return type

np.array

create_model(dictionary=None)

Create a model for the current signal

Parameters

dictionary ({None, dict}, optional) – A dictionary to be used to recreate a model. Usually generated using hyperspy.model.as_dictionary()

Returns

Return type

A Model class

crop_image(top=None, bottom=None, left=None, right=None, convert_units=False)

Crops an image in place.

Parameters
  • top ({int | float}) – If int the values are taken as indices. If float the values are converted to indices.

  • bottom ({int | float}) – If int the values are taken as indices. If float the values are converted to indices.

  • left ({int | float}) – If int the values are taken as indices. If float the values are converted to indices.

  • right ({int | float}) – If int the values are taken as indices. If float the values are converted to indices.

  • convert_units (bool) – Default is False If True, convert the signal units using the ‘convert_to_units’ method of the axes_manager. If False, does nothing.

See also

crop

estimate_shift2D(reference='current', correlation_threshold=None, chunk_size=30, roi=None, normalize_corr=False, sobel=True, medfilter=True, hanning=True, plot=False, dtype='float', show_progressbar=None, sub_pixel_factor=1)

Estimate the shifts in an image using phase correlation.

This method can only estimate the shift by comparing bi-dimensional features that should not change position between frames. To decrease the memory usage, the time of computation and the accuracy of the results it is convenient to select a region of interest by setting the roi argument.

Parameters
  • reference ({'current', 'cascade' ,'stat'}) – If ‘current’ (default) the image at the current coordinates is taken as reference. If ‘cascade’ each image is aligned with the previous one. If ‘stat’ the translation of every image with all the rest is estimated and by performing statistical analysis on the result the translation is estimated.

  • correlation_threshold ({None, 'auto', float}) – This parameter is only relevant when reference=’stat’. If float, the shift estimations with a maximum correlation value lower than the given value are not used to compute the estimated shifts. If ‘auto’ the threshold is calculated automatically as the minimum maximum correlation value of the automatically selected reference image.

  • chunk_size ({None, int}) – If int and reference=’stat’ the number of images used as reference are limited to the given value.

  • roi (tuple of ints or floats (left, right, top, bottom)) – Define the region of interest. If int(float) the position is given axis index(value). Note that ROIs can be used in place of a tuple.

  • normalize_corr (bool, default False) – If True, use phase correlation to align the images, otherwise use cross correlation.

  • sobel (bool, default True) – Apply a Sobel filter for edge enhancement

  • medfilter (bool, default True) – Apply a median filter for noise reduction

  • hanning (bool, default True) – Apply a 2D hanning filter

  • plot (bool or 'reuse') – If True plots the images after applying the filters and the phase correlation. If ‘reuse’, it will also plot the images, but it will only use one figure, and continuously update the images in that figure as it progresses through the stack.

  • dtype (str or dtype) – Typecode or data-type in which the calculations must be performed.

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

  • sub_pixel_factor (float) – Estimate shifts with a sub-pixel accuracy of 1/sub_pixel_factor parts of a pixel. Default is 1, i.e. no sub-pixel accuracy.

Returns

shifts – List of estimated shifts

Return type

list of array

Notes

The statistical analysis approach to the translation estimation when using reference='stat' roughly follows [Schaffer2004]. If you use it please cite their article.

References

Schaffer2004

Schaffer, Bernhard, Werner Grogger, and Gerald Kothleitner. “Automated Spatial Drift Correction for EFTEM Image Series.” Ultramicroscopy 102, no. 1 (December 2004): 27–36.

See also

find_peaks(method='local_max', interactive=True, current_index=False, show_progressbar=None, parallel=None, max_workers=None, display=True, toolkit=None, **kwargs)

Find peaks in a 2D signal.

Function to locate the positive peaks in an image using various, user specified, methods. Returns a structured array containing the peak positions.

Parameters
  • method (str) –

    Select peak finding algorithm to implement. Available methods are:

    • ’local_max’ - simple local maximum search using the skimage.feature.peak_local_max() function

    • ’max’ - simple local maximum search using the find_peaks_max().

    • ’minmax’ - finds peaks by comparing maximum filter results with minimum filter, calculates centers of mass. See the find_peaks_minmax() function.

    • ’zaefferer’ - based on gradient thresholding and refinement by local region of interest optimisation. See the find_peaks_zaefferer() function.

    • ’stat’ - based on statistical refinement and difference with respect to mean intensity. See the find_peaks_stat() function.

    • ’laplacian_of_gaussian’ - a blob finder using the laplacian of Gaussian matrices approach. See the find_peaks_log() function.

    • ’difference_of_gaussian’ - a blob finder using the difference of Gaussian matrices approach. See the find_peaks_log() function.

    • ’template_matching’ - A cross correlation peakfinder. This method requires providing a template with the template parameter, which is used as reference pattern to perform the template matching to the signal. It uses the skimage.feature.match_template() function and the peaks position are obtained by using minmax method on the template matching result.

  • interactive (bool) – If True, the method parameter can be adjusted interactively. If False, the results will be returned.

  • current_index (bool) – if True, the computation will be performed for the current index.

  • 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 multithreading. If None, the default from the preferences settings is used. The number of threads is controlled by the max_workers argument.

  • max_workers (None or int) – Maximum number of threads used when parallel=True. If None, defaults to min(32, os.cpu_count()).

  • 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.

  • **kwargs (dict) – Keywords parameters associated with above methods, see the documentation of each method for more details.

Notes

As a convenience, the ‘local_max’ method accepts the ‘distance’ and ‘threshold’ argument, which will be map to the ‘min_distance’ and ‘threshold_abs’ of the skimage.feature.peak_local_max() function.

Returns

peaks – Array of shape _navigation_shape_in_array in which each cell contains an array with dimensions (npeaks, 2) that contains the x, y pixel coordinates of peaks found in each image sorted first along y and then along x.

Return type

BaseSignal or numpy.ndarray if current_index=True

plot(navigator='auto', plot_markers=True, autoscale='v', saturated_pixels=None, norm='auto', vmin=None, vmax=None, gamma=1.0, linthresh=0.01, linscale=0.1, scalebar=True, scalebar_color='white', axes_ticks=None, axes_off=False, axes_manager=None, no_nans=False, colorbar=True, centre_colormap='auto', min_aspect=0.1, navigator_kwds={}, **kwargs)

Plot the signal at the current coordinates.

For multidimensional datasets an optional figure, the “navigator”, with a cursor to navigate that data is raised. In any case it is possible to navigate the data using the sliders. Currently only signals with signal_dimension equal to 0, 1 and 2 can be plotted.

Parameters
  • navigator (str, None, or BaseSignal (or subclass). Allowed string values are 'auto', 'slider', and 'spectrum'.) –

    If 'auto':

    • If navigation_dimension > 0, a navigator is provided to explore the data.

    • If navigation_dimension is 1 and the signal is an image the navigator is a sum spectrum obtained by integrating over the signal axes (the image).

    • If navigation_dimension is 1 and the signal is a spectrum the navigator is an image obtained by stacking all the spectra in the dataset horizontally.

    • If navigation_dimension is > 1, the navigator is a sum image obtained by integrating the data over the signal axes.

    • Additionally, if navigation_dimension > 2, a window with one slider per axis is raised to navigate the data.

    • For example, if the dataset consists of 3 navigation axes X, Y, Z and one signal axis, E, the default navigator will be an image obtained by integrating the data over E at the current Z index and a window with sliders for the X, Y, and Z axes will be raised. Notice that changing the Z-axis index changes the navigator in this case.

    • For lazy signals, the navigator will be calculated using the compute_navigator() method.

    If 'slider':

    • If navigation dimension > 0 a window with one slider per axis is raised to navigate the data.

    If 'spectrum':

    • If navigation_dimension > 0 the navigator is always a spectrum obtained by integrating the data over all other axes.

    • Not supported for lazy signals, the 'auto' option will be used instead.

    If None, no navigator will be provided.

    Alternatively a BaseSignal (or subclass) instance can be provided. The navigation or signal shape must match the navigation shape of the signal to plot or the navigation_shape + signal_shape must be equal to the navigator_shape of the current object (for a dynamic navigator). If the signal dtype is RGB or RGBA this parameter has no effect and the value is always set to 'slider'.

  • axes_manager (None or AxesManager) – If None, the signal’s axes_manager attribute is used.

  • plot_markers (bool, default True) – Plot markers added using s.add_marker(marker, permanent=True). Note, a large number of markers might lead to very slow plotting.

  • navigator_kwds (dict) – Only for image navigator, additional keyword arguments for matplotlib.pyplot.imshow().

  • colorbar (bool, optional) – If true, a colorbar is plotted for non-RGB images.

  • autoscale (str) – The string must contain any combination of the ‘x’, ‘y’ and ‘v’ characters. If ‘x’ or ‘y’ are in the string, the corresponding axis limits are set to cover the full range of the data at a given position. If ‘v’ (for values) is in the string, the contrast of the image will be set automatically according to vmin and vmax when the data or navigation indices change. Default is ‘v’.

  • saturated_pixels (scalar) –

    The percentage of pixels that are left out of the bounds. For example, the low and high bounds of a value of 1 are the 0.5% and 99.5% percentiles. It must be in the [0, 100] range. If None (default value), the value from the preferences is used.

    Deprecated since version 1.6.0: saturated_pixels will be removed in HyperSpy 2.0.0, it is replaced by vmin, vmax and autoscale.

  • norm ({“auto”, “linear”, “power”, “log”, “symlog” or a subclass of matplotlib.colors.Normalise}) – Set the norm of the image to display. If “auto”, a linear scale is used except if when power_spectrum=True in case of complex data type. “symlog” can be used to display negative value on a negative scale - read matplotlib.colors.SymLogNorm and the linthresh and linscale parameter for more details.

  • vmin ({scalar, str}, optional) – vmin and vmax are used to normalise the displayed data. It can be a float or a string. If string, it should be formatted as ‘xth’, where ‘x’ must be an float in the [0, 100] range. ‘x’ is used to compute the x-th percentile of the data. See numpy.percentile() for more information.

  • vmax ({scalar, str}, optional) – vmin and vmax are used to normalise the displayed data. It can be a float or a string. If string, it should be formatted as ‘xth’, where ‘x’ must be an float in the [0, 100] range. ‘x’ is used to compute the x-th percentile of the data. See numpy.percentile() for more information.

  • gamma (float) – Parameter used in the power-law normalisation when the parameter norm=”power”. Read matplotlib.colors.PowerNorm for more details. Default value is 1.0.

  • linthresh (float) – When used with norm=”symlog”, define the range within which the plot is linear (to avoid having the plot go to infinity around zero). Default value is 0.01.

  • linscale (float) – This allows the linear range (-linthresh to linthresh) to be stretched relative to the logarithmic range. Its value is the number of powers of base to use for each half of the linear range. See matplotlib.colors.SymLogNorm for more details. Defaulf value is 0.1.

  • scalebar (bool, optional) – If True and the units and scale of the x and y axes are the same a scale bar is plotted.

  • scalebar_color (str, optional) – A valid MPL color string; will be used as the scalebar color.

  • axes_ticks ({None, bool}, optional) – If True, plot the axes ticks. If None axes_ticks are only plotted when the scale bar is not plotted. If False the axes ticks are never plotted.

  • axes_off ({bool}) – Default is False.

  • no_nans (bool, optional) – If True, set nans to zero for plotting.

  • centre_colormap ({"auto", True, False}) – If True the centre of the color scheme is set to zero. This is specially useful when using diverging color schemes. If “auto” (default), diverging color schemes are automatically centred.

  • min_aspect (float) – Set the minimum aspect ratio of the image and the figure. To keep the image in the aspect limit the pixels are made rectangular.

  • **kwargs (dict) – Only when plotting an image: additional (optional) keyword arguments for matplotlib.pyplot.imshow().

hyperspy._signals.signal2d.estimate_image_shift(ref, image, roi=None, sobel=True, medfilter=True, hanning=True, plot=False, dtype='float', normalize_corr=False, sub_pixel_factor=1, return_maxval=True)

Estimate the shift in a image using phase correlation

This method can only estimate the shift by comparing bidimensional features that should not change the position in the given axis. To decrease the memory usage, the time of computation and the accuracy of the results it is convenient to select a region of interest by setting the roi keyword.

Parameters
  • ref (2D numpy.ndarray) – Reference image

  • image (2D numpy.ndarray) – Image to register

  • roi (tuple of ints (top, bottom, left, right)) – Define the region of interest

  • sobel (bool) – apply a sobel filter for edge enhancement

  • medfilter (bool) – apply a median filter for noise reduction

  • hanning (bool) – Apply a 2d hanning filter

  • plot (bool or matplotlib.Figure) – If True, plots the images after applying the filters and the phase correlation. If a figure instance, the images will be plotted to the given figure.

  • reference ('current' or 'cascade') – If ‘current’ (default) the image at the current coordinates is taken as reference. If ‘cascade’ each image is aligned with the previous one.

  • dtype (str or dtype) – Typecode or data-type in which the calculations must be performed.

  • normalize_corr (bool) – If True use phase correlation instead of standard correlation

  • sub_pixel_factor (float) – Estimate shifts with a sub-pixel accuracy of 1/sub_pixel_factor parts of a pixel. Default is 1, i.e. no sub-pixel accuracy.

Returns

  • shifts (np.array) – containing the estimate shifts

  • max_value (float) – The maximum value of the correlation

Notes

The statistical analysis approach to the translation estimation when using reference=’stat’ roughly follows * . If you use it please cite their article.

References

*

Bernhard Schaffer, Werner Grogger and Gerald Kothleitner. “Automated Spatial Drift Correction for EFTEM Image Series.” Ultramicroscopy 102, no. 1 (December 2004): 27–36.

hyperspy._signals.signal2d.fft_correlation(in1, in2, normalize=False, real_only=False)

Correlation of two N-dimensional arrays using FFT.

Adapted from scipy’s fftconvolve.

Parameters
  • in1 (array) – Input arrays to convolve.

  • in2 (array) – Input arrays to convolve.

  • normalize (bool, default False) – If True performs phase correlation.

  • real_only (bool, default False) – If True, and in1 and in2 are real-valued inputs, uses rfft instead of fft for approx. 2x speed-up.

hyperspy._signals.signal2d.hanning2d(M, N)

A 2D hanning window created by outer product.

hyperspy._signals.signal2d.triu_indices_minus_diag(n)

Returns the indices for the upper-triangle of an (n, n) array excluding its diagonal

Parameters

n (int) – The length of the square array