API documentation

Classes

PixelatedSTEM

class pixstem.pixelated_stem_class.PixelatedSTEM(*args, **kw)[source]
add_ellipse_array_as_markers(ellipse_array, inlier_array=None, peak_array=None, nr=20, color_ellipse='blue', linewidth=1, linestyle='solid', color_inlier='blue', color_outlier='red', point_size=20)[source]

Add a ellipse parameters array to a signal as HyperSpy markers.

Useful to visualize the ellipse results.

Parameters:
  • ellipse_array (NumPy array) –
  • inlier_array (NumPy array, optional) –
  • peak_array (NumPy array, optional) –
  • nr (scalar, optional) – Default 20
  • color_ellipse (string, optional) – Default ‘blue’
  • linewidth (scalar, optional) – Default 1
  • linestyle (string, optional) – Default ‘solid’
  • color_inlier (string, optional) – Default ‘blue’
  • color_outlier (string, optional) – Default ‘red’
  • point_size (scalar, optional) –

Examples

>>> s, parray = ps.dummy_data.get_simple_ellipse_signal_peak_array()
>>> import pixstem.ransac_ellipse_tools as ret
>>> ellipse_array, inlier_array = ret.get_ellipse_model_ransac(
...     parray, xf=95, yf=95, rf_lim=20, semi_len_min=40,
...     semi_len_max=100, semi_len_ratio_lim=5, max_trails=50)
>>> s.add_ellipse_array_as_markers(
...     ellipse_array, inlier_array=inlier_array, peak_array=parray)
>>> s.plot()
add_peak_array_as_markers(peak_array, color='red', size=20, bool_array=None, bool_invert=False)[source]

Add a peak array to the signal as HyperSpy markers.

Parameters:
  • peak_array (NumPy 4D array) –
  • color (string, optional) – Default ‘red’
  • size (scalar, optional) – Default 20
  • bool_array (NumPy array) – Must be the same size as peak_array
  • bool_invert (bool) –

Examples

>>> s, parray = ps.dummy_data.get_simple_ellipse_signal_peak_array()
>>> s.add_peak_array_as_markers(parray)
>>> s.plot()
angular_mask(angle0, angle1, centre_x_array=None, centre_y_array=None)[source]

Get a bool array with True values between angle0 and angle1. Will use the (0, 0) point as given by the signal as the centre, giving an “angular” slice. Useful for analysing anisotropy in diffraction patterns.

Parameters:
  • angle1 (angle0,) –
  • centre_y_array (centre_x_array,) – Has to have the same shape as the navigation axis of the signal.
Returns:

mask_array – The True values will be the region between angle0 and angle1. The array will have the same dimensions as the signal.

Return type:

NumPy array

Examples

>>> import pixstem.dummy_data as dd
>>> s = dd.get_holz_simple_test_signal()
>>> s.axes_manager.signal_axes[0].offset = -25
>>> s.axes_manager.signal_axes[1].offset = -25
>>> mask_array = s.angular_mask(0.5*np.pi, np.pi)
angular_slice_radial_average(angleN=20, centre_x=None, centre_y=None, slice_overlap=None, show_progressbar=True)[source]

Do radial average of different angular slices. Useful for analysing anisotropy in round diffraction features, such as diffraction rings from polycrystalline materials or higher order Laue zone rings.

Parameters:
  • angleN (int, default 20) – Number of angular slices. If angleN=4, each slice will be 90 degrees. The average will start in the top left corner (0, 0) when plotting using s.plot(), and go clockwise.
  • centre_y (centre_x,) – If given as int, all the diffraction patterns will have the same centre position. Each diffraction pattern can also have different centre position, by passing a NumPy array with the same dimensions as the navigation axes. Note: in either case both x and y values must be given. If one is missing, both will be set from the signal (0., 0.) positions. If no values are given, the (0., 0.) positions in the signal will be used.
  • slice_overlap (float, optional) – Amount of overlap between the slices, given in fractions of angle slice (0 to 1). For angleN=4, each slice will be 90 degrees. If slice_overlap=0.5, each slice will overlap by 45 degrees on each side. The range of the slices will then be: (-45, 135), (45, 225), (135, 315) and (225, 45). Default off: meaning there is no overlap between the slices.
  • show_progressbar (bool) – Default True
Returns:

signal – With one more navigation dimensions (the angular slices) compared to the input signal.

Return type:

HyperSpy 1D signal

Examples

>>> import pixstem.dummy_data as dd
>>> s = dd.get_holz_simple_test_signal()
>>> s_com = s.center_of_mass(show_progressbar=False)
>>> s_ar = s.angular_slice_radial_average(
...     angleN=10, centre_x=s_com.inav[0].data,
...     centre_y=s_com.inav[1].data, slice_overlap=0.2,
...     show_progressbar=False)
>>> s_ar.plot() # doctest: +SKIP
angular_slice_radial_integration()[source]
as_lazy(*args, **kwargs)[source]

Create a copy of the given Signal as a LazySignal.

Parameters:copy_variance (bool) – Whether or not to copy the variance from the original Signal to the new lazy version
Returns:res – The same signal, converted to be lazy
Return type:LazySignal
center_of_mass(threshold=None, mask=None, lazy_result=False, show_progressbar=True, chunk_calculations=None)[source]

Get the centre of the STEM diffraction pattern using center of mass. Threshold can be set to only use the most intense parts of the pattern. A mask can be used to exclude parts of the diffraction pattern.

Parameters:
  • threshold (number, optional) – The thresholding will be done at mean times this threshold value.
  • mask (tuple (x, y, r), optional) – Round mask centered on x and y, with radius r.
  • lazy_result (bool, optional) – If True, will not compute the data directly, but return a lazy signal. Default False
  • show_progressbar (bool, optional) – Default True
  • chunk_calculations (tuple, optional) – Chunking values when running the calculations.
Returns:

s_com – DPCSignal with beam shifts along the navigation dimension and spatial dimensions as the signal dimension(s).

Return type:

DPCSignal

Examples

With mask centered at x=105, y=120 and 30 pixel radius

>>> import pixstem.dummy_data as dd
>>> s = dd.get_disk_shift_simple_test_signal()
>>> mask = (25, 25, 10)
>>> s_com = s.center_of_mass(mask=mask, show_progressbar=False)
>>> s_color = s_com.get_color_signal()

Also threshold

>>> s_com = s.center_of_mass(threshold=1.5, show_progressbar=False)

Get a lazy signal, then calculate afterwards

>>> s_com = s.center_of_mass(lazy_result=True, show_progressbar=False)
>>> s_com.compute(progressbar=False)
correct_bad_pixels(bad_pixel_array, lazy_result=True, show_progressbar=True)[source]

Correct bad pixels by getting mean value of neighbors.

Note: this method is currently not very optimized with regards to memory use, so currently be careful when using it on large datasets.

Parameters:
  • bad_pixel_array (array-like) –
  • lazy_result (bool) – Default True.
  • show_progressbar (bool) – Default True
Returns:

signal_corrected

Return type:

PixelatedSTEM or LazyPixelatedSTEM

Examples

>>> s = ps.dummy_data.get_hot_pixel_signal()
>>> s_hot_pixels = s.find_hot_pixels(
...     show_progressbar=False, lazy_result=True)
>>> s_corr = s.correct_bad_pixels(s_hot_pixels)

Dead pixels

>>> s = ps.dummy_data.get_dead_pixel_signal()
>>> s_dead_pixels = s.find_dead_pixels(
...     show_progressbar=False, lazy_result=True)
>>> s_corr = s.correct_bad_pixels(s_dead_pixels)

Combine both dead pixels and hot pixels

>>> s_bad_pixels = s_hot_pixels + s_dead_pixels
>>> s_corr = s.correct_bad_pixels(s_bad_pixels)
fem_analysis(centre_x=None, centre_y=None, show_progressbar=True)[source]

Perform analysis of fluctuation electron microscopy (FEM) data.

This is outlined in: T. L. Daulton, et al., Ultramicroscopy 110 (2010) 1279-1289. doi:10.1016/j.ultramic.2010.05.010

Parameters:
  • centre_y (centre_x,) – All the diffraction patterns assumed to have the same centre position.
  • show_progressbar (bool, optional) – Default True
Returns:

results – Results of FEM data analysis, including the normalized variance of the annular mean (V-Omegak), mean of normalized variances of rings (V-rk), normalized variance of ring ensemble (Vrek), the normalized variance image (Omega-Vi), and annular mean of the variance image (Omega-Vk).

Return type:

Python dictionary

Examples

>>> s = ps.dummy_data.get_fem_signal()
>>> fem_results = s.fem_analysis(
...     centre_x=50, centre_y=50,
...     show_progressbar=False)
>>> fem_results['V-Omegak'].plot()
find_dead_pixels(dead_pixel_value=0, mask_array=None, lazy_result=False, show_progressbar=True)[source]

Find dead pixels in the diffraction images.

Parameters:
  • dead_pixel_value (scalar) – Default 0
  • mask_array (Boolean Numpy array) –
  • lazy_result (bool) – If True, return a lazy signal. If False, compute the result and return a non-lazy signal. Default False.
  • show_progressbar (bool) –
Returns:

s_dead_pixels – With dead pixels as True, rest as False.

Return type:

HyperSpy 2D signal

Examples

>>> s = ps.dummy_data.get_dead_pixel_signal()
>>> s_dead_pixels = s.find_dead_pixels(show_progressbar=False)

