API reference¶
This reference manual details the public classes, modules and functions in kikuchipy, as generated from their docstrings. Many of the docstrings contain examples, however, see the user guide for how to use kikuchipy.
Caution
kikuchipy is in an alpha stage, and there will likely be breaking changes with each release.
Signals¶
EBSD¶
All methods listed here are also available to
LazyEBSD
objects.
Enhance the local contrast in an EBSD scan inplace using adaptive histogram equalization. |
|
|
Average patterns in an EBSD scan inplace with its neighbours within a window. |
|
Filter an EBSD scan inplace in the frequency domain. |
|
Get the model signal generated with the selected number of principal components from a decomposition. |
|
Get the dynamic background per EBSD pattern in a scan. |
|
Compute the image quality map of patterns in an EBSD scan. |
|
Get a virtual backscatter electron (VBSE) image formed from intensities within a region of interest (ROI) on the detector. |
|
Normalize image intensities in inplace to a mean of zero with a given standard deviation. |
|
Plot an interactive virtual backscatter electron (VBSE) image formed from intensities within a specified and adjustable region of interest (ROI) on the detector. |
|
Rebin the signal into a smaller or larger shape, based on linear interpolation. |
|
Remove the dynamic background in an EBSD scan inplace. |
|
Remove the static background in an EBSD scan inplace. |
|
Rescale image intensities inplace. |
|
Write the signal to file in the specified format. |
|
Set detector pixel size in microns. |
|
Set experimental parameters in signal metadata. |
|
Set parameters for one phase in signal metadata. |
|
Set the step size in microns. |
-
class
kikuchipy.signals.
EBSD
(*args, **kwargs)[source]¶ Bases:
kikuchipy.signals._common_image.CommonImage
,hyperspy._signals.signal2d.Signal2D
Scan of Electron Backscatter Diffraction (EBSD) patterns.
This class extends HyperSpy’s Signal2D class for EBSD patterns, with common intensity processing methods and some analysis methods.
Methods inherited from HyperSpy can be found in the HyperSpy user guide.
See the docstring of
hyperspy.signal.BaseSignal
for a list of attributes.-
adaptive_histogram_equalization
(kernel_size=None, clip_limit=0, nbins=128)[source]¶ Enhance the local contrast in an EBSD scan inplace using adaptive histogram equalization.
This method uses
skimage.exposure.equalize_adapthist()
.- Parameters
kernel_size (
Union
[Tuple
[int
,int
],List
[int
],None
]) – Shape of contextual regions for adaptive histogram equalization, default is 1/4 of image height and 1/4 of image width.clip_limit (
Union
[int
,float
]) – Clipping limit, normalized between 0 and 1 (higher values give more contrast). Default is 0.nbins (
int
) – Number of gray bins for histogram (“data range”), default is 128.
Examples
To best understand how adaptive histogram equalization works, we plot the histogram of the same image before and after equalization:
>>> import numpy as np >>> import matplotlib.pyplot as plt >>> s2 = s.inav[0, 0] >>> s2.adaptive_histogram_equalization() >>> imin = np.iinfo(s.data.dtype_out).min >>> imax = np.iinfo(s.data.dtype_out).max + 1 >>> hist, _ = np.histogram( ... s.inav[0, 0].data, bins=imax, range=(imin, imax)) >>> hist2, _ = np.histogram( ... s2.inav[0, 0].data, bins=imax, range=(imin, imax)) >>> fig, ax = plt.subplots(nrows=2, ncols=2) >>> ax[0, 0].imshow(s.inav[0, 0].data) >>> ax[1, 0].plot(hist) >>> ax[0, 1].imshow(s2.inav[0, 0].data) >>> ax[1, 1].plot(hist2)
Notes
It is recommended to perform adaptive histogram equalization only after static and dynamic background corrections, otherwise some unwanted darkening towards the edges might occur.
The default window size might not fit all image sizes, so it may be necessary to search for the optimal window size.
-
average_neighbour_patterns
(window='circular', window_shape=3, 3, **kwargs)[source]¶ Average patterns in an EBSD scan inplace with its neighbours within a window.
The amount of averaging is specified by the window coefficients. All patterns are averaged with the same window. Map borders are extended with zeros.
Averaging is accomplished by correlating the window with the extended array of patterns using
scipy.ndimage.correlate()
.- Parameters
window (
Union
[str
,ndarray
,Array
,Window
]) – Name of averaging window or an array. Available types are listed inscipy.signal.windows.get_window()
, in addition to a “circular” window (default) filled with ones in which corner coefficients are set to zero. A window element is considered to be in a corner if its radial distance to the origin (window centre) is shorter or equal to the half width of the window’s longest axis. A 1D or 2Dnumpy.ndarray
,dask.array.Array
orWindow
can also be passed.window_shape (
Tuple
[int
, …]) – Shape of averaging window. Not used if a custom window orWindow
object is passed to window. This can be either 1D or 2D, and can be asymmetrical. Default is (3, 3).**kwargs – Keyword arguments passed to the available window type listed in
scipy.signal.windows.get_window()
. If none are passed, the default values of that particular window are used.
Examples
>>> import numpy as np >>> import kikuchipy as kp >>> s = kp.signals.EBSD(np.ones((4, 4, 1, 1))) >>> k = 1 >>> for i in range(4): ... for j in range(4): ... s.inav[j, i].data *= k ... k += 1 >>> s.data[:, :, 0, 0] array([[ 1., 2., 3., 4.], [ 5., 6., 7., 8.], [ 9., 10., 11., 12.], [13., 14., 15., 16.]]) >>> s2 = s.deepcopy() # For use below >>> s.average_neighbour_patterns( ... window="circular", shape=(3, 3)) >>> s.data[:, :, 0, 0] array([[ 2.66666667, 3. , 4. , 5. ], [ 5.25 , 6. , 7. , 7.75 ], [ 9.25 , 10. , 11. , 11.75 ], [12. , 13. , 14. , 14.33333333]])
A window object can also be passed
>>> w = kp.filters.Window(window="gaussian", std=2) >>> w Window (3, 3) gaussian [[0.7788 0.8825 0.7788] [0.8825 1. 0.8825] [0.7788 0.8825 0.7788]] >>> s2.average_neighbour_patterns(w) >>> s2.data[:, :, 0, 0] array([[ 3.34395304, 3.87516243, 4.87516264, 5.40637176], [ 5.46879047, 5.99999985, 6.99999991, 7.53120901], [ 9.46879095, 9.99999959, 11.00000015, 11.53120913], [11.59362845, 12.12483732, 13.1248368 , 13.65604717]])
This window can subsequently be plotted and saved
>>> figure, image, colorbar = w.plot() >>> figure.savefig('averaging_window.png')
-
fft_filter
(transfer_function, function_domain, shift=False)[source]¶ Filter an EBSD scan inplace in the frequency domain.
Patterns are transformed via the Fast Fourier Transform (FFT) to the frequency domain, where their spectrum is multiplied by the transfer_function, and the filtered spectrum is subsequently transformed to the spatial domain via the inverse FFT (IFFT). Filtered patterns are rescaled to input data type range.
Note that if function_domain is “spatial”, only real valued FFT and IFFT is used.
- Parameters
transfer_function (
Union
[ndarray
,Window
]) – Filter to apply to patterns. This can either be a transfer function in the frequency domain of pattern shape or a kernel in the spatial domain. What is passed is determined from function_domain.function_domain (
str
) – Options are “frequency” and “spatial”, indicating, respectively, whether the filter function passed to filter_function is a transfer function in the frequency domain or a kernel in the spatial domain.shift (
bool
) – Whether to shift the zero-frequency component to the centre. Default is False. This is only used when function_domain=”frequency”.
Examples
Applying a Gaussian low pass filter with a cutoff frequency of 20 to an EBSD object
s
:>>> pattern_shape = s.axes_manager.signal_shape[::-1] >>> w = kp.filters.Window( ... "lowpass", cutoff=20, shape=pattern_shape) >>> s.fft_filter( ... transfer_function=w, ... function_domain="frequency", ... shift=True, ... )
See also
-
get_decomposition_model
(components=None, dtype_out=<class 'numpy.float32'>)[source]¶ Get the model signal generated with the selected number of principal components from a decomposition.
Calls HyperSpy’s
hyperspy.learn.mva.MVA.get_decomposition_model()
. Learning results are preconditioned before this call, doing the following: (1) setnumpy.dtype
to desired dtype_out, (2) remove unwanted components, and (3) rechunk, ifdask.array.Array
, to suitable chunks.- Parameters
components (
Union
[None
,int
,List
[int
]]) – If None (default), rebuilds the signal from all components. If int, rebuilds signal from components in range 0-given int. If list of ints, rebuilds signal from only components in given list.dtype_out (
dtype
) – Data type to cast learning results to (default isnumpy.float32
). Note that HyperSpy casts them tonumpy.float64
.
- Returns
s_model
- Return type
-
get_dynamic_background
(filter_domain='frequency', std=None, truncate=4.0, dtype_out=None, **kwargs)[source]¶ Get the dynamic background per EBSD pattern in a scan.
- Parameters
filter_domain (
str
) – Whether to apply a Gaussian convolution filter in the “frequency” (default) or “spatial” domain.std (
Union
[None
,int
,float
]) – Standard deviation of the Gaussian window. If None (default), it is set to width/8.truncate (
Union
[int
,float
]) – Truncate the Gaussian filter at this many standard deviations. Default is 4.0.dtype_out (
Optional
[dtype
]) – Data type of the background patterns. If None (default), it is set to the same data type as the input pattern.kwargs – Keyword arguments passed to the Gaussian blurring function determined from filter_domain.
- Returns
background_signal – Signal with the large scale variations across the detector.
- Return type
-
get_image_quality
(normalize=True)[source]¶ Compute the image quality map of patterns in an EBSD scan.
The image quality is calculated based on the procedure defined by Krieger Lassen [Lassen1994].
- Parameters
normalize (
bool
) – Whether to normalize patterns to a mean of zero and standard deviation of 1 before calculating the image quality. Default is True.- Returns
image_quality_map – Image quality map of same shape as signal navigation axes.
- Return type
References
- Lassen1994(1,2,3)
Lassen, “Automated Determination of Crystal Orientations from Electron Backscattering Patterns,” Institute of Mathematical Modelling, (1994).
Examples
>>> iq = s.get_image_quality(normalize=True) # Default >>> plt.imshow(iq)
-
get_virtual_bse_intensity
(roi, out_signal_axes=None)[source]¶ Get a virtual backscatter electron (VBSE) image formed from intensities within a region of interest (ROI) on the detector.
Adapted from
pyxem.signals.common_diffraction.CommonDiffraction.get_integrated_intensity()
.- Parameters
- Returns
virtual_image – VBSE image formed from detector intensities within an ROI on the detector.
- Return type
Examples
>>> import hyperspy.api as hs >>> roi = hs.roi.RectangularROI( ... left=0, right=5, top=0, bottom=5) >>> vbse_image = s.get_virtual_bse_intensity(roi)
-
normalize_intensity
(num_std=1, divide_by_square_root=False, dtype_out=None)[source]¶ Normalize image intensities in inplace to a mean of zero with a given standard deviation.
- Parameters
num_std (
int
) – Number of standard deviations of the output intensities. Default is 1.divide_by_square_root (
bool
) – Whether to divide output intensities by the square root of the signal dimension size. Default is False.dtype_out (
Optional
[dtype
]) – Data type of normalized images. If None (default), the input images’ data type is used.
Notes
Data type should always be changed to floating point, e.g.
np.float32
withchange_dtype()
, before normalizing the intensities.Examples
>>> np.mean(s.data) 146.0670987654321 >>> s.change_dtype(np.float32) # Or passing dtype_out=np.float32 >>> s.normalize_intensity() >>> np.mean(s.data) 2.6373216e-08
Notes
Rescaling RGB images is not possible. Use RGB channel normalization when creating the image instead.
-
plot_virtual_bse_intensity
(roi, out_signal_axes=None, **kwargs)[source]¶ Plot an interactive virtual backscatter electron (VBSE) image formed from intensities within a specified and adjustable region of interest (ROI) on the detector.
Adapted from
pyxem.signals.common_diffraction.CommonDiffraction.plot_integrated_intensity()
.- Parameters
roi (
BaseInteractiveROI
) – Any interactive ROI detailed in HyperSpy.out_signal_axes (
Union
[None
,Iterable
[int
],Iterable
[str
]]) – Which navigation axes to use as signal axes in the virtual image. If None (default), the first two navigation axes are used.**kwargs – Keyword arguments passed to the plot method of the virtual image.
Examples
>>> import hyperspy.api as hs >>> roi = hs.roi.RectangularROI( ... left=0, right=5, top=0, bottom=5) >>> s.plot_virtual_bse_intensity(roi)
-
rebin
(new_shape=None, scale=None, crop=True, out=None)[source]¶ Rebin the signal into a smaller or larger shape, based on linear interpolation. Specify either new_shape or scale.
- Parameters
new_shape (list (of floats or integer) or None) – For each dimension specify the new_shape. This will internally be converted into a scale parameter.
scale (list (of floats or integer) or None) – For each dimension, specify the new:old pixel ratio, e.g. a ratio of 1 is no binning and a ratio of 2 means that each pixel in the new spectrum is twice the size of the pixels in the old spectrum. The length of the list should match the dimension of the Signal’s underlying data array. Note : Only one of `scale` or `new_shape` should be specified, otherwise the function will not run
crop (bool) –
Whether or not to crop the resulting rebinned data (default is
True
). When binning by a non-integer number of pixels it is likely that the final row in each dimension will contain fewer than the full quota to fill one pixel.e.g. a 5*5 array binned by 2.1 will produce two rows containing 2.1 pixels and one row containing only 0.8 pixels. Selection of
crop=True
orcrop=False
determines whether or not this “black” line is cropped from the final binned array or not.
Please note that if
crop=False
is used, the final row in each dimension may appear black if a fractional number of pixels are left over. It can be removed but has been left to preserve total counts before and after binning.out (
BaseSignal
(or subclasses) orNone
) – IfNone
, a new Signal is created with the result of the operation and returned (default). If a Signal is passed, it is used to receive the output of the operation, and nothing is returned.
- Returns
s – The resulting cropped signal.
- Return type
BaseSignal
(or subclass)
Examples
>>> spectrum = hs.signals.EDSTEMSpectrum(np.ones([4, 4, 10])) >>> spectrum.data[1, 2, 9] = 5 >>> print(spectrum) <EDXTEMSpectrum, title: dimensions: (4, 4|10)> >>> print ('Sum = ', sum(sum(sum(spectrum.data)))) Sum = 164.0 >>> scale = [2, 2, 5] >>> test = spectrum.rebin(scale) >>> print(test) <EDSTEMSpectrum, title: dimensions (2, 2|2)> >>> print('Sum = ', sum(sum(sum(test.data)))) Sum = 164.0
-
remove_dynamic_background
(operation='subtract', filter_domain='frequency', std=None, truncate=4.0, **kwargs)[source]¶ Remove the dynamic background in an EBSD scan inplace.
The removal is performed by subtracting or dividing by a Gaussian blurred version of each pattern. Resulting pattern intensities are rescaled to fill the input patterns’ data type range.
- Parameters
operation (
str
) – Whether to “subtract” (default) or “divide” by the dynamic background pattern.filter_domain (
str
) – Whether to obtain the dynamic background by applying a Gaussian convolution filter in the “frequency” (default) or “spatial” domain.std (
Union
[None
,int
,float
]) – Standard deviation of the Gaussian window. If None (default), it is set to width/8.truncate (
Union
[int
,float
]) – Truncate the Gaussian window at this many standard deviations. Default is 4.0.kwargs – Keyword arguments passed to the Gaussian blurring function determined from filter_domain.
See also
kikuchipy.signals.EBSD.remove_static_background()
,kikuchipy.signals.EBSD.get_dynamic_background()
,kikuchipy.pattern.remove_dynamic_background()
,kikuchipy.pattern.get_dynamic_background()
Examples
Traditional background correction includes static and dynamic corrections, loosing relative intensities between patterns after dynamic corrections (whether relative is set to True or False in
remove_static_background()
):>>> s.remove_static_background(operation="subtract") >>> s.remove_dynamic_background( ... operation="subtract", # Default ... filter_domain="frequency", # Default ... truncate=4.0, # Default ... std=5, ... )
-
remove_static_background
(operation='subtract', relative=True, static_bg=None, scale_bg=False)[source]¶ Remove the static background in an EBSD scan inplace.
The removal is performed by subtracting or dividing by a static background pattern. Resulting pattern intensities are rescaled keeping relative intensities or not and stretched to fill the available grey levels in the patterns’ data type range.
- Parameters
operation (
str
) – Whether to “subtract” (default) or “divide” by the static background pattern.relative (
bool
) – Keep relative intensities between patterns. Default is True.static_bg (
Union
[None
,ndarray
,Array
]) – Static background pattern. If None is passed (default) we try to read it from the signal metadata.scale_bg (
bool
) – Whether to scale the static background pattern to each individual pattern’s data range before removal. Must be False if relative is True. Default is False.
Examples
We assume that a static background pattern with the same shape and data type (e.g. 8-bit unsigned integer,
uint8
) as the patterns is available in signal metadata:>>> import kikuchipy as kp >>> ebsd_node = kp.signals.util.metadata_nodes("ebsd") >>> s.metadata.get_item(ebsd_node + '.static_background') [[84 87 90 ... 27 29 30] [87 90 93 ... 27 28 30] [92 94 97 ... 39 28 29] ... [80 82 84 ... 36 30 26] [79 80 82 ... 28 26 26] [76 78 80 ... 26 26 25]]
The static background can be removed by subtracting or dividing this background from each pattern while keeping relative intensities between patterns (or not):
>>> s.remove_static_background( ... operation='subtract', relative=True)
If the metadata has no background pattern, this must be passed in the static_bg parameter as a numpy or dask array.
-
rescale_intensity
(relative=False, in_range=None, out_range=None, dtype_out=None, percentiles=None)[source]¶ Rescale image intensities inplace.
Output min./max. intensity is determined from out_range or the data type range of the
numpy.dtype
passed to dtype_out if out_range is None.This method is based on
skimage.exposure.rescale_intensity()
.- Parameters
relative (
bool
) – Whether to keep relative intensities between images (default is False). If True, in_range must be None, because in_range is in this case set to the global min./max. intensity.in_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of input images. If None (default), stretching is performed when in_range is set to a narrower in_range is set to pattern min./max intensity. Contrast intensity range than the input patterns. Must be None if relative is True or percentiles are passed.out_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of output images. If None (default), out_range is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of rescaled images, default is input images’ data type.percentiles (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Disregard intensities outside these percentiles. Calculated per image. Must be None if in_range or relative is passed. Default is None.
See also
kikuchipy.signals.CommonImage.normalize_intensity()
,kikuchipy.pattern.rescale_intensity()
,skimage.exposure.rescale_intensity()
Examples
Image intensities are stretched to fill the available grey levels in the input images’ data type range or any
numpy.dtype
range passed to dtype_out, either keeping relative intensities between images or not:>>> print(s.data.dtype_out, s.data.min(), s.data.max(), ... s.inav[0, 0].data.min(), s.inav[0, 0].data.max()) uint8 20 254 24 233 >>> s2 = s.deepcopy() >>> s.rescale_intensity(dtype_out=np.uint16) >>> print(s.data.dtype_out, s.data.min(), s.data.max(), ... s.inav[0, 0].data.min(), s.inav[0, 0].data.max()) uint16 0 65535 0 65535 >>> s2.rescale_intensity(relative=True) >>> print(s2.data.dtype_out, s2.data.min(), s2.data.max(), ... s2.inav[0, 0].data.min(), s2.inav[0, 0].data.max()) uint8 0 255 4 232
Contrast stretching can be performed by passing percentiles:
>>> s.rescale_intensity(percentiles=(1, 99))
Here, the darkest and brightest pixels within the 1% percentile are set to the ends of the data type range, e.g. 0 and 255 respectively for images of
uint8
data type.Notes
Rescaling RGB images is not possible. Use RGB channel normalization when creating the image instead.
-
save
(filename=None, overwrite=None, extension=None, **kwargs)[source]¶ Write the signal to file in the specified format.
The function gets the format from the extension: h5, hdf5 or h5ebsd for kikuchipy’s specification of the the h5ebsd format, dat for the NORDIF binary format or hspy for HyperSpy’s HDF5 specification. If no extension is provided the signal is written to a file in kikuchipy’s h5ebsd format. Each format accepts a different set of parameters.
For details see the specific format documentation under “See Also” below.
This method is a modified version of HyperSpy’s function
hyperspy.signal.BaseSignal.save()
.- Parameters
filename (
Optional
[str
]) – If None (default) and tmp_parameters.filename and tmp_parameters.folder in signal metadata are defined, the filename and path will be taken from there. A valid extension can be provided e.g. “data.h5”, see extension.overwrite (
Optional
[bool
]) – If None and the file exists, it will query the user. If True (False) it (does not) overwrite the file if it exists.extension (
Optional
[str
]) – Extension of the file that defines the file format. Options are “h5”/”hdf5”/”h5ebsd”/”dat”/”hspy”. “h5”/”hdf5”/”h5ebsd” are equivalent. If None, the extension is determined from the following list in this order: i) the filename, ii) tmp_parameters.extension or iii) “h5” (kikuchipy’s h5ebsd format).**kwargs – Keyword arguments passed to writer.
-
set_detector_calibration
(delta)[source]¶ Set detector pixel size in microns. The offset is set to the the detector centre.
Examples
>>> s.axes_manager['dx'].scale # Default value 1.0 >>> s.set_detector_calibration(delta=70.) >>> s.axes_manager['dx'].scale 70.0
-
set_experimental_parameters
(detector=None, azimuth_angle=None, elevation_angle=None, sample_tilt=None, working_distance=None, binning=None, exposure_time=None, grid_type=None, gain=None, frame_number=None, frame_rate=None, scan_time=None, beam_energy=None, xpc=None, ypc=None, zpc=None, static_background=None, manufacturer=None, version=None, microscope=None, magnification=None)[source]¶ Set experimental parameters in signal metadata.
- Parameters
azimuth_angle (float, optional) – Azimuth angle of the detector in degrees. If the azimuth is zero, the detector is perpendicular to the tilt axis.
beam_energy (float, optional) – Energy of the electron beam in kV.
binning (int, optional) – Camera binning.
detector (str, optional) – Detector manufacturer and model.
elevation_angle (float, optional) – Elevation angle of the detector in degrees. If the elevation is zero, the detector is perpendicular to the incident beam.
exposure_time (float, optional) – Camera exposure time in µs.
frame_number (float, optional) – Number of patterns integrated during acquisition.
frame_rate (float, optional) – Frames per s.
gain (float, optional) – Camera gain, typically in dB.
grid_type (str, optional) – Scan grid type, only square grid is supported.
manufacturer (str, optional) – Manufacturer of software used to collect patterns.
microscope (str, optional) – Microscope used to collect patterns.
magnification (int, optional) – Microscope magnification at which patterns were collected.
sample_tilt (float, optional) – Sample tilt angle from horizontal in degrees.
scan_time (float, optional) – Scan time in s.
static_background (numpy.ndarray, optional) – Static background pattern.
version (str, optional) – Version of software used to collect patterns.
working_distance (float, optional) – Working distance in mm.
xpc (float, optional) – Pattern centre horizontal coordinate with respect to detector centre, as viewed from the detector to the sample.
ypc (float, optional) – Pattern centre vertical coordinate with respect to detector centre, as viewed from the detector to the sample.
zpc (float, optional) – Specimen to scintillator distance.
Examples
>>> import kikuchipy as kp >>> ebsd_node = metadata_nodes("ebsd") >>> s.metadata.get_item(ebsd_node + '.xpc') 1.0 >>> s.set_experimental_parameters(xpc=0.50726) >>> s.metadata.get_item(ebsd_node + '.xpc') 0.50726
-
set_phase_parameters
(number=1, atom_coordinates=None, formula=None, info=None, lattice_constants=None, laue_group=None, material_name=None, point_group=None, setting=None, source=None, space_group=None, symmetry=None)[source]¶ Set parameters for one phase in signal metadata.
A phase node with default values is created if none is present in the metadata when this method is called.
- Parameters
number (int, optional) – Phase number.
atom_coordinates (dict, optional) – Dictionary of dictionaries with one or more of the atoms in the unit cell, on the form {‘1’: {‘atom’: ‘Ni’, ‘coordinates’: [0, 0, 0], ‘site_occupation’: 1, ‘debye_waller_factor’: 0}, ‘2’: {‘atom’: ‘O’,… etc. debye_waller_factor in units of nm^2, and site_occupation in range [0, 1].
formula (str, optional) – Phase formula, e.g. ‘Fe2’ or ‘Ni’.
info (str, optional) – Whatever phase info the user finds relevant.
lattice_constants (numpy.ndarray or list of floats, optional) – Six lattice constants a, b, c, alpha, beta, gamma.
laue_group (str, optional) – Phase Laue group.
material_name (str, optional) – Name of material.
point_group (str, optional) – Phase point group.
setting (int, optional) – Space group’s origin setting.
source (str, optional) – Literature reference for phase data.
space_group (int, optional) – Number between 1 and 230.
symmetry (int, optional) – Phase symmetry.
Examples
>>> s.metadata.Sample.Phases.Number_1.atom_coordinates.Number_1 ├── atom = ├── coordinates = array([0., 0., 0.]) ├── debye_waller_factor = 0.0 └── site_occupation = 0.0 >>> s.set_phase_parameters( ... number=1, atom_coordinates={ ... '1': {'atom': 'Ni', 'coordinates': [0, 0, 0], ... 'site_occupation': 1, ... 'debye_waller_factor': 0.0035}}) >>> s.metadata.Sample.Phases.Number_1.atom_coordinates.Number_1 ├── atom = Ni ├── coordinates = array([0., 0., 0.]) ├── debye_waller_factor = 0.0035 └── site_occupation = 1
-
These methods are exclusive to LazyEBSD objects.
|
Write the model signal generated from the selected number of principal components directly to an .hspy file. |
-
class
kikuchipy.signals.
LazyEBSD
(*args, **kwargs)[source]¶ Bases:
kikuchipy.signals.ebsd.EBSD
,hyperspy._signals.signal2d.LazySignal2D
Lazy implementation of the
EBSD
class.This class extends HyperSpy’s LazySignal2D class for EBSD patterns.
Methods inherited from HyperSpy can be found in the HyperSpy user guide.
See docstring of
EBSD
for attributes and methods.-
get_decomposition_model_write
(components=None, dtype_learn=<class 'numpy.float32'>, mbytes_chunk=100, dir_out=None, fname_out=None)[source]¶ Write the model signal generated from the selected number of principal components directly to an .hspy file.
The model signal intensities are rescaled to the original signals’ data type range, keeping relative intensities.
- Parameters
components (
Union
[None
,int
,List
[int
]]) – If None (default), rebuilds the signal from all components. If int, rebuilds signal from components in range 0-given int. If list of ints, rebuilds signal from only components in given list.dtype_learn (
dtype
) – Data type to set learning results to (default isnumpy.float32
) before multiplication.mbytes_chunk (
int
) – Size of learning results chunks in MB, default is 100 MB as suggested in the Dask documentation.dir_out (
Optional
[str
]) – Directory to place output signal in.
Notes
Multiplying the learning results’ factors and loadings in memory to create the model signal cannot sometimes be done due to too large matrices. Here, instead, learning results are written to file, read into dask arrays and multiplied using
dask.array.matmul()
, out of core.
-
EBSDMasterPattern¶
All methods listed here are also available to
LazyEBSDMasterPattern
objects.
|
Normalize image intensities in inplace to a mean of zero with a given standard deviation. |
|
Rescale image intensities inplace. |
|
Set simulated parameters in signal metadata. |
|
Set parameters for one phase in signal metadata. |
-
class
kikuchipy.signals.
EBSDMasterPattern
(*args, **kwargs)[source]¶ Bases:
kikuchipy.signals._common_image.CommonImage
,hyperspy._signals.signal2d.Signal2D
Simulated Electron Backscatter Diffraction (EBSD) master pattern.
This class extends HyperSpy’s Signal2D class for EBSD master patterns.
Methods inherited from HyperSpy can be found in the HyperSpy user guide.
See the docstring of
hyperspy.signal.BaseSignal
for a list of attributes.-
normalize_intensity
(num_std=1, divide_by_square_root=False, dtype_out=None)[source]¶ Normalize image intensities in inplace to a mean of zero with a given standard deviation.
- Parameters
num_std (
int
) – Number of standard deviations of the output intensities. Default is 1.divide_by_square_root (
bool
) – Whether to divide output intensities by the square root of the signal dimension size. Default is False.dtype_out (
Optional
[dtype
]) – Data type of normalized images. If None (default), the input images’ data type is used.
Notes
Data type should always be changed to floating point, e.g.
np.float32
withchange_dtype()
, before normalizing the intensities.Examples
>>> np.mean(s.data) 146.0670987654321 >>> s.change_dtype(np.float32) # Or passing dtype_out=np.float32 >>> s.normalize_intensity() >>> np.mean(s.data) 2.6373216e-08
Notes
Rescaling RGB images is not possible. Use RGB channel normalization when creating the image instead.
-
rescale_intensity
(relative=False, in_range=None, out_range=None, dtype_out=None, percentiles=None)[source]¶ Rescale image intensities inplace.
Output min./max. intensity is determined from out_range or the data type range of the
numpy.dtype
passed to dtype_out if out_range is None.This method is based on
skimage.exposure.rescale_intensity()
.- Parameters
relative (
bool
) – Whether to keep relative intensities between images (default is False). If True, in_range must be None, because in_range is in this case set to the global min./max. intensity.in_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of input images. If None (default), stretching is performed when in_range is set to a narrower in_range is set to pattern min./max intensity. Contrast intensity range than the input patterns. Must be None if relative is True or percentiles are passed.out_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of output images. If None (default), out_range is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of rescaled images, default is input images’ data type.percentiles (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Disregard intensities outside these percentiles. Calculated per image. Must be None if in_range or relative is passed. Default is None.
See also
kikuchipy.signals.CommonImage.normalize_intensity()
,kikuchipy.pattern.rescale_intensity()
,skimage.exposure.rescale_intensity()
Examples
Image intensities are stretched to fill the available grey levels in the input images’ data type range or any
numpy.dtype
range passed to dtype_out, either keeping relative intensities between images or not:>>> print(s.data.dtype_out, s.data.min(), s.data.max(), ... s.inav[0, 0].data.min(), s.inav[0, 0].data.max()) uint8 20 254 24 233 >>> s2 = s.deepcopy() >>> s.rescale_intensity(dtype_out=np.uint16) >>> print(s.data.dtype_out, s.data.min(), s.data.max(), ... s.inav[0, 0].data.min(), s.inav[0, 0].data.max()) uint16 0 65535 0 65535 >>> s2.rescale_intensity(relative=True) >>> print(s2.data.dtype_out, s2.data.min(), s2.data.max(), ... s2.inav[0, 0].data.min(), s2.inav[0, 0].data.max()) uint8 0 255 4 232
Contrast stretching can be performed by passing percentiles:
>>> s.rescale_intensity(percentiles=(1, 99))
Here, the darkest and brightest pixels within the 1% percentile are set to the ends of the data type range, e.g. 0 and 255 respectively for images of
uint8
data type.Notes
Rescaling RGB images is not possible. Use RGB channel normalization when creating the image instead.
-
set_phase_parameters
(number=1, atom_coordinates=None, formula=None, info=None, lattice_constants=None, laue_group=None, material_name=None, point_group=None, setting=None, source=None, space_group=None, symmetry=None)[source]¶ Set parameters for one phase in signal metadata.
A phase node with default values is created if none is present in the metadata when this method is called.
- Parameters
number (
int
) – Phase number.atom_coordinates (
Optional
[dict
]) – Dictionary of dictionaries with one or more of the atoms in the unit cell, on the form {‘1’: {‘atom’: ‘Ni’, ‘coordinates’: [0, 0, 0], ‘site_occupation’: 1, ‘debye_waller_factor’: 0}, ‘2’: {‘atom’: ‘O’,… etc. debye_waller_factor in units of nm^2, and site_occupation in range [0, 1].formula (
Optional
[str
]) – Phase formula, e.g. ‘Fe2’ or ‘Ni’.info (
Optional
[str
]) – Whatever phase info the user finds relevant.lattice_constants (
Union
[None
,ndarray
,List
[float
],List
[int
]]) – Six lattice constants a, b, c, alpha, beta, gamma.source (
Optional
[str
]) – Literature reference for phase data.
See also
Examples
>>> s.metadata.Sample.Phases.Number_1.atom_coordinates.Number_1 ├── atom = ├── coordinates = array([0., 0., 0.]) ├── debye_waller_factor = 0.0 └── site_occupation = 0.0 >>> s.set_phase_parameters( ... number=1, atom_coordinates={ ... '1': {'atom': 'Ni', 'coordinates': [0, 0, 0], ... 'site_occupation': 1, ... 'debye_waller_factor': 0.0035}}) >>> s.metadata.Sample.Phases.Number_1.atom_coordinates.Number_1 ├── atom = Ni ├── coordinates = array([0., 0., 0.]) ├── debye_waller_factor = 0.0035 └── site_occupation = 1
-
set_simulation_parameters
(complete_cutoff=None, depth_step=None, energy_step=None, hemisphere=None, incident_beam_energy=None, max_depth=None, min_beam_energy=None, mode=None, number_of_electrons=None, pixels_along_x=None, projection=None, sample_tilt=None, smallest_interplanar_spacing=None, strong_beam_cutoff=None, weak_beam_cutoff=None)[source]¶ Set simulated parameters in signal metadata.
- Parameters
complete_cutoff (
Union
[None
,int
,float
]) – Bethe parameter c3.depth_step (
Union
[None
,int
,float
]) – Material penetration depth step size, in nm.energy_step (
Union
[None
,int
,float
]) – Energy bin size, in keV.hemisphere (
Optional
[str
]) – Which hemisphere(s) the data contains.incident_beam_energy (
Union
[None
,int
,float
]) – Incident beam energy, in keV.max_depth (
Union
[None
,int
,float
]) – Maximum material penetration depth, in nm.min_beam_energy (
Union
[None
,int
,float
]) – Minimum electron energy to consider, in keV.mode (
Optional
[str
]) – Simulation mode, e.g. Continuous slowing down approximation (CSDA) used by EMsoft.number_of_electrons (
Optional
[int
]) – Total number of incident electrons.pixels_along_x (
Optional
[int
]) – Pixels along horizontal direction.projection (
Optional
[str
]) – Which projection the pattern is in.sample_tilt (
Union
[None
,int
,float
]) – Sample tilte angle from horizontal, in degrees.smallest_interplanar_spacing (
Union
[None
,int
,float
]) – Smallest interplanar spacing, d-spacing, taken into account in the computation of the electrostatic lattice potential, in nm.strong_beam_cutoff (
Union
[None
,int
,float
]) – Bethe parameter c1.weak_beam_cutoff (
Union
[None
,int
,float
]) – Bethe parameter c2.
See also
Examples
>>> import kikuchipy as kp >>> ebsd_mp_node = kp.signals.util.metadata_nodes( ... "ebsd_master_pattern") >>> s.metadata.get_item(ebsd_mp_node + '.incident_beam_energy') 15.0 >>> s.set_simulated_parameters(incident_beam_energy=20.5) >>> s.metadata.get_item(ebsd_mp_node + '.incident_beam_energy') 20.5
-
There are no methods exclusive to LazyEBSDMasterPattern objects.
-
class
kikuchipy.signals.
LazyEBSDMasterPattern
(*args, **kwargs)[source]¶ Bases:
kikuchipy.signals.ebsd_master_pattern.EBSDMasterPattern
,hyperspy._signals.signal2d.LazySignal2D
Lazy implementation of the
EBSDMasterPattern
class.This class extends HyperSpy’s LazySignal2D class for EBSD master patterns.
Methods inherited from HyperSpy can be found in the HyperSpy user guide.
See docstring of
EBSDMasterPattern
for attributes and methods.
VirtualBSEImage¶
|
Normalize image intensities in inplace to a mean of zero with a given standard deviation. |
|
Rescale image intensities inplace. |
-
class
kikuchipy.signals.
VirtualBSEImage
(*args, **kwargs)[source]¶ Virtual backscatter electron (BSE) image(s).
This class extends HyperSpy’s Signal2D class for virtual BSE images.
Methods inherited from HyperSpy can be found in the HyperSpy user guide.
See the docstring of
hyperspy.signal.BaseSignal
for a list of attributes.-
normalize_intensity
(num_std=1, divide_by_square_root=False, dtype_out=None)[source]¶ Normalize image intensities in inplace to a mean of zero with a given standard deviation.
- Parameters
num_std (
int
) – Number of standard deviations of the output intensities. Default is 1.divide_by_square_root (
bool
) – Whether to divide output intensities by the square root of the signal dimension size. Default is False.dtype_out (
Optional
[dtype
]) – Data type of normalized images. If None (default), the input images’ data type is used.
Notes
Data type should always be changed to floating point, e.g.
np.float32
withchange_dtype()
, before normalizing the intensities.Examples
>>> np.mean(s.data) 146.0670987654321 >>> s.change_dtype(np.float32) # Or passing dtype_out=np.float32 >>> s.normalize_intensity() >>> np.mean(s.data) 2.6373216e-08
Notes
Rescaling RGB images is not possible. Use RGB channel normalization when creating the image instead.
-
rescale_intensity
(relative=False, in_range=None, out_range=None, dtype_out=None, percentiles=None)[source]¶ Rescale image intensities inplace.
Output min./max. intensity is determined from out_range or the data type range of the
numpy.dtype
passed to dtype_out if out_range is None.This method is based on
skimage.exposure.rescale_intensity()
.- Parameters
relative (
bool
) – Whether to keep relative intensities between images (default is False). If True, in_range must be None, because in_range is in this case set to the global min./max. intensity.in_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of input images. If None (default), stretching is performed when in_range is set to a narrower in_range is set to pattern min./max intensity. Contrast intensity range than the input patterns. Must be None if relative is True or percentiles are passed.out_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of output images. If None (default), out_range is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of rescaled images, default is input images’ data type.percentiles (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Disregard intensities outside these percentiles. Calculated per image. Must be None if in_range or relative is passed. Default is None.
See also
kikuchipy.signals.CommonImage.normalize_intensity()
,kikuchipy.pattern.rescale_intensity()
,skimage.exposure.rescale_intensity()
Examples
Image intensities are stretched to fill the available grey levels in the input images’ data type range or any
numpy.dtype
range passed to dtype_out, either keeping relative intensities between images or not:>>> print(s.data.dtype_out, s.data.min(), s.data.max(), ... s.inav[0, 0].data.min(), s.inav[0, 0].data.max()) uint8 20 254 24 233 >>> s2 = s.deepcopy() >>> s.rescale_intensity(dtype_out=np.uint16) >>> print(s.data.dtype_out, s.data.min(), s.data.max(), ... s.inav[0, 0].data.min(), s.inav[0, 0].data.max()) uint16 0 65535 0 65535 >>> s2.rescale_intensity(relative=True) >>> print(s2.data.dtype_out, s2.data.min(), s2.data.max(), ... s2.inav[0, 0].data.min(), s2.inav[0, 0].data.max()) uint8 0 255 4 232
Contrast stretching can be performed by passing percentiles:
>>> s.rescale_intensity(percentiles=(1, 99))
Here, the darkest and brightest pixels within the 1% percentile are set to the ends of the data type range, e.g. 0 and 255 respectively for images of
uint8
data type.Notes
Rescaling RGB images is not possible. Use RGB channel normalization when creating the image instead.
-
Utilities¶
-
kikuchipy.signals.util.
ebsd_master_pattern_metadata
()[source]¶ Return a dictionary in HyperSpy’s DictionaryTreeBrowser format with the default kikuchipy EBSD master pattern metadata.
The parameters are chosen based on the contents in EMsoft’s EBSD master pattern HDF5 file.
See
set_simulation_parameters()
for an explanation of the parameters.- Returns
md
- Return type
-
kikuchipy.signals.util.
ebsd_metadata
()[source]¶ Return a dictionary in HyperSpy’s DictionaryTreeBrowser format with the default kikuchipy EBSD metadata.
See
set_experimental_parameters()
for an explanation of the parameters.- Returns
md
- Return type
Generators¶
Generators producing signals. This workflow with generators is adopted from pyxem.
VirtualBSEGenerator¶
|
Return an in-memory signal with a stack of virtual backscatter electron (BSE) images by integrating the intensities within regions of interest (ROI) defined by the detector grid_shape. |
|
Return an in-memory RGB virtual BSE image from three regions of interest (ROIs) on the EBSD detector, with a potential “alpha channel” in which all three arrays are multiplied by a fourth. |
|
Plot a pattern with the detector grid superimposed, potentially coloring the edges of three grid tiles red, green and blue. |
|
Return a rectangular region of interest (ROI) on the EBSD detector from one or multiple generator grid tile indices as row(s) and column(s). |
-
class
kikuchipy.generators.
VirtualBSEGenerator
(signal)[source]¶ Generates virtual backscatter electron (BSE) images for a specified electron backscatter diffraction (EBSD) signal and a set of EBSD detector areas.
-
signal
¶
-
get_images_from_grid
(dtype_out=<class 'numpy.float32'>)[source]¶ Return an in-memory signal with a stack of virtual backscatter electron (BSE) images by integrating the intensities within regions of interest (ROI) defined by the detector grid_shape.
- Parameters
dtype_out (
dtype
) – Output data type, default is float32.- Returns
vbse_images – In-memory signal with virtual BSE images.
- Return type
Examples
>>> s <EBSD, title: Pattern, dimensions: (200, 149|60, 60)> >>> vbse_gen = VirtualBSEGenerator(s) >>> vbse_gen.grid_shape = (5, 5) >>> vbse = vbse_gen.get_images_from_grid() >>> vbse <VirtualBSEImage, title: , dimensions: (5, 5|200, 149)>
-
get_rgb_image
(r, g, b, percentiles=None, normalize=True, alpha=None, dtype_out=<class 'numpy.uint8'>, **kwargs)[source]¶ Return an in-memory RGB virtual BSE image from three regions of interest (ROIs) on the EBSD detector, with a potential “alpha channel” in which all three arrays are multiplied by a fourth.
- Parameters
r (
Union
[BaseInteractiveROI
,Tuple
,List
[BaseInteractiveROI
],List
[Tuple
]]) – One ROI or a list of ROIs, or one tuple or a list of tuples with detector grid indices specifying one or more ROI(s). Intensities within the specified ROI(s) are summed up to form the red color channel.g (
Union
[BaseInteractiveROI
,Tuple
,List
[BaseInteractiveROI
],List
[Tuple
]]) – One ROI or a list of ROIs, or one tuple or a list of tuples with detector grid indices specifying one or more ROI(s). Intensities within the specified ROI(s) are summed up to form the green color channel.b (
Union
[BaseInteractiveROI
,Tuple
,List
[BaseInteractiveROI
],List
[Tuple
]]) – One ROI or a list of ROIs, or one tuple or a list of tuples with detector grid indices specifying one or more ROI(s). Intensities within the specified ROI(s) are summed up to form the blue color channel.normalize (
bool
) – Whether to normalize the individual images (channels) before RGB image creation.alpha (
Union
[None
,ndarray
,VirtualBSEImage
]) – “Alpha channel”. If None (default), no “alpha channel” is added to the image.percentiles (
Optional
[Tuple
]) – Whether to apply contrast stretching with a given percentile tuple with percentages, e.g. (0.5, 99.5), after creating the RGB image. If None (default), no contrast stretching is performed.dtype_out (
Union
[uint8
,uint16
]) – Output data type, either np.uint16 or np.uint8 (default).kwargs – Keyword arguments passed to :func:` ~kikuchipy.generators.util.virtual_bse.get_rgb_image`.
- Returns
vbse_rgb_image – Virtual RGB image in memory.
- Return type
See also
kikuchipy.signals.EBSD.plot_virtual_bse_intensity()
,kikuchipy.signals.EBSD.get_virtual_bse_intensity()
,kikuchipy.generators.util.get_rgb_image()
Notes
HyperSpy only allows for RGB signal dimensions with data types unsigned 8 or 16 bit.
-
plot_grid
(pattern_idx=None, rgb_channels=None, visible_indices=True, **kwargs)[source]¶ Plot a pattern with the detector grid superimposed, potentially coloring the edges of three grid tiles red, green and blue.
- Parameters
pattern_idx (
Optional
[Tuple
[int
, …]]) – A tuple of integers defining the pattern to superimpose the grid on. If None (default), the first pattern is used.rgb_channels (
Union
[None
,List
[Tuple
],List
[List
[Tuple
]]]) – A list of tuple indices defining three or more detector grid tiles which edges to color red, green and blue. If None (default), no tiles’ edges are colored.visible_indices (
bool
) – Whether to show grid indices. Default is True.kwargs – Keyword arguments passed to
matplotlib.pyplot.axhline()
and axvline, used by HyperSpy to draw lines.
- Returns
pattern – A single pattern with the markers added.
- Return type
-
roi_from_grid
(index)[source]¶ Return a rectangular region of interest (ROI) on the EBSD detector from one or multiple generator grid tile indices as row(s) and column(s).
-
Utilities¶
-
kikuchipy.generators.util.
get_rgb_image
(channels, percentiles=None, normalize=True, alpha=None, dtype_out=<class 'numpy.uint8'>, **kwargs)[source]¶ Return an RGB image from three numpy arrays, with a potential alpha channel.
- Parameters
channels (
List
[ndarray
]) – A list of np.ndarray for the red, green and blue channel, respectively.normalize (
bool
) – Whether to normalize the individual channels before RGB image creation.alpha (
Optional
[ndarray
]) – Potential alpha channel. If None (default), no alpha channel is added to the image.percentiles (
Optional
[Tuple
]) – Whether to apply contrast stretching with a given percentile tuple with percentages, e.g. (0.5, 99.5), after creating the RGB image. If None (default), no contrast stretching is performed.dtype_out (
Union
[uint8
,uint16
]) – Output data type, either np.uint16 or np.uint8 (default).kwargs – Keyword arguments passed to :func:` ~kikuchipy.generators.util.virtual_bse.normalize_image`.
- Returns
rgb_image – RGB image.
- Return type
np.ndarray
-
kikuchipy.generators.util.
normalize_image
(image, add_bright=0, contrast=1.0, dtype_out=<class 'numpy.uint8'>)[source]¶ Normalize an image’s intensities to a mean of 0 and a standard deviation of 1, with the possibility to also scale by a contrast factor and shift the brightness values.
Clips intensities to uint8 data type range, [0, 255].
Adapted from the aloe/xcdskd package.
Pattern¶
Single pattern processing¶
This module mainly includes functions operating on single EBSD patterns as
numpy.ndarray
.
|
Compute the discrete Fast Fourier Transform (FFT) of an EBSD pattern. |
|
Filter an EBSD patterns in the frequency domain. |
|
Get the frequency vectors in a Fourier Transform spectrum. |
|
Compute the FFT spectrum of a Fourier transformed EBSD pattern. |
|
Get the dynamic background in an EBSD pattern. |
|
Return the image quality of an EBSD pattern. |
|
Compute the inverse Fast Fourier Transform (IFFT) of an FFT of an EBSD pattern. |
|
Normalize image intensities to a mean of zero and a given standard deviation. |
|
Remove the dynamic background in an EBSD pattern. |
|
Rescale intensities in an EBSD pattern. |
-
kikuchipy.pattern.
fft
(pattern, apodization_window=None, shift=False, real_fft_only=False, **kwargs)[source]¶ Compute the discrete Fast Fourier Transform (FFT) of an EBSD pattern.
Very light wrapper around routines in
scipy.fft
. The routines are wrapped instead of used directly to accommodate easy setting of shift and real_fft_only.- Parameters
pattern (
ndarray
) – EBSD pattern.apodization_window (
Union
[None
,ndarray
,Window
]) – An apodization window to apply before the FFT in order to suppress streaks.shift (
bool
) – Whether to shift the zero-frequency component to the centre of the spectrum (default is False).real_fft_only (
bool
) – If True, the discrete FFT is computed for real input usingscipy.fft.rfft2()
. If False (default), it is computed usingscipy.fft.fft2()
.kwargs – Keyword arguments pass to
scipy.fft.fft2()
orscipy.fft.rfft2()
.
- Returns
out – The result of the 2D FFT.
- Return type
-
kikuchipy.pattern.
fft_filter
(pattern, transfer_function, apodization_window=None, shift=False)[source]¶ Filter an EBSD patterns in the frequency domain.
- Parameters
pattern (
ndarray
) – EBSD pattern.transfer_function (
Union
[ndarray
,Window
]) – Filter transfer function in the frequency domain.apodization_window (
Union
[None
,ndarray
,Window
]) – An apodization window to apply before the FFT in order to suppress streaks.shift (
bool
) – Whether to shift the zero-frequency component to the centre of the spectrum. Default is False.
- Returns
filtered_pattern – Filtered EBSD pattern.
- Return type
-
kikuchipy.pattern.
fft_frequency_vectors
(shape)[source]¶ Get the frequency vectors in a Fourier Transform spectrum.
- Parameters
- Returns
frequency_vectors – Frequency vectors.
- Return type
-
kikuchipy.pattern.
fft_spectrum
(fft_pattern)[source]¶ Compute the FFT spectrum of a Fourier transformed EBSD pattern.
- Parameters
fft_pattern (
ndarray
) – Fourier transformed EBSD pattern.- Returns
fft_spectrum – 2D FFT spectrum of the EBSD pattern.
- Return type
-
kikuchipy.pattern.
get_dynamic_background
(pattern, filter_domain='frequency', std=None, truncate=4.0)[source]¶ Get the dynamic background in an EBSD pattern.
The background is obtained either in the frequency domain, by a low pass Fast Fourier Transform (FFT) Gaussian filter, or in the spatial domain by a Gaussian filter.
Data type is preserved.
- Parameters
pattern (
ndarray
) – EBSD pattern.filter_domain (
str
) – Whether to obtain the dynamic background by applying a Gaussian convolution filter in the “frequency” (default) or “spatial” domain.std (
Union
[None
,int
,float
]) – Standard deviation of the Gaussian window. If None (default), a deviation of pattern width/8 is chosen.truncate (
Union
[int
,float
]) – Truncate the Gaussian window at this many standard deviations. Default is 4.0.
- Returns
dynamic_bg – The dynamic background.
- Return type
-
kikuchipy.pattern.
get_image_quality
(pattern, normalize=True, frequency_vectors=None, inertia_max=None)[source]¶ Return the image quality of an EBSD pattern.
The image quality is calculated based on the procedure defined by Krieger Lassen [Lassen1994].
- Parameters
pattern (
ndarray
) – EBSD pattern.normalize (
bool
) – Whether to normalize the pattern to a mean of zero and standard deviation of 1 before calculating the image quality (default is True).frequency_vectors (
Optional
[ndarray
]) – Integer 2D array assigning each FFT spectrum frequency component a weight. If None (default), these are calculated fromfft_frequency_vectors()
. This only depends on the pattern shape.inertia_max (
Union
[None
,int
,float
]) – Maximum possible inertia of the FFT power spectrum of the image. If None (default), this is calculated from the frequency_vectors, which in this case must be passed. This only depends on the pattern shape.
- Returns
image_quality – Image quality of the pattern.
- Return type
-
kikuchipy.pattern.
ifft
(fft_pattern, shift=False, real_fft_only=False, **kwargs)[source]¶ Compute the inverse Fast Fourier Transform (IFFT) of an FFT of an EBSD pattern.
Very light wrapper around routines in
scipy.fft
. The routines are wrapped instead of used directly to accommodate easy setting of shift and real_fft_only.- Parameters
fft_pattern (
ndarray
) – FFT of EBSD pattern.shift (
bool
) – Whether to shift the zero-frequency component back to the corners of the spectrum (default is False).real_fft_only (
bool
) – If True, the discrete IFFT is computed for real input usingscipy.fft.irfft2()
. If False (default), it is computed usingscipy.fft.ifft2()
.kwargs – Keyword arguments pass to
scipy.fft.ifft()
.
- Returns
pattern – Real part of the IFFT of the EBSD pattern.
- Return type
-
kikuchipy.pattern.
normalize_intensity
(pattern, num_std=1, divide_by_square_root=False)[source]¶ Normalize image intensities to a mean of zero and a given standard deviation.
Data type is preserved.
- Parameters
- Returns
normalized_pattern – Normalized pattern.
- Return type
Notes
Data type should always be changed to floating point, e.g.
np.float32
withnumpy.ndarray.astype()
, before normalizing the intensities.
-
kikuchipy.pattern.
remove_dynamic_background
(pattern, operation='subtract', filter_domain='frequency', std=None, truncate=4.0, dtype_out=None)[source]¶ Remove the dynamic background in an EBSD pattern.
The removal is performed by subtracting or dividing by a Gaussian blurred version of the pattern. The blurred version is obtained either in the frequency domain, by a low pass Fast Fourier Transform (FFT) Gaussian filter, or in the spatial domain by a Gaussian filter. Returned pattern intensities are rescaled to fill the input data type range.
- Parameters
pattern (
ndarray
) – EBSD pattern.operation (
str
) – Whether to “subtract” (default) or “divide” by the dynamic background pattern.filter_domain (
str
) – Whether to obtain the dynamic background by applying a Gaussian convolution filter in the “frequency” (default) or “spatial” domain.std (
Union
[None
,int
,float
]) – Standard deviation of the Gaussian window. If None (default), it is set to width/8.truncate (
Union
[int
,float
]) – Truncate the Gaussian window at this many standard deviations. Default is 4.0.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of corrected pattern. If None (default), it is set to input patterns’ data type.
- Returns
corrected_pattern – Pattern with the dynamic background removed.
- Return type
-
kikuchipy.pattern.
rescale_intensity
(pattern, in_range=None, out_range=None, dtype_out=None)[source]¶ Rescale intensities in an EBSD pattern.
Pattern max./min. intensity is determined from out_range or the data type range of
numpy.dtype
passed to dtype_out.This method is based on
skimage.exposure.rescale_intensity()
.- Parameters
pattern (
ndarray
) – EBSD pattern.in_range (
Optional
[Tuple
[Union
[int
,float
], …]]) – Min./max. intensity values of the input pattern. If None (default), it is set to the pattern’s min./max intensity.out_range (
Optional
[Tuple
[Union
[int
,float
], …]]) – Min./max. intensity values of the rescaled pattern. If None (default), it is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Optional
[dtype
]) – Data type of the rescaled pattern. If None (default), it is set to the same data type as the input pattern.
- Returns
rescaled_pattern – Rescaled pattern.
- Return type
Chunk processing¶
This module includes functions for operating on numpy.ndarray
or
dask.array.Array
chunks of EBSD patterns.
|
Local contrast enhancement of a chunk of EBSD patterns with adaptive histogram equalization. |
|
Average a chunk of patterns with its neighbours within a window. |
|
Filter a chunk of EBSD patterns in the frequency domain. |
|
Obtain the dynamic background in a chunk of EBSD patterns. |
|
Compute the image quality in a chunk of EBSD patterns. |
|
Normalize intensities in a chunk of EBSD patterns to a mean of zero with a given standard deviation. |
|
Correct the dynamic background in a chunk of EBSD patterns. |
|
Remove the static background in a chunk of EBSD patterns. |
|
Rescale pattern intensities in a chunk of EBSD patterns. |
-
kikuchipy.pattern.chunk.
adaptive_histogram_equalization
(patterns, kernel_size, clip_limit=0, nbins=128)[source]¶ Local contrast enhancement of a chunk of EBSD patterns with adaptive histogram equalization.
This method makes use of
skimage.exposure.equalize_adapthist()
.- Parameters
kernel_size (
Union
[Tuple
[int
,int
],List
[int
]]) – Shape of contextual regions for adaptive histogram equalization.clip_limit (
Union
[int
,float
]) – Clipping limit, normalized between 0 and 1 (higher values give more contrast). Default is 0.nbins (
int
) – Number of gray bins for histogram. Default is 128.
- Returns
equalized_patterns – Patterns with enhanced contrast.
- Return type
-
kikuchipy.pattern.chunk.
average_neighbour_patterns
(patterns, window_sums, window, dtype_out=None)[source]¶ Average a chunk of patterns with its neighbours within a window.
- Parameters
patterns (
ndarray
) – Patterns to average, with some overlap with surrounding chunks.window_sums (
ndarray
) – Sum of window data for each image.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of averaged patterns. If None (default), it is set to the same data type as the input patterns.
- Returns
averaged_patterns – Averaged patterns.
- Return type
-
kikuchipy.pattern.chunk.
fft_filter
(patterns, filter_func, transfer_function, dtype_out=None, **kwargs)[source]¶ Filter a chunk of EBSD patterns in the frequency domain.
Patterns are transformed via the Fast Fourier Transform (FFT) to the frequency domain, where their spectrum is multiplied by a filter transfer_function, and the filtered spectrum is subsequently transformed to the spatial domain via the inverse FFT (IFFT).
Filtered patterns are rescaled to the data type range of dtype_out.
- Parameters
patterns (
ndarray
) – EBSD patterns.filter_func (
Union
[fft_filter
,_fft_filter
]) – Function to apply transfer_function with.transfer_function (
Union
[ndarray
,Window
]) – Filter transfer function in the frequency domain.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of output patterns. If None (default), it is set to the input patterns’ data type.kwargs – Keyword arguments passed to the filter_func.
- Returns
filtered_patterns – Filtered EBSD patterns.
- Return type
-
kikuchipy.pattern.chunk.
get_dynamic_background
(patterns, filter_func, dtype_out=None, **kwargs)[source]¶ Obtain the dynamic background in a chunk of EBSD patterns.
- Parameters
filter_func (
Union
[gaussian_filter
,fft_filter
]) – Function where a Gaussian convolution filter is applied, in the frequency or spatial domain. Eitherscipy.ndimage.gaussian_filter()
orkikuchipy.util.barnes_fftfilter.fft_filter()
.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of background patterns. If None (default), it is set to input patterns’ data type.kwargs – Keyword arguments passed to the Gaussian blurring function passed to filter_func.
- Returns
background – Large scale variations in the input EBSD patterns.
- Return type
-
kikuchipy.pattern.chunk.
get_image_quality
(patterns, frequency_vectors=None, inertia_max=None, normalize=True)[source]¶ Compute the image quality in a chunk of EBSD patterns.
The image quality is calculated based on the procedure defined by Krieger Lassen [Lassen1994].
- Parameters
frequency_vectors (
Optional
[ndarray
]) – Integer 2D array with values corresponding to the weight given each FFT spectrum frequency component. If None (default), these are calculated fromfft_frequency_vectors()
.inertia_max (
Union
[None
,int
,float
]) – Maximum inertia of the FFT power spectrum of the image. If None (default), this is calculated from the frequency_vectors.normalize (
bool
) – Whether to normalize patterns to a mean of zero and standard deviation of 1 before calculating the image quality. Default is True.
- Returns
image_quality_chunk – Image quality of patterns.
- Return type
-
kikuchipy.pattern.chunk.
normalize_intensity
(patterns, num_std=1, divide_by_square_root=False, dtype_out=None)[source]¶ Normalize intensities in a chunk of EBSD patterns to a mean of zero with a given standard deviation.
- Parameters
patterns (
Union
[ndarray
,Array
]) – Patterns to normalize the intensity in.num_std (
int
) – Number of standard deviations of the output intensities. Default is 1.divide_by_square_root (
bool
) – Whether to divide output intensities by the square root of the pattern size. Default is False.dtype_out (
Optional
[dtype
]) – Data type of normalized patterns. If None (default), the input patterns’ data type is used.
- Returns
normalized_patterns – Normalized patterns.
- Return type
Notes
Data type should always be changed to floating point, e.g.
np.float32
withnumpy.ndarray.astype()
, before normalizing the intensities.
-
kikuchipy.pattern.chunk.
remove_dynamic_background
(patterns, filter_func, operation_func, out_range=None, dtype_out=None, **kwargs)[source]¶ Correct the dynamic background in a chunk of EBSD patterns.
The correction is performed by subtracting or dividing by a Gaussian blurred version of each pattern. Returned pattern intensities are rescaled to fill the input data type range.
- Parameters
filter_func (
Union
[gaussian_filter
,fft_filter
]) – Function where a Gaussian convolution filter is applied, in the frequency or spatial domain. Eitherscipy.ndimage.gaussian_filter()
orkikuchipy.util.barnes_fftfilter.fft_filter()
.operation_func (
Union
[<ufunc ‘subtract’>, <ufunc ‘true_divide’>]) – Function to subtract or divide by the dynamic background pattern.out_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity values of the output patterns. If None (default), out_range is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of corrected patterns. If None (default), it is set to input patterns’ data type.kwargs – Keyword arguments passed to the Gaussian blurring function passed to filter_func.
- Returns
corrected_patterns – Dynamic background corrected patterns.
- Return type
See also
kikuchipy.signals.ebsd.EBSD.remove_dynamic_background()
,kikuchipy.util.pattern.remove_dynamic_background()
-
kikuchipy.pattern.chunk.
remove_static_background
(patterns, static_bg, operation_func, scale_bg=False, in_range=None, out_range=None, dtype_out=None)[source]¶ Remove the static background in a chunk of EBSD patterns.
Removal is performed by subtracting or dividing by a static background pattern. Resulting pattern intensities are rescaled keeping relative intensities or not and stretched to fill the available grey levels in the patterns’ data type range.
- Parameters
static_bg (
Union
[ndarray
,Array
]) – Static background pattern. If None is passed (default) we try to read it from the signal metadata.operation_func (
Union
[<ufunc ‘subtract’>, <ufunc ‘true_divide’>]) – Function to subtract or divide by the dynamic background pattern.scale_bg (
bool
) – Whether to scale the static background pattern to each individual pattern’s data range before removal (default is False).in_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity values of input and output patterns. If None (default), it is set to the overall pattern min./max, losing relative intensities between patterns.out_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity values of the output patterns. If None (default), out_range is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of corrected patterns. If None (default), it is set to input patterns’ data type.
- Returns
corrected_patterns – Patterns with the static background removed.
- Return type
-
kikuchipy.pattern.chunk.
rescale_intensity
(patterns, in_range=None, out_range=None, dtype_out=None, percentiles=None)[source]¶ Rescale pattern intensities in a chunk of EBSD patterns.
Chunk max./min. intensity is determined from out_range or the data type range of
numpy.dtype
passed to dtype_out.- Parameters
in_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of input patterns. If None (default), in_range is set to pattern min./max. Contrast stretching is performed when in_range is set to a narrower intensity range than the input patterns.out_range (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Min./max. intensity of output patterns. If None (default), out_range is set to dtype_out min./max according to skimage.util.dtype.dtype_range.dtype_out (
Union
[None
,dtype
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Data type of rescaled patterns. If None (default), it is set to the same data type as the input patterns.percentiles (
Union
[None
,Tuple
[int
,int
],Tuple
[float
,float
]]) – Disregard intensities outside these percentiles. Calculated per pattern. Must be None if in_range is passed (default is None).
- Returns
rescaled_patterns – Rescaled patterns.
- Return type
Correlate¶
Utilities for computing similarities between EBSD patterns.
-
kikuchipy.pattern.correlate.
normalized_correlation_coefficient
(pattern, template, zero_normalised=True)[source]¶ Calculate the normalized or zero-normalized correlation coefficient between a pattern and a template following [Gonzalez2017].
- Parameters
- Returns
coefficient – Correlation coefficient in range [-1, 1] if zero normalised, otherwise [0, 1].
- Return type
References
- Gonzalez2017
Gonzalez, R. E. Woods, “Digital Image Processing,” 4th edition, Pearson Education Limited, 2017.
Filters¶
-
kikuchipy.filters.
distance_to_origin
(shape, origin=None)[source]¶ Return the distance to the window origin in pixels.
-
kikuchipy.filters.
highpass_fft_filter
(shape, cutoff, cutoff_width=None)[source]¶ Return a frequency domain high-pass filter transfer function in 2D.
Used in [Wilkinson2006].
- Parameters
- Returns
w – 2D transfer function.
- Return type
Notes
The high-pass filter transfer function is defined as
\[\begin{split}w(r) = e^{-\left(\frac{c - r}{\sqrt{2}w_c/2}\right)^2}, w(r) = \begin{cases} 0, & r < c - 2w_c\\ 1, & r > c, \end{cases}\end{split}\]where \(r\) is the radial distance to the window centre, \(c\) is the cut-off frequency, and \(w_c\) is the width of the cut-off region.
Examples
>>> import kikuchipy as kp >>> w1 = kp.filters.Window( ... "highpass", cutoff=1, cutoff_width=0.5, shape=(96, 96)) >>> w2 = kp.filters.highpass_fft_filter( ... shape=(96, 96), cutoff=1, cutoff_width=0.5) >>> np.allclose(w1, w2) True
-
kikuchipy.filters.
lowpass_fft_filter
(shape, cutoff, cutoff_width=None)[source]¶ Return a frequency domain low-pass filter transfer function in 2D.
Used in [Wilkinson2006].
- Parameters
- Returns
w – 2D transfer function.
- Return type
Notes
The low-pass filter transfer function is defined as
\[\begin{split}w(r) = e^{-\left(\frac{r - c}{\sqrt{2}w_c/2}\right)^2}, w(r) = \begin{cases} 0, & r > c + 2w_c \\ 1, & r < c, \end{cases}\end{split}\]where \(r\) is the radial distance to the window centre, \(c\) is the cut-off frequency, and \(w_c\) is the width of the cut-off region.
Examples
>>> import kikuchipy as kp >>> w1 = kp.filters.Window( ... "lowpass", cutoff=30, cutoff_width=15, shape=(96, 96)) >>> w2 = kp.filters.lowpass_fft_filter( shape=(96, 96), cutoff=30, cutoff_width=15) >>> np.allclose(w1, w2) True
-
kikuchipy.filters.
modified_hann
(Nx)[source]¶ Return a 1D modified Hann window with the maximum value normalized to 1.
Used in [Wilkinson2006].
- Parameters
Nx (
int
) – Number of points in the window.- Returns
w – 1D Hann window.
- Return type
Notes
The modified Hann window is defined as
\[w(x) = \cos\left(\frac{\pi x}{N_x}\right),\]with \(x\) relative to the window centre.
References
- Wilkinson2006(1,2,3)
A. J. Wilkinson, G. Meaden, D. J. Dingley, “High resolution mapping of strains and rotations using electron backscatter diffraction,” Materials Science and Technology 22(11), (2006), doi: https://doi.org/10.1179/174328406X130966.
Examples
>>> import kikuchipy as kp >>> w1 = kp.filters.modified_hann(Nx=30) >>> w2 = kp.filters.Window("modified_hann", shape=(30,)) >>> np.allclose(w1, w2) True
-
class
kikuchipy.filters.
Window
[source]¶ A window/kernel/mask/filter of a given shape with some values.
This class is a subclass of
numpy.ndarray
with some additional convenience methods.It can be used to create a transfer function for filtering in the frequency domain, create an averaging window for averaging patterns with their nearest neighbours, and so on.
- Parameters
window ("circular", "rectangular", "gaussian", str, numpy.ndarray, or dask.array.Array, optional) – Window type to create. Available types are listed in
scipy.signal.windows.get_window()
and includes “rectangular” and “gaussian”, in addition to a “circular” window (default) filled with ones in which corner data are set to zero, a “modified_hann” window and “lowpass” and “highpass” FFT windows. A window element is considered to be in a corner if its radial distance to the origin (window centre) is shorter or equal to the half width of the windows’s longest axis. A 1D or 2Dnumpy.ndarray
ordask.array.Array
can also be passed.shape (sequence of int, optional) – Shape of the window. Not used if a custom window is passed to window. This can be either 1D or 2D, and can be asymmetrical. Default is (3, 3).
**kwargs – Keyword arguments passed to the window type. If none are passed, the default values of that particular window are used.
Examples
>>> import kikuchipy as kp >>> import numpy as np
The following passed parameters are the default
>>> w = kp.filters.Window(window="circular", shape=(3, 3)) >>> w Window (3, 3) circular [[0. 1. 0.] [1. 1. 1.] [0. 1. 0.]]
A window can be made circular
>>> w = kp.filters.Window(window="rectangular") >>> w Window (3, 3) rectangular [[1. 1. 1.] [1. 1. 1.] [1. 1. 1.]] >>> w.make_circular() >>> w Window (3, 3) circular [[0. 1. 0.] [1. 1. 1.] [0. 1. 0.]]
A custom window can be created
>>> w = kp.filters.Window(np.arange(6).reshape(3, 2)) >>> w Window (3, 2) custom [[0, 1] [2, 3] [4, 5]]
To create a Gaussian window with a standard deviation of 2, obtained from
scipy.signal.windows.gaussian()
>>> w = kp.filters.Window(window="gaussian", std=2) >>> w Window (3, 3) gaussian [[0.77880078, 0.8824969 , 0.77880078] [0.8824969 , 1. , 0.8824969 ] [0.77880078, 0.8824969 , 0.77880078]]
See also
scipy.signal.windows.get_window()
,kikuchipy.filters.modified_hann
,kikuchipy.filters.highpass_fft_filter
,kikuchipy.filters.lowpass_fft_filter
-
circular
: bool = False¶
-
make_circular
()[source]¶ Make window circular.
The data of window elements who’s radial distance to the window origin is shorter or equal to the half width of the window’s longest axis are set to zero. This has no effect if the window has only one axis.
-
name
: str = None¶
-
plot
(grid=True, show_values=True, cmap='viridis', textcolors=None, cmap_label='Value')[source]¶ Plot window values with indices relative to the origin.
- Parameters
grid (
bool
) – Whether to separate each value with a white spacing in a grid. Default is True.show_values (
bool
) – Whether to show values as text in centre of element. Default is True.cmap (
str
) – A color map to color data with, available inmatplotlib.colors.ListedColormap
. Default is “viridis”.textcolors (
Optional
[List
[str
]]) – A list of two color specifications. The first is used for values below a threshold, the second for those above. If None (default), this is set to [“white”, “black”].cmap_label (
str
) – Color map label. Default is “Value”.
- Return type
- Returns
fig
image
colorbar
Examples
A plot of window data with indices relative to the origin, showing element values and x/y ticks, can be produced and written to file
>>> figure, image, colorbar = w.plot( ... cmap="inferno", grid=True, show_values=True) >>> figure.savefig('my_kernel.png')
Input/output¶
-
kikuchipy.io._io.
load
(filename, lazy=False, **kwargs)[source]¶ Load an
EBSD
orEBSDMasterPattern
object from a supported file.This function is a modified version of
hyperspy.io.load()
.- Parameters
filename (
str
) – Name of file to load.lazy (
bool
) – Open the data lazily without actually reading the data from disk until required. Allows opening arbitrary sized datasets. Default is False.kwargs – Keyword arguments passed to the corresponding kikuchipy reader. See their individual documentation for available options.
These plugin functions import patterns and parameters from file formats into
EBSD
or
EBSDMasterPattern
(or
LazyEBSD
or
LazyEBSDMasterPattern
if loading lazily) objects.
h5ebsd¶
-
kikuchipy.io.plugins.h5ebsd.
brukerheader2dicts
(scan_group, md)[source]¶ Return scan metadata dictionaries from a Bruker h5ebsd file.
- Parameters
scan_group (Group) – HDF group of scan data and header.
md (hyperspy.misc.utils.DictionaryTreeBrowser) – Dictionary with empty fields from kikuchipy’s metadata.
- Return type
Tuple
[DictionaryTreeBrowser
,DictionaryTreeBrowser
,DictionaryTreeBrowser
]- Returns
md – kikuchipy
metadata
elements available in Bruker file.omd – All metadata available in Bruker file.
scan_size – Scan, image, step and detector pixel size available in Bruker file.
-
kikuchipy.io.plugins.h5ebsd.
check_h5ebsd
(file)[source]¶ Check if HDF file is an h5ebsd file by searching for datasets containing manufacturer, version and scans in the top group.
- Parameters
file (File) – File where manufacturer, version and scan datasets should reside in the top group.
-
kikuchipy.io.plugins.h5ebsd.
dict2h5ebsdgroup
(dictionary, group, **kwargs)[source]¶ Write a dictionary from
metadata
to datasets in a new group in an opened HDF file in the h5ebsd format.- Parameters
dictionary (
dict
) –Metadata
, with keys as dataset names.group (Group) – HDF group to write dictionary to.
kwargs – Keyword arguments passed to
Group.require_dataset()
.
-
kikuchipy.io.plugins.h5ebsd.
edaxheader2dicts
(scan_group, md)[source]¶ Return scan metadata dictionaries from an EDAX TSL h5ebsd file.
- Parameters
scan_group (Group) – HDF group of scan data and header.
md (
DictionaryTreeBrowser
) – Dictionary with empty fields from kikuchipy’s metadata.
- Return type
Tuple
[DictionaryTreeBrowser
,DictionaryTreeBrowser
,DictionaryTreeBrowser
]- Returns
md – kikuchipy
metadata
elements available in EDAX file.omd – All metadata available in EDAX file.
scan_size – Scan, image, step and detector pixel size available in EDAX file.
-
kikuchipy.io.plugins.h5ebsd.
file_reader
(filename, scan_group_names=None, lazy=False, **kwargs)[source]¶ Read electron backscatter diffraction patterns from an h5ebsd file [Jackson2014]. A valid h5ebsd file has at least one top group with the subgroup ‘EBSD’ with the subgroups ‘Data’ (patterns etc.) and ‘Header’ (
metadata
etc.).- Parameters
filename (
str
) – Full file path of the HDF file.scan_group_names (
Union
[None
,str
,List
[str
]]) – Name or a list of names of HDF5 top group(s) containing the scan(s) to return. If None, the first scan in the file is returned.lazy (
bool
) – Open the data lazily without actually reading the data from disk until required. Allows opening arbitrary sized datasets. Default is False.kwargs – Key word arguments passed to
File
.
- Returns
scan_dict_list – Data, axes, metadata and original metadata.
- Return type
list of dicts
References
- Jackson2014
M. A. Jackson, M. A. Groeber, M. D. Uchic, D. J. Rowenhorst and M. De Graef, “h5ebsd: an archival data format for electron back-scatter diffraction data sets,” Integrating Materials and Manufacturing Innovation 3 (2014), doi: https://doi.org/10.1186/2193-9772-3-4.
-
kikuchipy.io.plugins.h5ebsd.
file_writer
(filename, signal, add_scan=None, scan_number=1, **kwargs)[source]¶ Write an
EBSD
orLazyEBSD
signal to an existing, but not open, or new h5ebsd file.Only writing to kikuchipy’s h5ebsd format is supported.
- Parameters
filename (
str
) – Full path of HDF file.signal (kikuchipy.signals.EBSD or kikuchipy.signals.LazyEBSD) – Signal instance.
add_scan (
Optional
[bool
]) – Add signal to an existing, but not open, h5ebsd file. If it does not exist it is created and the signal is written to it.scan_number (
int
) – Scan number in name of HDF dataset when writing to an existing, but not open, h5ebsd file.kwargs – Keyword arguments passed to
Group.require_dataset()
.
-
kikuchipy.io.plugins.h5ebsd.
get_desired_scan_groups
(file, scan_group_names=None)[source]¶ Get the desired HDF5 groups with scans within them.
- Parameters
- Returns
A list of the desired scan group(s) in the file.
- Return type
scan_groups
-
kikuchipy.io.plugins.h5ebsd.
get_scan_groups
(file)[source]¶ Return a list of the scan group names from an h5ebsd file.
-
kikuchipy.io.plugins.h5ebsd.
h5ebsd2signaldict
(scan_group, manufacturer, version, lazy=False)[source]¶ Return a dictionary with
signal
,metadata
andoriginal_metadata
from an h5ebsd scan.- Parameters
- Returns
scan – Dictionary with patterns,
metadata
andoriginal_metadata
.- Return type
-
kikuchipy.io.plugins.h5ebsd.
h5ebsdheader2dicts
(scan_group, manufacturer, version, lazy=False)[source]¶ Return three dictionaries in HyperSpy’s
hyperspy.misc.utils.DictionaryTreeBrowser
format, one with the h5ebsd scan header parameters as kikuchipy metadata, another with all datasets in the header as original metadata, and the last with info about scan size, image size and detector pixel size.- Parameters
- Return type
Tuple
[DictionaryTreeBrowser
,DictionaryTreeBrowser
,DictionaryTreeBrowser
]- Returns
md – kikuchipy
metadata
elements available in file.omd – All metadata available in file.
scan_size – Scan, image, step and detector pixel size available in file.
-
kikuchipy.io.plugins.h5ebsd.
hdf5group2dict
(group, dictionary=None, recursive=False, data_dset_names=None, **kwargs)[source]¶ Return a dictionary with values from datasets in a group in an opened HDF5 file.
- Parameters
- Returns
dictionary – Dataset values in group (and subgroups if recursive=True).
- Return type
-
kikuchipy.io.plugins.h5ebsd.
kikuchipyheader2dicts
(scan_group, md)[source]¶ Return scan metadata dictionaries from a kikuchipy h5ebsd file.
- Parameters
scan_group (Group) – HDF group of scan data and header.
md (
DictionaryTreeBrowser
) – Dictionary with empty fields from kikuchipy’s metadata.
- Return type
Tuple
[DictionaryTreeBrowser
,DictionaryTreeBrowser
,DictionaryTreeBrowser
]- Returns
md – kikuchipy
metadata
elements available in kikuchipy file.omd – All metadata available in kikuchipy file.
scan_size – Scan, image, step and detector pixel size available in kikuchipy file.
NORDIF¶
-
kikuchipy.io.plugins.nordif.
file_reader
(filename, mmap_mode=None, scan_size=None, pattern_size=None, setting_file=None, lazy=False)[source]¶ Read electron backscatter patterns from a NORDIF data file.
- Parameters
filename (
str
) – File path to NORDIF data file.scan_size (
Union
[None
,int
,Tuple
[int
, …]]) – Scan size in number of patterns in width and height.pattern_size (
Optional
[Tuple
[int
, …]]) – Pattern size in detector pixels in width and height.setting_file (
Optional
[str
]) – File path to NORDIF setting file (default is Setting.txt in same directory asfilename
).lazy (
bool
) – Open the data lazily without actually reading the data from disk until required. Allows opening arbitrary sized datasets. Default is False.
- Returns
scan – Data, axes, metadata and original metadata.
- Return type
list of dicts
-
kikuchipy.io.plugins.nordif.
file_writer
(filename, signal)[source]¶ Write an
EBSD
orLazyEBSD
object to a NORDIF binary file.- Parameters
filename (
str
) – Full path of HDF file.signal (kikuchipy.signals.EBSD or kikuchipy.signals.LazyEBSD) – Signal instance.
-
kikuchipy.io.plugins.nordif.
get_settings_from_file
(filename)[source]¶ Return metadata with parameters from NORDIF setting file.
- Parameters
filename (
str
) – File path of NORDIF setting file.- Return type
Tuple
[DictionaryTreeBrowser
,DictionaryTreeBrowser
,DictionaryTreeBrowser
]- Returns
md – Metadata complying with HyperSpy’s metadata structure.
omd – Metadata that does not fit into HyperSpy’s metadata structure.
scan_size – Information on image size, scan size and scan steps.
EMsoft EBSD master pattern HDF5¶
-
kikuchipy.io.plugins.emsoft_ebsd_master_pattern.
file_reader
(filename, energy_range=None, projection='spherical', hemisphere='north', lazy=False, **kwargs)[source]¶ Read electron backscatter diffraction master patterns from EMsoft’s HDF5 file format [Callahan2013].
- Parameters
filename (
str
) – Full file path of the HDF file.energy_range (
Optional
[range
]) – Range of beam energies for patterns to read. If None is passed (default), all available energies are read.projection (
str
) – Projection(s) to read. Options are “spherical” (default) or “lambert”.hemisphere (
str
) – Projection hemisphere(s) to read. Options are “north” (default), “south” or “both”. If “both”, these will be stacked in the vertical navigation axis.lazy (
bool
) – Open the data lazily without actually reading the data from disk until requested. Allows opening datasets larger than available memory. Default is False.kwargs – Keyword arguments passed to h5py.File.
- Returns
signal_dict_list – Data, axes, metadata and original metadata.
- Return type
list of dicts
References
- Callahan2013
P. G. Callahan and M. De Graef, “Dynamical Electron Backscatter Diffraction Patterns. Part I: Pattern Simulations,” Microscopy and Microanalysis 19 (2013), doi: https://doi.org/10.1017/S1431927613001840.