Using a mask array

>>> import numpy as np
>>> mask_array = np.zeros((128, 128), dtype=np.bool)
>>> mask_array[:, 100:] = True
>>> s = ps.dummy_data.get_dead_pixel_signal()
>>> s_dead_pixels = s.find_dead_pixels(
...     mask_array=mask_array, show_progressbar=False)

Getting a lazy signal as output

>>> s_dead_pixels = s.find_dead_pixels(
...     lazy_result=True, show_progressbar=False)
find_hot_pixels(threshold_multiplier=500, mask_array=None, lazy_result=True, show_progressbar=True)[source]

Find hot pixels in the diffraction images.

Note: this method will be default return a lazy signal, since the size of the returned signal is the same shape as the original signal. So for large datasets actually calculating computing the results can use a lot of memory.

In addition, this signal is currently not very optimized with regards to memory use, so be careful when using this method for large datasets.

Parameters:
  • threshold_multiplier (scalar) – Default 500
  • mask_array (Boolean NumPy array) –
  • lazy_result (bool) – If True, return a lazy signal. If False, compute the result and return a non-lazy signal. Default True.
  • show_progressbar (bool) –

Examples

>>> s = ps.dummy_data.get_hot_pixel_signal()
>>> s_hot_pixels = s.find_hot_pixels(show_progressbar=False)

Using a mask array

>>> import numpy as np
>>> mask_array = np.zeros((128, 128), dtype=np.bool)
>>> mask_array[:, 100:] = True
>>> s = ps.dummy_data.get_hot_pixel_signal()
>>> s_hot_pixels = s.find_hot_pixels(
...     mask_array=mask_array, show_progressbar=False)

Getting a non-lazy signal as output

>>> s_hot_pixels = s.find_hot_pixels(
...     lazy_result=False, show_progressbar=False)
find_peaks(method='dog', lazy_result=True, show_progressbar=True, **kwargs)[source]

Find peaks in the signal dimensions.

Can use either skimage’s blob_dog or blob_log.

Parameters:
  • method (string, optional) – ‘dog’: difference of Gaussians. ‘log’: Laplacian of Gaussian. Default ‘dog’.
  • min_sigma (float, optional) – Default 0.98.
  • max_sigma (float, optional) – Default 55.
  • sigma_ratio (float, optional) – For method ‘dog’. Default 1.76.
  • num_sigma (float, optional) – For method ‘log’. Default 10.
  • threshold (float, optional) – Default 0.36.
  • overlap (float, optional) – Default 0.81.
  • normalize_value (float, optional) – All the values in the signal will be divided by this value. If no value is specified, the max value in each individual image will be used.
  • max_r (float) – Maximum radius compared from the center of the diffraction pattern
  • lazy_result (bool, optional) – Default True
  • show_progressbar (bool, optional) – Default True
Returns:

peak_array – Same size as the two last dimensions in data, in the form [[y0, x0], [y1, x1], …]. The peak positions themselves are stored in 2D NumPy arrays inside each position in peak_array. This is done instead of making a 4D NumPy array, since the number of found peaks can vary in each position.

Return type:

dask 2D object array

Example

>>> s = ps.dummy_data.get_cbed_signal()
>>> peak_array = s.find_peaks()
>>> peak_array_computed = peak_array.compute(show_progressbar=False)
>>> peak02 = peak_array_computed[0, 2]
>>> s.add_peak_array_as_markers(peak_array_computed)
>>> s.plot()

Change parameters

>>> peak_array = s.find_peaks(
...     method='dog', min_sigma=1.2, max_sigma=27, sigma_ratio=2.2,
...     threshold=0.6, overlap=0.6, lazy_result=False,
...     show_progressbar=False)

Using Laplacian of Gaussian

>>> s = ps.dummy_data.get_cbed_signal()
>>> peak_array = s.find_peaks(
...     method='log', min_sigma=5, max_sigma=55, num_sigma=10,
...     threshold=0.2, overlap=0.86, lazy_result=False,
...     show_progressbar=False)
>>> s.add_peak_array_as_markers(peak_array)
>>> s.plot()
flip_diffraction_x()[source]

Flip the dataset along the diffraction x-axis.

The function returns a new signal, but the data itself is a view of the original signal. So changing the returned signal will also change the original signal (and visa versa). To avoid changing the original signal, use the deepcopy method afterwards, but note that this requires double the amount of memory. See below for an example of this.

Returns:flipped_signal
Return type:PixelatedSTEM signal

Example

>>> s = ps.dummy_data.get_holz_simple_test_signal()
>>> s_flip = s.flip_diffraction_x()

To avoid changing the original object afterwards

>>> s_flip = s.flip_diffraction_x().deepcopy()
flip_diffraction_y()[source]

Flip the dataset along the diffraction y-axis.

The function returns a new signal, but the data itself is a view of the original signal. So changing the returned signal will also change the original signal (and visa versa). To avoid changing the original signal, use the deepcopy method afterwards, but note that this requires double the amount of memory. See below for an example of this.

Returns:flipped_signal
Return type:PixelatedSTEM signal

Example

>>> s = ps.dummy_data.get_holz_simple_test_signal()
>>> s_flip = s.flip_diffraction_y()

To avoid changing the original object afterwards

>>> s_flip = s.flip_diffraction_y().deepcopy()
intensity_peaks(peak_array, disk_r=4, lazy_result=True, show_progressbar=True)[source]

Get intensity of a peak in the diffraction data.

The intensity is calculated by taking the mean of the pixel values inside radius disk_r from the peak position.

Parameters:
  • peak_array (Numpy or Dask array) – Must have the same navigation shape as this signal.
  • disk_r (int) – Radius of the disc chosen to take the mean value of
  • lazy_result (bool, default True) – If True, will return a LazyPixelatedSTEM object. If False, will compute the result and return a PixelatedSTEM object.
  • show_progressbar (bool, default True) –
Returns:

intensity_array – Same navigation shape as this signal, with peak position in x and y coordinates and the mean intensity.

Return type:

Numpy or Dask array

Examples

>>> s = ps.dummy_data.get_cbed_signal()
>>> peak_array = s.find_peaks()
>>> intensity_array = s.intensity_peaks(peak_array, disk_r=6)
>>> intensity_array_computed = intensity_array.compute()
peak_position_refinement_com(peak_array, square_size=10, lazy_result=True, show_progressbar=True)[source]

Refines the peak position using the center of mass.

Parameters:
  • peak_array (Numpy or Dask array) – Object with x and y coordinates of the peak positions. Must have the same dimensions as this signal’s navigation dimensions.
  • square_size (int) – Even integer, sub image from which the center of mass is calculated. Default 5.
  • lazy_result (bool, default True) – If True, will return a LazyPixelatedSTEM object. If False, will compute the result and return a PixelatedSTEM object.
  • show_progressbar (bool, default True) –
Returns:

output_array – Same size as the two last dimensions in data, in the form [[y0, x0], [y1, x1], …]. The peak positions themselves are stored in 2D NumPy arrays inside each position in peak_array. This is done instead of making a 4D NumPy array, since the number of found peaks can vary in each position.

Return type:

dask 2D object array

Examples

>>> s = ps.dummy_data.get_cbed_signal()
>>> peak_array = s.find_peaks()
>>> refined_peak_array = s.peak_position_refinement_com(peak_array, 20)
>>> refined_peak_array_com = refined_peak_array.compute(
...     show_progressbar=False)
>>> s.add_peak_array_as_markers(refined_peak_array_com)
>>> s.plot()
plot(*args, **kwargs)[source]

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)) –
  • string values are 'auto', 'slider', and 'spectrum'. (Allowed) –
  • 'auto' (If) –
    • 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.

    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.

    If None, no navigator will be provided.

    Alternatively a BaseSignal (or subclass) instance can be provided. The signal_dimension must be 1 (for a spectrum navigator) or 2 (for a image navigator) and navigation_shape must be 0 (for a static navigator) or 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.
:param .. deprecated:: 1.6.0:
saturated_pixels will be removed in HyperSpy 2.0.0, it is replaced
by vmin, vmax and autoscale.
Parameters:
  • 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.
  • vmax (vmin,) – 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 – Only when plotting an image: additional (optional) keyword arguments for matplotlib.pyplot.imshow().
radial_average(centre_x=None, centre_y=None, mask_array=None, normalize=True, parallel=True, show_progressbar=True)[source]

Radially average a pixelated STEM diffraction signal.

Done by integrating over the azimuthal dimension, giving a profile of intensity as a function of scattering angle.

Parameters:
  • centre_y (centre_x,) – If given as int, all the diffraction patterns will have the same centre position. Each diffraction pattern can also have different centre position, by passing a NumPy array with the same dimensions as the navigation axes. Note: in either case both x and y values must be given. If one is missing, both will be set from the signal (0., 0.) positions. If no values are given, the (0., 0.) positions in the signal will be used.
  • mask_array (Boolean NumPy array, optional) – Mask with the same shape as the signal.
  • normalize (bool, default True) – If true, the returned radial profile will be normalized by the number of bins used for each average.
  • parallel (bool, default True) – If True, run the processing on several cores. In most cases this should be True, but for debugging False can be useful.
  • show_progressbar (bool) – Default True
Returns:

Return type:

HyperSpy signal, one less signal dimension than the input signal.

Examples

>>> import pixstem.dummy_data as dd
>>> s = dd.get_holz_simple_test_signal()
>>> s_r = s.radial_average(centre_x=25, centre_y=25,
...     show_progressbar=False)
>>> s_r.plot()

Using center_of_mass to find bright field disk position

>>> s = dd.get_disk_shift_simple_test_signal()
>>> s_com = s.center_of_mass(threshold=2, show_progressbar=False)
>>> s_r = s.radial_average(
...     centre_x=s_com.inav[0].data, centre_y=s_com.inav[1].data,
...     show_progressbar=False)
>>> s_r.plot()
radial_integration()[source]
rotate_diffraction(angle, parallel=True, show_progressbar=True)[source]

Rotate the diffraction dimensions.

Parameters:
  • angle (scalar) – Clockwise rotation in degrees.
  • parallel (bool) – Default True
  • show_progressbar (bool) – Default True
Returns:

rotated_signal

Return type:

PixelatedSTEM class

Examples

>>> s = ps.dummy_data.get_holz_simple_test_signal()
>>> s_rot = s.rotate_diffraction(30, show_progressbar=False)
shift_diffraction(shift_x, shift_y, interpolation_order=1, parallel=True, inplace=False, show_progressbar=True)[source]

Shift the diffraction patterns in a pixelated STEM signal.

The points outside the boundaries are set to zero.

Parameters:
  • shift_y (shift_x,) – If given as int, all the diffraction patterns will have the same shifts. Each diffraction pattern can also have different shifts, by passing a NumPy array with the same dimensions as the navigation axes.
  • interpolation_order (int) – When shifting, a spline interpolation is used. This parameter sets the order of this spline. Must be between 0 and 5. Note that in some low-signal and high noise datasets, using a non-zero order might lead to artifacts. See the docstring in scipy.ndimage.shift for more information. Default 1.
  • parallel (bool) – If True, run the processing on several cores. In most cases this should be True, but for debugging False can be useful. Default True
  • inplace (bool) – If True (default), the data is replaced by the result. Useful when working with very large datasets, as this avoids doubling the amount of memory needed. If False, a new signal with the results is returned.
  • show_progressbar (bool) – Default True.
Returns:

shifted_signal

Return type:

PixelatedSTEM signal

Examples

>>> s = ps.dummy_data.get_disk_shift_simple_test_signal()
>>> s_c = s.center_of_mass(threshold=3., show_progressbar=False)
>>> s_c -= 25 # To shift the center disk to the middle (25, 25)
>>> s_shift = s.shift_diffraction(
...     s_c.inav[0].data, s_c.inav[1].data,
...     show_progressbar=False)
>>> s_shift.plot()

Using a different interpolation order

>>> s_shift = s.shift_diffraction(
...     s_c.inav[0].data, s_c.inav[1].data, interpolation_order=3,
...     show_progressbar=False)
subtract_diffraction_background(method='median kernel', lazy_result=True, show_progressbar=True, **kwargs)[source]

Background subtraction of the diffraction data.

There are three different methods for doing this: - Difference of Gaussians - Median kernel - Radial median

Parameters:
  • method (string) – ‘difference of gaussians’, ‘median kernel’ and ‘radial median’. Default ‘median kernel’.
  • lazy_result (bool, default True) – If True, will return a LazyPixelatedSTEM object. If False, will compute the result and return a PixelatedSTEM object.
  • show_progressbar (bool, default True) –
  • sigma_min (float, optional) – Standard deviation for the minimum Gaussian convolution (difference of Gaussians only)
  • sigma_max (float, optional) – Standard deviation for the maximum Gaussian convolution (difference of Gaussians only)
  • footprint (int, optional) – Size of the window that is convoluted with the array to determine the median. Should be large enough that it is about 3x as big as the size of the peaks (median kernel only).
  • centre_x (int, optional) – Centre x position of the coordinate system on which to map to radial coordinates (radial median only).
  • centre_y (int, optional) – Centre y position of the coordinate system on which to map to radial coordinates (radial median only).
Returns:

s

Return type:

PixelatedSTEM or LazyPixelatedSTEM signal

Examples

>>> s = ps.dummy_data.get_cbed_signal()
>>> s_r = s.subtract_diffraction_background(method='median kernel',
...     footprint=20, lazy_result=False, show_progressbar=False)
>>> s_r.plot()
template_match_disk(disk_r=4, lazy_result=True, show_progressbar=True)[source]

Template match the signal dimensions with a disk.

Used to find diffraction disks in convergent beam electron diffraction data.

Parameters:
  • disk_r (scalar, optional) – Radius of the disk. Default 4.
  • lazy_result (bool, default True) – If True, will return a LazyPixelatedSTEM object. If False, will compute the result and return a PixelatedSTEM object.
  • show_progressbar (bool, default True) –
Returns:

template_match

Return type:

PixelatedSTEM object

Examples

>>> s = ps.dummy_data.get_cbed_signal()
>>> s_template = s.template_match_disk(
...     disk_r=5, show_progressbar=False)
>>> s.plot()
template_match_ring(r_inner=5, r_outer=7, lazy_result=True, show_progressbar=True)[source]

Template match the signal dimensions with a ring.

Used to find diffraction rings in convergent beam electron diffraction data.

Parameters:
  • r_outer (r_inner,) – Inner and outer radius of the rings.
  • lazy_result (bool, default True) – If True, will return a LazyPixelatedSTEM object. If False, will compute the result and return a PixelatedSTEM object.
  • show_progressbar (bool, default True) –
Returns:

template_match

Return type:

PixelatedSTEM object

Examples

>>> s = ps.dummy_data.get_cbed_signal()
>>> s_template = s.template_match_ring(show_progressbar=False)
>>> s.plot()
template_match_with_binary_image(binary_image, lazy_result=True, show_progressbar=True)[source]

Template match the signal dimensions with a binary image.

Used to find diffraction disks in convergent beam electron diffraction data.

Might also work with non-binary images, but this haven’t been extensively tested.

Parameters:
  • binary_image (2-D NumPy array) –
  • lazy_result (bool, default True) – If True, will return a LazyPixelatedSTEM object. If False, will compute the result and return a PixelatedSTEM object.
  • show_progressbar (bool, default True) –
Returns:

template_match

Return type:

PixelatedSTEM object

Examples

>>> s = ps.dummy_data.get_cbed_signal()
>>> binary_image = np.random.randint(0, 2, (6, 6))
>>> s_template = s.template_match_with_binary_image(
...     binary_image, show_progressbar=False)
>>> s.plot()
threshold_and_mask(threshold=None, mask=None, show_progressbar=True)[source]

Get a thresholded and masked of the signal.

Useful for figuring out optimal settings for the center_of_mass method.

Parameters:
  • threshold (number, optional) – The thresholding will be done at mean times this threshold value.
  • mask (tuple (x, y, r)) – Round mask centered on x and y, with radius r.
  • show_progressbar (bool) – Default True
Returns:

s_out

Return type:

PixelatedSTEM signal

Examples

>>> import pixstem.dummy_data as dd
>>> s = dd.get_disk_shift_simple_test_signal()
>>> mask = (25, 25, 10)
>>> s_out = s.threshold_and_mask(
...     mask=mask, threshold=2, show_progressbar=False)
>>> s_out.plot()
virtual_annular_dark_field(cx, cy, r_inner, r, lazy_result=False, show_progressbar=True)[source]

Get a virtual annular dark field signal.

Parameters:
  • cy (cx,) – x- and y-centre positions.
  • r_inner (float) – Inner radius.
  • r (float) – Outer radius.
  • lazy_result (bool, optional) – If True, will not compute the data directly, but return a lazy signal. Default False
  • show_progressbar (bool, default True) –
Returns:

virtual_adf_signal

Return type:

HyperSpy 2D signal

Examples

>>> s = ps.dummy_data.get_holz_heterostructure_test_signal()
>>> s_adf = s.virtual_annular_dark_field(
...     40, 40, 20, 40, show_progressbar=False)
>>> s_adf.plot()

Get a lazy signal, then compute

>>> s_adf = s.virtual_annular_dark_field(
...     40, 40, 20, 40, lazy_result=True, show_progressbar=False)
>>> s_adf.compute(progressbar=False)
>>> s_adf.plot()
virtual_bright_field(cx=None, cy=None, r=None, lazy_result=False, show_progressbar=True)[source]

Get a virtual bright field signal.

Can be sum the whole diffraction plane, or a circle subset. If any of the parameters are None, it will sum the whole diffraction plane.

Parameters:
  • cy (cx,) – x- and y-centre positions.
  • r (float, optional) – Outer radius.
  • lazy_result (bool, optional) – If True, will not compute the data directly, but return a lazy signal. Default False
  • show_progressbar (bool, optional) – Default True.
Returns:

virtual_bf_signal

Return type:

HyperSpy 2D signal

Examples

>>> s = ps.dummy_data.get_holz_heterostructure_test_signal()
>>> s_bf = s.virtual_bright_field(show_progressbar=False)
>>> s_bf.plot()

Sum a subset of the diffraction pattern

>>> s_bf = s.virtual_bright_field(40, 40, 10, show_progressbar=False)
>>> s_bf.plot()

Get a lazy signal, then compute

>>> s_bf = s.virtual_bright_field(
...     lazy_result=True, show_progressbar=False)
>>> s_bf.compute(progressbar=False)

DPCBaseSignal

class pixstem.pixelated_stem_class.DPCBaseSignal(*args, **kwargs)[source]

Signal for processing differential phase contrast (DPC) acquired using scanning transmission electron microscopy (STEM).

The signal assumes the data is 3 dimensions, where the two signal dimensions are the probe positions, and the navigation dimension is the x and y disk shifts.

The first navigation index (s.inav[0]) is assumed to the be x-shift and the second navigation is the y-shift (s.inav[1]).

DPCSignal1D

class pixstem.pixelated_stem_class.DPCSignal1D(*args, **kwargs)[source]

Signal for processing differential phase contrast (DPC) acquired using scanning transmission electron microscopy (STEM).

The signal assumes the data is 2 dimensions, where the signal dimension is the probe position, and the navigation dimension is the x and y disk shifts.

The first navigation index (s.inav[0]) is assumed to the be x-shift and the second navigation is the y-shift (s.inav[1]).

get_bivariate_histogram(histogram_range=None, masked=None, bins=200, spatial_std=3)[source]

Useful for finding the distribution of magnetic vectors(?).

Parameters:
  • histogram_range (tuple, optional) – Set the minimum and maximum of the histogram range. Default is setting it automatically.
  • masked (1-D NumPy bool array, optional) – Mask parts of the data. The array must be the same size as the signal. The True values are masked. Default is not masking anything.
  • bins (integer, default 200) – Number of bins in the histogram
  • spatial_std (number, optional) – If histogram_range is not given, this value will be used to set the automatic histogram range. Default value is 3.
Returns:

s_hist

Return type:

Signal2D

DPCSignal2D

class pixstem.pixelated_stem_class.DPCSignal2D(*args, **kw)[source]

Signal for processing differential phase contrast (DPC) acquired using scanning transmission electron microscopy (STEM).

The signal assumes the data is 3 dimensions, where the two signal dimensions are the probe positions, and the navigation dimension is the x and y disk shifts.

The first navigation index (s.inav[0]) is assumed to the be x-shift and the second navigation is the y-shift (s.inav[1]).

correct_ramp(corner_size=0.05, only_offset=False, out=None)[source]

Subtracts a plane from the signal, useful for removing the effects of d-scan in a STEM beam shift dataset.

The plane is calculated by fitting a plane to the corner values of the signal. This will only work well when the property one wants to measure is zero in these corners.

Parameters:
  • corner_size (number, optional) – The size of the corners, as a percentage of the image’s axis. If corner_size is 0.05 (5%), and the image is 500 x 1000, the size of the corners will be (500*0.05) x (1000*0.05) = 25 x 50. Default 0.05
  • only_offset (bool, optional) – If True, will subtract a “flat” plane, i.e. it will subtract the mean value of the corners. Default False
  • out (optional, DPCSignal2D signal) –
Returns:

corrected_signal

Return type:

Signal2D

Examples

>>> s = ps.dummy_data.get_square_dpc_signal(add_ramp=True)
>>> s_corr = s.correct_ramp()
>>> s_corr.plot()

Only correct offset

>>> s_corr = s.correct_ramp(only_offset=True)
>>> s_corr.plot()
flip_axis_90_degrees(flips=1)[source]

Flip both the spatial and beam deflection axis

Will rotate both the image and the beam deflections by 90 degrees.

Parameters:flips (int, default 1) – Number of flips. The default (1) gives 90 degrees rotation. 2 gives 180, 3 gives 270, …

Examples

>>> s = ps.dummy_data.get_stripe_pattern_dpc_signal()
>>> s
<DPCSignal2D, title: , dimensions: (2|50, 100)>
>>> s_rot = s.flip_axis_90_degrees()
>>> s_rot
<DPCSignal2D, title: , dimensions: (2|100, 50)>

Do several flips

>>> s_rot = s.flip_axis_90_degrees(2)
>>> s_rot
<DPCSignal2D, title: , dimensions: (2|50, 100)>
>>> s_rot = s.flip_axis_90_degrees(3)
>>> s_rot
<DPCSignal2D, title: , dimensions: (2|100, 50)>
gaussian_blur(sigma=2, output=None)[source]

Blur the x- and y-beam shifts.

Useful for reducing the effects of structural diffraction effects.

Parameters:
  • sigma (scalar, default 2) –
  • output (HyperSpy signal) –
Returns:

blurred_signal

Return type:

HyperSpy 2D Signal

Examples

>>> s = ps.dummy_data.get_square_dpc_signal(add_ramp=False)
>>> s_blur = s.gaussian_blur()

Different sigma

>>> s_blur = s.gaussian_blur(sigma=1.2)

Using the signal itself as output

>>> s.gaussian_blur(output=s)
>>> s.plot()
get_bivariate_histogram(histogram_range=None, masked=None, bins=200, spatial_std=3)[source]

Useful for finding the distribution of magnetic vectors(?).

Parameters:
  • histogram_range (tuple, optional) – Set the minimum and maximum of the histogram range. Default is setting it automatically.
  • masked (2-D NumPy bool array, optional) – Mask parts of the data. The array must be the same size as the signal. The True values are masked. Default is not masking anything.
  • bins (integer, default 200) – Number of bins in the histogram
  • spatial_std (number, optional) – If histogram_range is not given, this value will be used to set the automatic histogram range. Default value is 3.
Returns:

s_hist

Return type:

HyperSpy Signal2D

Examples

>>> s = ps.dummy_data.get_stripe_pattern_dpc_signal()
>>> s_hist = s.get_bivariate_histogram()
>>> s_hist.plot()
get_color_image_with_indicator(phase_rotation=0, indicator_rotation=0, only_phase=False, autolim=True, autolim_sigma=4, scalebar_size=None, ax=None, ax_indicator=None)[source]

Make a matplotlib figure showing DPC contrast.

Parameters:
  • phase_rotation (float, default 0) – Changes the phase of the plotted data. Useful for correcting scan rotation.
  • indicator_rotation (float, default 0) – Changes the color wheel rotation.
  • only_phase (bool, default False) – If False, will plot both the magnitude and phase. If True, will only plot the phase.
  • autolim (bool, default True) –
  • autolim_sigma (float, default 4) –
  • scalebar_size (int, optional) –
  • ax (Matplotlib subplot, optional) –
  • ax_indicator (Matplotlib subplot, optional) – If None, generate a new subplot for the indicator. If False, do not include an indicator

Examples

>>> s = ps.dummy_data.get_simple_dpc_signal()
>>> fig = s.get_color_image_with_indicator()
>>> fig.savefig("simple_dpc_test_signal.png")

Only plotting the phase

>>> fig = s.get_color_image_with_indicator(only_phase=True)
>>> fig.savefig("simple_dpc_test_signal.png")

Matplotlib subplot as input

>>> import matplotlib.pyplot as plt
>>> fig, ax = plt.subplots()
>>> ax_indicator = fig.add_subplot(331)
>>> fig_return = s.get_color_image_with_indicator(
...     scalebar_size=10, ax=ax, ax_indicator=ax_indicator)
get_color_signal(rotation=None, autolim=True, autolim_sigma=4)[source]

Get DPC image visualized using continuous color scale.

Converts the x and y beam shifts into an RGB array, showing the magnitude and direction of the beam shifts.

Useful for visualizing magnetic domain structures.

Parameters:
  • rotation (float, optional) – In degrees. Useful for correcting the mismatch between scan direction and diffraction pattern rotation.
  • autolim (bool, default True) –
  • autolim_sigma (float, default 4) –
Returns:

color_signal

Return type:

HyperSpy 2D RGB signal

Examples

>>> s = ps.dummy_data.get_simple_dpc_signal()
>>> s_color = s.get_color_signal()
>>> s_color.plot()

Rotate the beam shift by 30 degrees

>>> s_color = s.get_color_signal(rotation=30)

See also

get_color_signal()
Signal showing both phase and magnitude
get_phase_signal()
Signal showing the phase
get_magnitude_signal(autolim=True, autolim_sigma=4)[source]

Get DPC magnitude image visualized as greyscale.

Converts the x and y beam shifts into a magnitude map, showing the magnitude of the beam shifts.

Useful for visualizing magnetic domain structures.

Parameters:
  • autolim (bool, default True) –
  • autolim_sigma (float, default 4) –
Returns:

magnitude_signal

Return type:

HyperSpy 2D signal

Examples

>>> s = ps.dummy_data.get_simple_dpc_signal()
>>> s_magnitude = s.get_magnitude_signal()
>>> s_magnitude.plot()

See also

get_color_signal()
Signal showing both phase and magnitude
get_phase_signal()
Signal showing the phase
get_phase_signal(rotation=None)[source]

Get DPC phase image visualized using continuous color scale.

Converts the x and y beam shifts into an RGB array, showing the direction of the beam shifts.

Useful for visualizing magnetic domain structures.

Parameters:
  • rotation (float, optional) – In degrees. Useful for correcting the mismatch between scan direction and diffraction pattern rotation.
  • autolim (bool, default True) –
  • autolim_sigma (float, default 4) –
Returns:

phase_signal

Return type:

HyperSpy 2D RGB signal

Examples

>>> s = ps.dummy_data.get_simple_dpc_signal()
>>> s_color = s.get_phase_signal(rotation=20)
>>> s_color.plot()

See also

get_color_signal()
Signal showing both phase and magnitude
get_magnitude_signal()
Signal showing the magnitude
rotate_beam_shifts(angle)[source]

Rotate the beam shift vector.

Parameters:angle (float) – Clockwise rotation in degrees
Returns:shift_rotated_signal
Return type:DPCSignal2D

Example

Rotate beam shifts by 10 degrees clockwise

>>> s = ps.dummy_data.get_simple_dpc_signal()
>>> s_new = s.rotate_beam_shifts(10)
>>> s_new.plot()
rotate_data(angle, reshape=False)[source]

Rotate the scan dimensions by angle.

Parameters:angle (float) – Clockwise rotation in degrees
Returns:rotated_signal
Return type:DPCSignal2D

Example

Rotate data by 10 degrees clockwise

>>> s = ps.dummy_data.get_simple_dpc_signal()
>>> s_rot = s.rotate_data(10)
>>> s_rot.plot()

Modules

IO tools

pixstem.io_tools.load_binary_merlin_signal(filename, probe_x=None, probe_y=None, chunks=(32, 32, 32, 32), flyback_pixels=1, datatype=None, lazy_result=True)[source]

Temporary function for loading Merlin binary data.

Parameters:
  • filename (string) –
  • probe_y (probe_x,) –
  • chunks (tuple) – Default (32, 32, 32, 32)
  • flyback_pixels (int) –
  • datatype (string) – 6 bit: “>u1”. 12 bit: “>u2”. 24 bit: “>u4”.
  • lazy_result (bool) –
  • function will be replaced at some point, so do not rely on it for (This) –
  • functions! (other) –
pixstem.io_tools.load_dpc_signal(filename)[source]

Load a differential phase contrast style signal.

This function can both files saved directly using HyperSpy, and saved using this library. The only requirement is that the signal has one navigation dimension, with this one dimension having a size of two. The first navigation index is the x-shift, while the second is the y-shift. The signal dimension contains the spatial dimension(s), i.e. the probe positions.

The return signal depends on the dimensions of the input file: - If two signal dimensions: DPCSignal2D - If one signal dimension: DPCSignal1D - If zero signal dimension: DPCBaseSignal

Parameters:filename (string) –
Returns:dpc_signal – The type of return signal depends on the signal dimensions of the input file.
Return type:DPCBaseSignal, DPCSignal1D, DPCSignal2D

Examples

>>> import numpy as np
>>> s = ps.DPCSignal2D(np.random.random((2, 90, 50)))
>>> s.save("test_dpc_signal2d.hspy", overwrite=True)
>>> s_dpc = ps.load_dpc_signal("test_dpc_signal2d.hspy")
>>> s_dpc
<DPCSignal2D, title: , dimensions: (2|50, 90)>
>>> s_dpc.plot()

Saving a HyperSpy signal

>>> import hyperspy.api as hs
>>> s = hs.signals.Signal1D(np.random.random((2, 10)))
>>> s.save("test_dpc_signal1d.hspy", overwrite=True)
>>> s_dpc_1d = ps.load_dpc_signal("test_dpc_signal1d.hspy")
>>> s_dpc_1d
<DPCSignal1D, title: , dimensions: (2|10)>
pixstem.io_tools.load_ps_signal(filename, lazy=False, chunk_size=None, navigation_signal=None)[source]
Parameters:
  • filename (string) –
  • lazy (bool, default False) –
  • chunk_size (tuple, optional) – Used if Lazy is True. Sets the chunk size of the signal. If it is not specified, the file chunking will be used. Higher number will potentially make the calculations be faster, but use more memory.
  • navigation_signal (Signal2D) –
pixstem.io_tools.signal_to_pixelated_stem(s)[source]

Make a PixelatedSTEM object from a HyperSpy signal.

This will retain both the axes information and the metadata. If the signal is lazy, the function will return LazyPixelatedSTEM.

Parameters:s (HyperSpy signal) – Should work for any HyperSpy signal.
Returns:pixelated_stem_signal
Return type:PixelatedSTEM or LazyPixelatedSTEM object

Examples

>>> import numpy as np
>>> import hyperspy.api as hs
>>> s = hs.signals.Signal2D(np.random.random((8, 11, 21, 13)))
>>> s.metadata.General.title = "test dataset"
>>> s
<Signal2D, title: test dataset, dimensions: (11, 8|13, 21)>
>>> from pixstem.io_tools import signal_to_pixelated_stem
>>> s_new = signal_to_pixelated_stem(s)
>>> s_new
<PixelatedSTEM, title: test dataset, dimensions: (11, 8|13, 21)>

Radial

pixstem.radial.fit_ellipses_to_signal(s, radial_signal_span_list, prepeak_range=None, angleN=20, show_progressbar=True)[source]
Parameters:
  • s (HyperSpy Signal2D) –
  • radial_signal_span_list (tuple) –
  • prepeak_range (tuple, optional) –
  • angleN (list or int, default 20) –
  • show_progressbar (bool, default True) –
Returns:

  • signal (HyperSpy 2D signal) – Fitted ellipse and fitting points are stored as HyperSpy markers in the metadata
  • xC, yC (floats) – Centre position
  • semi_len0, semi_len1 (floats) – Length of the two semi-axes.
  • rot (float) – Angle between semi_len0 and the positive x-axis. Since semi_len0, is not necessarily the longest semi-axis, the rotation will _not_ be between the major semi-axis and the positive x-axis. In radians, between 0 and pi. The rotation is clockwise, so at rot = 0.1 the ellipse will be pointing in the positive x-direction, and negative y-direction.
  • eccen (float) – Eccentricity of the ellipse

Examples

>>> import pixstem.make_diffraction_test_data as mdtd
>>> s = ps.PixelatedSTEM(np.zeros((200, 220)))
>>> s.axes_manager[0].offset, s.axes_manager[1].offset = -100, -110
>>> xx, yy = np.meshgrid(s.axes_manager[0].axis, s.axes_manager[1].axis)
>>> ellipse_ring = mdtd._get_elliptical_ring(xx, yy, 0, 0, 50, 70, 0.8)
>>> s.data += ellipse_ring
>>> from pixstem.radial import fit_ellipses_to_signal
>>> output = fit_ellipses_to_signal(
...     s, [(40, 80)], angleN=30, show_progressbar=False)
pixstem.radial.fit_single_ellipse_to_signal(s, radial_signal_span, prepeak_range=None, angleN=20, show_progressbar=True)[source]
Parameters:
  • s (HyperSpy Signal2D) –
  • radial_signal_span (tuple) –
  • prepeak_range (tuple, optional) –
  • angleN (int, default 20) –
  • show_progressbar (bool, default True) –
Returns:

  • signal (HyperSpy 2D signal) – Fitted ellipse and fitting points are stored as HyperSpy markers in the metadata
  • xC, yC (floats) – Centre position
  • semi_len0, semi_len1 (floats) – Length of minor and major semi-axis
  • rot (float) – Angle between semi_len0 and the positive x-axis. Since semi_len0, is not necessarily the longest semi-axis, the rotation will _not_ be between the major semi-axis and the positive x-axis. In radians, between 0 and pi. The rotation is clockwise, so at rot = 0.1 the ellipse will be pointing in the positive x-direction, and negative y-direction.
  • eccen (float) – Eccentricity of the ellipse
  • signal, xC, yC, semi0, semi1, rot, ecc

Examples

>>> import pixstem.make_diffraction_test_data as mdtd
>>> s = ps.PixelatedSTEM(np.zeros((200, 220)))
>>> s.axes_manager[0].offset, s.axes_manager[1].offset = -100, -110
>>> xx, yy = np.meshgrid(s.axes_manager[0].axis, s.axes_manager[1].axis)
>>> ellipse_ring = mdtd._get_elliptical_ring(xx, yy, 0, 0, 50, 70, 0.8)
>>> s.data += ellipse_ring
>>> from pixstem.radial import fit_single_ellipse_to_signal
>>> output = fit_single_ellipse_to_signal(
...     s, (40, 80), angleN=30, show_progressbar=False)
pixstem.radial.get_angle_image_comparison(s0, s1, angleN=12, mask_radius=None)[source]

Compare two images by overlaying one on the other in slices.

This function takes two images, extracts different angular slices and combines them into one image.

Useful for comparing two diffraction images, to see if the rings have the same radius.

Parameters:
  • s1 (s0,) – Both signals need to have the same shape, and no navigation dimensions.
  • angleN (int, default 12) – Number of angular slices.
  • mask_radius (int, optional) – Mask the centre of the image. The default is not to mask anything. Useful to mask the most intense parts of the diffraction pattern, so the less intense parts can be visualized.
Returns:

comparison_signal

Return type:

HyperSpy 2D

Examples

>>> from pixstem.make_diffraction_test_data import MakeTestData
>>> test_data0 = MakeTestData(300, 300)
>>> test_data0.add_ring(150, 150, 40)
>>> test_data1 = MakeTestData(300, 300)
>>> test_data1.add_ring(150, 150, 60)
>>> s0 = test_data0.signal
>>> s1 = test_data1.signal
>>> s0.axes_manager[0].offset, s0.axes_manager[1].offset = -150, -150
>>> s1.axes_manager[0].offset, s1.axes_manager[1].offset = -150, -150
>>> import pixstem.radial as ra
>>> s = ra.get_angle_image_comparison(s0, s1)
>>> s.plot()

Mask the inner parts

>>> s = ra.get_angle_image_comparison(s0, s1, mask_radius=10)
pixstem.radial.get_centre_position_list(s, steps, step_size)[source]

Returns a zip of x and y coordinates based on the offset of the center- point, and number of steps in each direction and a step-size given in pixels. Scale and offset taken from axes_manager.

pixstem.radial.get_coordinate_of_min(s)[source]

Returns the x and y values of the minimum in a signal.

pixstem.radial.get_optimal_centre_position(s, radial_signal_span, steps=3, step_size=1, angleN=8, show_progressbar=True)[source]

Find centre position of a ring by using angle sliced radial average.

Takes signal s, radial span of feature used to determine the centre position. Radially averages the feature in angleN segments, for each possible centre position. Models this averaged signal as a Gaussian and returns array of with the standard deviations of these Gaussians.

Note, the approximate centre position must be set using the offset parameter in the signal. The offset parameter is the negative of the position shown with s.plot(). So if the centre position in the image is x=52 and y=55: s.axes_manager.signal_axes[0].offset = -52 s.axes_manager.signal_axes[1].offset = -55 This can be checked by replotting the signal, to see if the centre position has the values x=0 and y=0.

Parameters:
  • s (HyperSpy 2D signal) – Only supports signals with no navigation dimensions.
  • radial_signal_span (tuple) – Range for finding the circular feature.
  • steps (int, default 3) – Number of steps in x/y direction to look for the optimal centre position. If the offset is (50, 55), and step_size 1: the positions x=(45, 46, 47, 48, …, 55) and y=(50, 51, 52, …, 60), will be checked.
  • step_size (int, default 1) –
  • angleN (int, default 8) – See angular_slice_radial_average for information about angleN.
Returns:

Optimal centre position

Return type:

HyperSpy 2D signal

Examples

>>> import pixstem.dummy_data as dd
>>> s = dd.get_single_ring_diffraction_signal()
>>> s.axes_manager.signal_axes[0].offset = -105
>>> s.axes_manager.signal_axes[1].offset = -67
>>> import pixstem.radial as ra
>>> s_ocp = ra.get_optimal_centre_position(s, (35, 45), steps=2)
>>> centre_pos = ra.get_coordinate_of_min(s_ocp)
pixstem.radial.get_radius_vs_angle(signal, radial_signal_span, angleN=15, centre_x=None, centre_y=None, prepeak_range=None, show_progressbar=True)[source]

Get radius of a ring as a function of angle.

This is done by radially averaging angular slices of the image, followed by firstly fitting the background around the specified ring, then fitting a Gaussian to the remaining intensity in the ring (which gets reduced to a peak due to the radial average).

Useful for finding if a ring is circular or elliptical. The centre position has be to set using the offset parameter in signal.axes_manager.

Parameters:
  • signal (HyperSpy 2D signal) – Diffraction image with calibrated centre point. This can be set manually, by looking at the diffraction pattern, or use a function like radial.get_optimal_centre_position
  • radial_signal_span (tuple) –
  • angleN (default 15) –
  • centre_y (centre_x,) –
  • prepeak_range (tuple, optional) –
  • show_progressbar (default True) –
Returns:

s_centre

Return type:

HyperSpy 1D signal

Examples

>>> import pixstem.dummy_data as dd
>>> s = dd.get_single_ring_diffraction_signal()
>>> s.axes_manager.signal_axes[0].offset = -105
>>> s.axes_manager.signal_axes[1].offset = -67
>>> import pixstem.radial as ra
>>> s_centre = ra.get_radius_vs_angle(s, (35, 45), show_progressbar=False)
pixstem.radial.refine_signal_centre_position(s, radial_signal_span, refine_step_size=None, **kwargs)[source]

Refine centre position of a diffraction signal by using round feature.

This function simply calls get_optimal_centre_position with a smaller and smaller step_size, giving a more accurate centre position.

The calculated centre position is stored in the offset value in the input signal.

Parameters:
  • s (HyperSpy 2D signal) – Approximate centre position must be set beforehand, see get_optimal_centre_position for more information.
  • radial_signal_span (tuple) –
  • refine_step_size (list, default [1, 0.5, 0.25]) –
  • **kwargs – Passed to get_optimal_centre_position

See also

get_optimal_centre_position()
finds the centre positions for a step_size

Example

>>> import pixstem.dummy_data as dd
>>> import pixstem.radial as ra
>>> s = dd.get_single_ring_diffraction_signal()
>>> s.axes_manager[0].offset, s.axes_manager[1].offset = -103, -69
>>> ra.refine_signal_centre_position(s, (32., 48.), angleN=4)
>>> x, y = s.axes_manager[0].offset, s.axes_manager[1].offset

Fluctuation electron microscopy analysis

pixstem.fem_tools.fem_calc(s, centre_x=None, centre_y=None, show_progressbar=True)[source]

Perform analysis of fluctuation electron microscopy (FEM) data as outlined in:

T. L. Daulton, et al., Ultramicroscopy 110 (2010) 1279-1289. doi:10.1016/j.ultramic.2010.05.010

Parameters:
  • s (PixelatedSTEM) – Signal on which FEM analysis was performed
  • centre_y (centre_x,) – All the diffraction patterns assumed to have the same centre position.
  • show_progressbar (bool) – Default True
Returns:

results – Results of FEM data analysis, including the normalized variance of the annular mean (V-Omegak), mean of normalized variances of rings (V-rk), normalized variance of ring ensemble (Vrek), the normalized variance image (Omega-Vi), and annular mean of the variance image (Omega-Vk).

Return type:

Python dictionary

Examples

>>> import pixstem.dummy_data as dd
>>> import pixstem.fem_tools as femt
>>> s = dd.get_fem_signal()
>>> fem_results = femt.fem_calc(
...     s,
...     centre_x=128,
...     centre_y=128,
...     show_progressbar=False)
>>> fem_results['V-Omegak'].plot()
pixstem.fem_tools.load_fem(rootname)[source]
Load individual HDF5 files containing previously saved FEM parameters to
dictionary.
Parameters:rootname (str) – Root name for output files. ‘_FEM_Parameter.hdf5’ will be appended
Returns:results – Series of HyperSpy signals previously saved using saveFEM
Return type:Dictionary

Examples

>>> import pixstem.fem_tools as femt
>>> s = ps.dummy_data.get_fem_signal()
>>> fem_results = s.fem_analysis(
...     centre_x=50,
...     centre_y=50,
...     show_progressbar=False)
>>> femt.save_fem(fem_results,'TestData')
>>> new_results = femt.load_fem('TestData')
pixstem.fem_tools.plot_fem(s, results, lowcutoff=10, highcutoff=120, k_cal=None)[source]

Produce a summary plot of all calculated FEM parameters

Parameters:
  • s (PixelatedSTEM) – Signal on which FEM analysis was performed
  • results (Dictionary) – Series of HyperSpy Signals containing results of FEM analysis performed on s
  • lowcutoff (integer) – Position of low-q cutoff for plots
  • highcutoff (integer) – Position of high-q cutoff for plots
  • k_cal (float or None) – Reciprocal space unit of length per pixel in inverse Angstroms
Returns:

fig

Return type:

Matplotlib Figure

Examples

>>> s = ps.dummy_data.get_fem_signal()
>>> import pixstem.fem_tools as femt
>>> fem_results = s.fem_analysis(
...     centre_x=50,
...     centre_y=50,
...     show_progressbar=False)
>>> fig = femt.plot_fem(s, fem_results, 10, 120)
pixstem.fem_tools.save_fem(results, rootname)[source]

Save dictionary of FEM results to individual HDF5 files.

Parameters:
  • results (Dictionary) – Results of FEM analysis performed on s
  • rootname (str) – Root name for output files. ‘_FEM_Parameter.hdf5’ will be appended

Examples

>>> import pixstem.fem_tools as femt
>>> s = ps.dummy_data.get_fem_signal()
>>> fem_results = s.fem_analysis(
...     centre_x=50,
...     centre_y=50,
...     show_progressbar=False)
>>> femt.save_fem(fem_results,'TestData')

Marker tools

pixstem.marker_tools.add_peak_array_to_signal_as_markers(signal, peak_array, color='red', size=20, bool_array=None, bool_invert=False)[source]

Add an array of points to a signal as HyperSpy markers.

Parameters:
  • signal (PixelatedSTEM or Signal2D) –
  • peak_array (4D NumPy array) –
  • color (string, optional) – Default ‘red’
  • size (scalar, optional) – Default 20
  • bool_array (NumPy array, optional) – Same shape as peaks_list.
  • bool_invert (bool, optional) – Default False.

Example

>>> s = ps.dummy_data.get_cbed_signal()
>>> peak_array = s.find_peaks(lazy_result=False, show_progressbar=False)
>>> import pixstem.marker_tools as mt
>>> mt.add_peak_array_to_signal_as_markers(s, peak_array)
>>> s.plot()

Dummy data

pixstem.dummy_data.get_cbed_signal()[source]

Get artificial pixelated STEM signal similar to CBED data.

Returns:cbed_signal
Return type:PixelatedSTEM signal

Example

>>> s = ps.dummy_data.get_cbed_signal()
>>> s.plot()
pixstem.dummy_data.get_dead_pixel_signal(lazy=False)[source]

Get 2D PixelatedSTEM signal with a disk in the middle.

Has 4 pixels with value equal to 0, to simulate dead pixels.

Example

>>> s = ps.dummy_data.get_dead_pixel_signal()

Lazy signal

>>> s_lazy = ps.dummy_data.get_dead_pixel_signal(lazy=True)
pixstem.dummy_data.get_disk_shift_simple_test_signal(lazy=False)[source]

Get HyperSpy 2D signal with 2D navigation dimensions for DPC testing.

Probe size x/y (20, 20), and image size x/y (50, 50). Disk moves from 22-28 x/y.

Parameters:lazy (bool, default False) –
Returns:disk_shift_signal
Return type:PixelatedSTEM signal

Examples

>>> s = ps.dummy_data.get_disk_shift_simple_test_signal()
>>> s.plot()

Load as lazy

>>> s = ps.dummy_data.get_disk_shift_simple_test_signal(lazy=True)
pixstem.dummy_data.get_fem_signal(lazy=False)[source]

Get a 2D signal that approximates a fluctuation electron microscopy (FEM) dataset.

Parameters:lazy (bool, default False) – If True, resulting signal will be lazy.
Returns:fem_signal
Return type:PixelatedSTEM

Examples

>>> s = ps.dummy_data.get_fem_signal()
>>> s.plot()
pixstem.dummy_data.get_generic_fem_signal(probe_x=2, probe_y=2, image_x=50, image_y=50, lazy=False)[source]

Get a 2D signal that approximates a fluctuation electron microscopy (FEM) dataset with user defined dimensions.

Parameters:
  • probe_x (int, default 2) – Horizontal dimension of the navigation axes
  • probe_y (int, default 2) – Vertical dimension of the navigation axes
  • image_x (int, default 2) – Horizontal dimension of the signal axes
  • image_y (int, default 2) – Vertical dimension of the signal axes
  • lazy (bool, default False) – If True, resulting signal will be lazy.
Returns:

fem_signal

Return type:

PixelatedSTEM

Examples

>>> s = ps.dummy_data.get_generic_fem_signal(probe_x=5, probe_y=10,
...     image_x=25, image_y=30, lazy=False)
>>> s.plot()
pixstem.dummy_data.get_holz_heterostructure_test_signal(lazy=False)[source]

Get HyperSpy 2D signal with 2D navigation dimensions for HOLZ testing.

The centre, radius and intensity of the ring varies as a function of probe position. The disk centre position varies as a function of probe position.

Parameters:lazy (bool, default False) –
Returns:holz_signal
Return type:HyperSpy 2D signal

Example

>>> s = ps.dummy_data.get_holz_heterostructure_test_signal()
>>> s.plot()

Load as lazy

>>> s = ps.dummy_data.get_holz_heterostructure_test_signal(lazy=True)
pixstem.dummy_data.get_holz_simple_test_signal(lazy=False)[source]

Get HyperSpy 2D signal with 2D navigation dimensions for HOLZ testing.

Probe size x/y (20, 20), and image size x/y (50, 50). Contains a disk and a ring. The disk stays at x, y = 25, 25, with radius 2. The ring has a radius of 20, and moves from x, y = 24-26, 24-26.

Parameters:lazy (bool, default False) –
Returns:holz_signal
Return type:PixelatedSTEM signal

Examples

>>> s = ps.dummy_data.get_holz_simple_test_signal()
>>> s.plot()

Load as lazy

>>> s = ps.dummy_data.get_holz_simple_test_signal(lazy=True)
pixstem.dummy_data.get_hot_pixel_signal(lazy=False)[source]

Get 2D PixelatedSTEM signal with a disk in the middle.

Has 4 pixels with value equal to 50000, to simulate hot pixels.

Example

>>> s = ps.dummy_data.get_hot_pixel_signal()

Lazy signal

>>> s_lazy = ps.dummy_data.get_hot_pixel_signal(lazy=True)
pixstem.dummy_data.get_nanobeam_electron_diffraction_signal()[source]

Get a signal emulating a NBED dataset.

Returns:signal
Return type:PixelatedSTEM

Example

>>> s = ps.dummy_data.get_nanobeam_electron_diffraction_signal()
>>> s.plot()
pixstem.dummy_data.get_simple_dpc_signal()[source]

Get a simple DPCSignal2D with a zero point in the centre.

Example

>>> s = ps.dummy_data.get_simple_dpc_signal()
pixstem.dummy_data.get_simple_ellipse_signal_peak_array()[source]

Get a signal and peak array of an ellipse.

Returns:signal, peak_array
Return type:HyperSpy Signal2D, NumPy array

Examples

>>> s, peak_array = ps.dummy_data.get_simple_ellipse_signal_peak_array()
>>> s.add_peak_array_as_markers(peak_array, color='blue', size=30)
pixstem.dummy_data.get_simple_fem_signal(lazy=False)[source]

Get a 2D signal that approximates a very small fluctuation electron microscopy (FEM) dataset.

Parameters:lazy (bool, default False) – If True, resulting signal will be lazy.
Returns:fem_signal
Return type:PixelatedSTEM

Examples

>>> s = ps.dummy_data.get_simple_fem_signal()
>>> s.plot()
pixstem.dummy_data.get_single_ring_diffraction_signal()[source]

Get HyperSpy 2D signal with a single ring with centre position.

The ring has a centre at x=105 and y=67, and radius=40.

pixstem.dummy_data.get_square_dpc_signal(add_ramp=False)[source]

Get a 2D DPC signal resembling a Landau domain.

Parameters:add_ramp (bool, default False) – If True, will add a ramp in the beam shift across the image to emulate the effects of d-scan.
Returns:square_dpc_signal
Return type:DPCSignal2D

Examples

>>> s = ps.dummy_data.get_square_dpc_signal()
>>> s.plot()

Adding a ramp

>>> s = ps.dummy_data.get_square_dpc_signal(add_ramp=True)
>>> s.plot()
pixstem.dummy_data.get_stripe_pattern_dpc_signal()[source]

Get a 2D DPC signal with a stripe pattern.

The stripe pattern only has an x-component, with alternating left/right directions. There is a small a net moment in the positive x-direction (leftwards).

Returns:stripe_dpc_signal
Return type:DPCSignal2D

Example

>>> s = ps.dummy_data.get_stripe_pattern_dpc_signal()

Make diffraction test data

class pixstem.make_diffraction_test_data.Circle(xx, yy, x0, y0, r, intensity, scale, lw=None)[source]
centre_on_image(xx, yy)[source]
get_centre_pixel(xx, yy, scale)[source]

This function sets the indices for the pixels on which the centre point is. Because the centre point can sometimes be exactly on the boundary of two pixels, the pixels are held in a list. One list for x (self.centre_x_pixels) and one for y (self.centre_x_pixels). If the centre is outside the image, the lists will be empty.

mask_outside_r(scale)[source]
set_uniform_intensity()[source]
update_axis(xx, yy)[source]
class pixstem.make_diffraction_test_data.DiffractionTestDataset(probe_x=10, probe_y=10, detector_x=256, detector_y=256, noise=0.5, dtype=<class 'numpy.float32'>)[source]
add_diffraction_image(diffraction_test_image, position_array=None)[source]

Add a diffraction image to all or a subset of the dataset.

See the class docstring for example on how to use this.

Parameters:
  • diffraction_test_image (PixStem DiffractionTestData object) –
  • position_array (NumPy array) – Boolean array, specifying which positions in the dataset the diffraction_test_image should be added to. Must have two dimensions, and the same shape as (probe_x, probe_y).
get_signal()[source]
plot()[source]
class pixstem.make_diffraction_test_data.DiffractionTestImage(disk_r=7, blur=2, image_x=256, image_y=256, rotation=0, diff_intensity_reduction=1.0, intensity_noise=0.5)[source]
add_background_lorentz(width=10, intensity=5)[source]
add_cubic_disks(vx, vy, intensity=1, n=1)[source]

Add disks in a cubic pattern around the centre point of the image.

Parameters:
  • vy (vx,) –
  • intensity (scalar) –
  • n (int) – Number of orders of diffraction disks. If n=1, 8 disks, if n=2, 24 disks.
add_disk(x, y, intensity=1)[source]
copy()[source]
get_diffraction_test_image(dtype=<class 'numpy.float32'>)[source]
get_signal()[source]
plot()[source]
class pixstem.make_diffraction_test_data.Disk(xx, yy, scale, x0, y0, r, intensity)[source]

Disk object, with outer edge of the ring at r

get_signal()[source]
set_centre_intensity()[source]

Sets the intensity of the centre pixels to I. Coordinates are self.z.circle[y, x], due to how numpy works.

update_axis(xx, yy)[source]
class pixstem.make_diffraction_test_data.EllipseDisk(xx, yy, x0, y0, semi_len0, semi_len1, rotation, intensity)[source]
get_signal()[source]
class pixstem.make_diffraction_test_data.EllipseRing(xx, yy, x0, y0, semi_len0, semi_len1, rotation, intensity, lw_r=1)[source]
get_signal()[source]
class pixstem.make_diffraction_test_data.MakeTestData(size_x=100, size_y=100, scale=1, default=False, blur=True, blur_sigma=1, downscale=True)[source]

MakeTestData is an object containing a generated test signal.

The default signal consists of a Disk and concentric Ring, with the Ring being less intensive than the centre Disk. Unlimited number of Rings and Disk can be added separately.

Parameters:
  • size_y (size_x,) – The range of the x and y axis goes from 0 to size_x, size_y
  • scale (float, int) – The step size of the x and y axis
  • default (bool, default False) – If True, the default object should be generated. If false, Ring and Disk must be added separately by self.add_ring(), self.add_disk()
  • blur (bool, default True) – If True, do a Gaussian blur of the disk.
  • blur_sigma (int, default 1) – Sigma of the Gaussian blurring, if blur is True.
  • downscale (bool, default True) – Note: currently using downscale and adding a disk, will lead to the center of the disk being shifted. Ergo: the disk_x and disk_y will not be correct when downscale is True.
signal

Test signal

Type:hyperspy.signals.Signal2D
z_list

List containing Ring and Disk objects added to the signal

Type:list
downscale_factor

The data is upscaled before Circle is added, and similarly downscaled to return to given dimensions. This improves the quality of Circle

Type:int

Examples

Default settings

>>> from pixstem.make_diffraction_test_data import MakeTestData
>>> test_data = MakeTestData()
>>> test_data.signal.plot()

More control

>>> test_data = MakeTestData(default=False)
>>> test_data.add_disk(x0=50, y0=50, r=10, intensity=30)
>>> test_data.add_ring(x0=45, y0=52, r=25, intensity=10)
>>> test_data.signal.plot()
add_disk(x0=50, y0=50, r=5, intensity=10)[source]
add_disk_ellipse(x0=50, y0=50, semi_len0=5, semi_len1=8, rotation=0.78, intensity=10)[source]
add_ring(x0=50, y0=50, r=20, intensity=10, lw_pix=0)[source]

Add a ring to the test data.

Parameters:
  • y0 (x0,) – Centre position of the ring
  • r (number, default 20) – Radius of the ring, defined as the distance from the centre to the middle of the line of the ring, which will be most intense after blurring.
  • intensity (number, default 10) – Pixel value of the ring. Note, this value will be lowered if blur or downscale is True
  • lw_pix (number, default 0) – Distance in pixels from radius to the outer and inner edge of the ring. Inner radius: r-lw, outer radius: r+lw. In total this gives a ring line width in pixels of 2*lw+1.
add_ring_ellipse(x0=50, y0=50, semi_len0=11, semi_len1=13, rotation=0.78, intensity=10, lw_r=2)[source]
blur()[source]
downscale()[source]
generate_mesh()[source]
make_signal()[source]
set_downscale_factor(factor)[source]
set_signal_zero()[source]
to_signal()[source]
update_signal()[source]
class pixstem.make_diffraction_test_data.Ring(xx, yy, scale, x0, y0, r, intensity, lr)[source]

Ring object, with outer edge of the ring at r+lr, and inner r-lr.

The radius of the ring is defined as in the middle of the line making up the ring.

get_signal()[source]
mask_inside_r(scale)[source]
update_axis(xx, yy)[source]
pixstem.make_diffraction_test_data.generate_4d_data(probe_size_x=10, probe_size_y=10, image_size_x=50, image_size_y=50, disk_x=25, disk_y=25, disk_r=5, disk_I=20, ring_x=25, ring_y=25, ring_r=20, ring_I=6, ring_lw=0, ring_e_x=None, ring_e_y=25, ring_e_semi_len0=15, ring_e_semi_len1=15, ring_e_r=0, ring_e_I=6, ring_e_lw=1, blur=True, blur_sigma=1, downscale=True, add_noise=False, noise_amplitude=1, lazy=False, lazy_chunks=None, show_progressbar=True)[source]

Generate a test dataset containing a disk and diffraction ring.

Useful for checking that radial average algorithms are working properly.

The centre, intensity and radius position of the ring and disk can vary as a function of probe position, through the disk_x, disk_y, disk_r, disk_I, ring_x, ring_y, ring_r and ring_I arguments. In addition, the line width of the ring can be varied with ring_lw.

There is also an elliptical ring, which can be added separately to the circular ring. This elliptical ring uses the ring_e_* arguments. It is disabled by default.

The ring can be deactivated by setting ring_x=None. The disk can be deactivated by setting disk_x=None. The elliptical ring can be deactivated by setting ring_e_x=None.

Parameters:
  • probe_size_y (probe_size_x,) – Size of the navigation dimension.
  • image_size_y (image_size_x,) – Size of the signal dimension.
  • disk_y (disk_x,) – Centre position of the disk. Either integer or NumPy 2-D array. See examples on how to make them the correct size. To deactivate the disk, set disk_x=None.
  • disk_r (int or NumPy 2D-array, default 5) – Radius of the disk. Either integer or NumPy 2-D array. See examples on how to make it the correct size.
  • disk_I (int or NumPy 2D-array, default 20) – Intensity of the disk, for each of the pixels. So if I=30, the each pixel in the disk will have a value of 30. Note, this value will change if blur=True or downscale=True.
  • ring_y (ring_x,) – Centre position of the ring. Either integer or NumPy 2-D array. See examples on how to make them the correct size. To deactivate the ring, set ring_x=None.
  • ring_r (int or NumPy 2D-array, default 20) – Radius of the ring. Either integer or NumPy 2-D array. See examples on how to make it the correct size.
  • ring_I (int or NumPy 2D-array, default 6) – Intensity of the ring, for each of the pixels. So if I=5, each pixel in the ring will have a value of 5. Note, this value will change if blur=True or downscale=True.
  • ring_lw (int or NumPy 2D-array, default 0) – Line width of the ring. If ring_lw=1, the line will be 3 pixels wide. If ring_lw=2, the line will be 5 pixels wide.
  • ring_e_y (ring_e_x,) – Centre position of the elliptical ring. Either integer or NumPy 2-D array. See examples on how to make them the correct size. To deactivate the ring, set ring_x=None (which is the default).
  • ring_e_semi_len1 (ring_e_semi_len0,) – Semi lengths of the elliptical ring. Either integer or NumPy 2-D arrays. See examples on how to make it the correct size.
  • ring_e_I (int or NumPy 2D-array, default 6) – Intensity of the elliptical ring, for each of the pixels. So if I=5, each pixel in the ring will have a value of 5. Note, this value will change if blur=True or downscale=True.
  • ring_e_r (int or NumPy 2D-array, default 0) – Rotation of the elliptical ring, in radians.
  • ring_e_lw (int or NumPy 2D-array, default 0) – Line width of the ring. If ring_lw=1, the line will be 3 pixels wide. If ring_lw=2, the line will be 5 pixels wide.
  • blur (bool, default True) – If True, do a Gaussian blur of the disk.
  • blur_sigma (int, default 1) – Sigma of the Gaussian blurring, if blur is True.
  • downscale (bool, default True) – If True, use upscaling (then downscaling) to anti-alise the disk.
  • add_noise (bool, default False) – Add Gaussian random noise.
  • noise_amplitude (float, default 1) – The amplitude of the noise, if add_noise is True.
  • lazy (bool, default False) – If True, the signal will be lazy
  • lazy_chunks (tuple, optional) – Used if lazy is True, default (10, 10, 10, 10).
Returns:

signal – Signal with 2 navigation dimensions and 2 signal dimensions.

Return type:

HyperSpy Signal2D

Examples

>>> import pixstem.make_diffraction_test_data as mdtd
>>> s = mdtd.generate_4d_data(show_progressbar=False)
>>> s.plot()

Using more arguments

>>> s = mdtd.generate_4d_data(probe_size_x=20, probe_size_y=30,
...         image_size_x=50, image_size_y=90,
...         disk_x=30, disk_y=70, disk_r=9, disk_I=30,
...         ring_x=35, ring_y=65, ring_r=20, ring_I=10,
...         blur=False, downscale=False, show_progressbar=False)

Adding some Gaussian random noise

>>> s = mdtd.generate_4d_data(add_noise=True, noise_amplitude=3,
...         show_progressbar=False)

Different centre positions for each probe position. Note the size=(20, 10), and probe_x=10, probe_y=20: size=(y, x).

>>> import numpy as np
>>> disk_x = np.random.randint(5, 35, size=(20, 10))
>>> disk_y = np.random.randint(5, 45, size=(20, 10))
>>> disk_I = np.random.randint(50, 100, size=(20, 10))
>>> ring_x = np.random.randint(5, 35, size=(20, 10))
>>> ring_y = np.random.randint(5, 45, size=(20, 10))
>>> ring_r = np.random.randint(10, 15, size=(20, 10))
>>> ring_I = np.random.randint(1, 30, size=(20, 10))
>>> ring_lw = np.random.randint(1, 5, size=(20, 10))
>>> s = mdtd.generate_4d_data(probe_size_x=10, probe_size_y=20,
...         image_size_x=40, image_size_y=50, disk_x=disk_x, disk_y=disk_y,
...         disk_I=disk_I, ring_x=ring_x, ring_y=ring_y, ring_r=ring_r,
...         ring_I=ring_I, ring_lw=ring_lw, show_progressbar=False)

Do not plot the disk

>>> s = mdtd.generate_4d_data(disk_x=None, show_progressbar=False)

Do not plot the ring

>>> s = mdtd.generate_4d_data(ring_x=None, show_progressbar=False)

Plot only an elliptical ring

>>> from numpy.random import randint, random
>>> s = mdtd.generate_4d_data(
...        probe_size_x=10, probe_size_y=10,
...        disk_x=None, ring_x=None,
...        ring_e_x=randint(20, 30, (10, 10)),
...        ring_e_y=randint(20, 30, (10, 10)),
...        ring_e_semi_len0=randint(10, 20, (10, 10)),
...        ring_e_semi_len1=randint(10, 20, (10, 10)),
...        ring_e_r=random((10, 10))*np.pi,
...        ring_e_lw=randint(1, 3, (10, 10)))