openNCEM Documentation

A collection of packages and tools for electron microscopy data analysis supported by the National Center for Electron Microscopy facility of the Molecular Foundry.

Installation

We recommend installing the ncempy package from PyPi:

>> pip install ncempy

If you wish to install the optional EDSTomo module then run:

pip install 'ncempy[edstomo]'

Components

The openNCEM collection comes with different components. The general functionality is condensed in libraries/packages, distinguished by programming language or framework. Tools are provided, addressing particular problems or tasks, leveraging the functionality provided in the libraries/packages.

ncempy package

Subpackages

ncempy.algo package
Submodules
ncempy.algo.align module
ncempy.algo.distortion module

Module to handle distortions in diffraction patterns.

ncempy.algo.distortion.filter_ring(points, center, rminmax)[source]

Filter points to be in a certain radial distance range from center.

Parameters:
  • points (np.ndarray) – Candidate points.
  • center (np.ndarray or tuple) – Center position.
  • rminmax (tuple) – Tuple of min and max radial distance.
Returns:

List of filtered points, two column array.

Return type:

np.ndarray

ncempy.algo.distortion.optimize_center(points, center, maxfev=1000, verbose=None)[source]

Optimize the center by minimizing the sum of square deviations from the mean radial distance.

Parameters:
  • points (np.ndarray) – The points to which the optimization is done (x,y coords in org image).
  • center (np.ndarray or tuple) – Initial center guess.
  • maxfev (int) – Max number of iterations forwarded to scipy.optimize.leastsq().
  • verbose (bool) – Set to get verbose output.
Returns:

The optimized center.

Return type:

np.ndarray

ncempy.algo.distortion.optimize_distortion(points, ns, maxfev=1000, verbose=False)[source]

Optimize distortions.

The orders in the list ns are first fitted subsequently and the result is refined in a final fit simultaneously fitting all orders.

Parameters:
  • points (np.ndarray) – Points to optimize to (in polar coords).
  • ns (tuple) – List of orders to correct for.
  • maxfev (int) – Max number of iterations forwarded to scipy.optimize.leastsq().
  • verbose (bool) – Set for verbose output.
Returns:

Optimized parameters according to ns.

Return type:

np.ndarray

ncempy.algo.distortion.points_topolar(points, center)[source]

Convert points to polar coordinate system.

Can be either in pixel or real dim, but should be the same for points and center.

Parameters:
  • points (np.ndarray) – Positions as two column array.
  • center (np.ndarray or tuple) – Origin of the polar coordinate system.
Returns:

Positions in polar coordinate system as two column array (r, theta).

Return type:

np.ndarray

ncempy.algo.distortion.rad_dis(theta, alpha, beta, order=2)[source]

Radial distortion due to ellipticity or higher order distortion.

Relative distortion, to be multiplied with radial distance.

Parameters:
  • theta (np.ndarray) – Angles at which to evaluate. Must be float.
  • alpha (float) – Orientation of major axis.
  • beta (float) – Strength of distortion (beta = (1-r_min/r_max)/(1+r_min/r_max).
  • order (int) – Order of distortion.
Returns:

Distortion factor.

Return type:

np.ndarray

ncempy.algo.distortion.residuals_center(param, data)[source]

Residual function for minimizing the deviations from the mean radial distance.

Parameters:
  • param (np.ndarray) – The center to optimize.
  • data (np.ndarray) – The points in x,y coordinates of the original image.
Returns:

Residuals.

Return type:

np.ndarray

ncempy.algo.distortion.residuals_dis(param, points, ns)[source]

Residual function for distortions.

Parameters:
  • param (np.ndarray) – Parameters for distortion.
  • points (np.ndarray) – Points to fit to.
  • ns (tuple) – List of orders to account for.
Returns:

Residuals.

Return type:

np.ndarray

ncempy.algo.fourier_operations module
ncempy.algo.fourier_operations.bandpass_filter(im, inner, outer, sigma)[source]

Apply a bandpass filter with smoothed edges in Fourier space.

Parameters:
  • im (ndarray, 2D) – The image to apply the bandpass filter to.
  • inner (float) – The lowest spatial frequency allowed in terms of the Nyquist frequency. Should be between 0 and 0.5.
  • outer (float) – The highest spatial frequency allowed in terms of the Nyquist frequency. Should be between 0 and 0.5.
  • sigma (float) – The standard deviation of the Gaussian smoothing to remove hard edges.
Returns:

The image with the bandpass filter applied.

Return type:

ndarray, 2D

ncempy.algo.fourier_operations.rotateImage(im, theta, pad=False)[source]

Use three shears in Fourier space to exactly (and reversibly) rotate an image.

Currently only works for square images.

Parameters:
  • im (ndarray, 2D) – The image to rotate
  • theta (float) – The angle to rotate by in radians
  • pad (bool) – Add padding to the image before rotating
ncempy.algo.fourier_operations.shearImage(im, dim, shear_factor)[source]

Exact shear of an image using the frequency space of the FFT of the image.

Currently only works for square images.

Parameters:
  • im (ndarray, 2D) – The image to shear.
  • dim (int) – The axis to shear along
  • shear_factor (float) – The amount of shear.
Returns:

The sheared image.

Return type:

ndarray

ncempy.algo.fourier_operations.shiftImage(im, shift)[source]

Exact shear of an image using the frequency space of the FFT of the image.

Parameters:
  • im (ndarray, 2D) – The image to shear.
  • shift (tuple of floats) – The number of pixels to shift in each direction
Returns:

The shifted image.

Return type:

ndarray

ncempy.algo.local_max module

Module to find local maxima in an image.

ncempy.algo.local_max.local_max(img, r, thresh)[source]

Find local maxima from comparing dilated and eroded images.

Calculates images with maximum and minimum within given radius. If the difference is larger than the threshold, the original pixel position with max value is detected as local maximum.

Parameters:
  • img (np.ndarray) – Input image.
  • r (int) – Radius for locality.
  • thresh (int or float) – Intensity difference threshold.
Returns:

Array of points.

Return type:

np.ndarray

ncempy.algo.local_max.points_todim(points, dims)[source]

Convert points from px coordinates to real dim.

Points are expected to be array indices for the first two dimensions in dims.

Parameters:
  • points (np.ndarray) – Points to convert.
  • dims (tuple) – Tuple of dimensions.
Returns:

Converted points.

Return type:

np.ndarray

ncempy.algo.math module

Module containing definitions of basic math functions, which can be used for fitting.

If you add functions here, do not forget to update the lookup table at the bottom as well.

ncempy.algo.math.const(x, param)[source]

Constant function.

f = param[0] .

Parameters:
  • x (np.ndarray) – Positions at which to evaluate function.
  • param (np.ndarray) – Necessary parameters.
Returns:

Values of function at x.

Return type:

(np.ndarray)

ncempy.algo.math.linear(x, param)[source]

Linear function.

f = param[0] * x + param[1] .

Parameters:
  • x (np.ndarray) – Positions at which to evaluate function.
  • param (np.ndarray) – Necessary parameters.
Returns:

Values of function at x.

Return type:

(np.ndarray)

ncempy.algo.math.lkp_funcs = {'const': (<function const>, 1), 'linear': (<function linear>, 2), 'powlaw': (<function powlaw>, 2), 'voigt': (<function voigt>, 4)}

Look-up table for functions implemented in this module. Functions are identified by strings, the entries give handles to the functions as well as the number of arguments necessary.

Type:(dict)
ncempy.algo.math.powlaw(x, param)[source]

Power law.

f = param[0] * x^param[1] .

Parameters:
  • x (np.ndarray) – Positions at which to evaluate function.
  • param (np.ndarray) – Necessary parameters.
Returns:

Values of function at x.

Return type:

(np.ndarray)

ncempy.algo.math.sum_functions(x, funcs, param)[source]

Sum functions in funcs over range x.

Parameters:
  • x (np.ndarray) – Positions at which to evaluate the functions.
  • funcs (list) – List of strings identifying function implemented in ncempy.algo.math.
  • param (np.ndarray) – Concatenated parameters for functions in funcs.
Returns:

Values of sum of functions at x.

Return type:

(np.ndarray)

ncempy.algo.math.voigt(x, param)[source]

Voigt peak function.

f = param[0]/(param[2] * sqrt(2*pi)) * Re( Faddeeva( ((x - param[1]) + i*param[3])/(param[2] * sqrt(2)) ) )

Parameters:
  • x (np.ndarray) – Positions at which to evaluate function.
  • param (np.ndarray) – Necessary parameters.
Returns:

Values of function at x.

Return type:

(np.ndarray)

ncempy.algo.moments module

Calculate the image raw moments.

param im:A 2D array of the image to calculate the moments for
type im:ndarray
param order:The order to calculate the moments to.
type order:int, default=3
returns:An ndarray of (order+1, order+1)
rtype:ndarray
ncempy.algo.multicorr_funcs module

Module to correlate two images, functionally written.

Todo

  • Cant use rfft2 currently. This gives one shift as 1/2 the value. How can this be improved to improve speed?
ncempy.algo.multicorr_funcs.dftUpsample(image_corr, upsample_factor, xy_shift)[source]

This performs a matrix multiply DFT around a small neighboring region of the initial correlation peak. By using the matrix multiply DFT to do the Fourier upsampling, the efficiency is greatly improved. This is adapted from the subfunction dftups found in the dftregistration function on the Matlab File Exchange.

https://www.mathworks.com/matlabcentral/fileexchange/18401-efficient-subpixel-image-registration-by-cross-correlation

The matrix multiplication DFT is from Manuel Guizar-Sicairos, Samuel T. Thurman, and James R. Fienup, “Efficient subpixel image registration algorithms,” Opt. Lett. 33, 156-158 (2008). http://www.sciencedirect.com/science/article/pii/S0045790612000778

Parameters:
  • image_corr (ndarray) – Correlation image between two images in Fourier space.
  • upsample_factor (int) – Scalar integer of how much to upsample.
  • xy_shift (list of 2 floats) – Single pixel shift between images previously computed. Used to center the matrix multiplication on the correlation peak.
Returns:

image_upsample – Upsampled image from region around correlation peak.

Return type:

ndarray

ncempy.algo.multicorr_funcs.imageShifter(g1, xy_shift)[source]

Multiply im by a plane wave that has the real space effect of shifting ifft2(G2) by [x, y] pixels.

Parameters

g1 : complex ndarray
The Fourier transform of an image.
xy_shift : list
A two element list of the shifts along each axis.
G2shift : complex ndarray
Fourier shifted image FFT
>> shiftIm0 = np.real(np.fft.ifft2(multicorr.imageShifter(np.fft.fft2(im0),[11.1,22.2]))) >> plt.imshow(shiftIm0)
ncempy.algo.multicorr_funcs.initial_correlation_image(g1, g2, method='cross', verbose=False)[source]

Generate correlation image at initial resolution using the method specified.

Parameters:
  • g1 (complex ndarray) – Fourier transform of reference image.
  • g2 (complex ndarray) – Fourier transform of the image to register (the kernel).
  • method (str, optional) – The correlation method to use. Must be ‘phase’ or ‘cross’ or ‘hybrid’ (default = ‘cross’)
  • verbose (bool, default is False) – Print output.
Returns:

imageCorr – Correlation array which has not yet been inverse Fourier transformed.

Return type:

ndarray complex

ncempy.algo.multicorr_funcs.upsampleFFT(image_init, upsample_factor)[source]

This does a Fourier upsample of the imageInit. imageInit is the Fourier transform of the correlation image. The function returns the real space correlation image that has been Fourier upsampled by the upsample_factor. An upsample factor of 2 is generally sufficient.

The way it works is that it embeds imageInit in a larger array of zeros, then does the inverse Fourier transform to return the Fourier upsampled image in real space.

Parameters:
  • image_init (ndarray complex) – The image to be Fourier upsampled. This should be in the Fourier domain.
  • upsample_factor (int) – THe upsample factor (usually 2).
Returns:

imageUpsampleReal – The inverse Fourier transform of imageInit upsampled by the upsampleFactor.

Return type:

ndarray complex

imageSize = imageInit.shape imageUpsample = np.zeros(tuple((i*upsampleFactor for i in imageSize))) + 0j imageUpsample[:imageSize[0], :imageSize[1]] = imageInit imageUpsample = np.roll(np.roll(imageUpsample, -int(imageSize[0]/2), 0), -int(imageSize[1]/2),1) imageUpsampleReal = np.real(np.fft.ifft2(imageUpsample)) return imageUpsampleReal

ncempy.algo.multicorr_funcs.upsampled_correlation(image_corr, upsample_factor, verbose=False)[source]

Upsamples the correlation image by a set integer factor upsample_factor. If upsample_factor == 2, then it is naively Fourier upsampled. If the upsample_factor is higher than 2, then it uses dftUpsample, which is a more efficient way to Fourier upsample the image.

Parameters:
  • image_corr (ndarray complex) – Fourier transformed correlation image returned by initial_correlation_image.
  • upsample_factor (int) – Upsampling factor.
  • verbose (bool) – Provide output for debugging
Returns:

xyShift – Shift in x and y of G2 with respect to G1.

Return type:

list

ncempy.algo.radial_profile module

Module to calculate radial profiles.

The settings dict used for the evaluation has the following entries:

  • lmax_r (int): Radius used for local maximum detection in [px].
  • lmax_thresh (int): Threshold for local maximum detection in [intensity].
  • lmax_cinit (int, int): Initial guess for center in [px].
  • lmax_range (float, float): r range to cinit used to filter local maxima [dims].
  • ns (int, …): Distortion orders to correct.
  • fit_rrange (float, float): r range used to fit radial profile [dims].
  • back_xs (float, …): Fixpoints for fitting background [dims].
  • back_xswidth (float): Range around fixpoints to take into account [dims].
  • back_init (float): Initial guess for subtracting background.
  • fit_funcs (str, …): List of functions to model radial profile.
  • fit_init (float, …): Initial guess for fitting.

optional (set to None to use defaults):

  • plt_imgminmax (float, float): Relative range to plot the img.
  • rad_rmax (float): Max r to be used in radial profile [dims].
  • rad_dr (float): Stepsize for r-axis in radial profile [dims].
  • rad_sigma (float): Sigma for Gaussian used as kernel density estimator [dims].
  • mask (np.ndarray): Binary image as img, 0 for pixels to exclude.
  • fit_maxfev (int): Maxfev forwarded to scipy optimize.
ncempy.algo.radial_profile.calc_polarcoords(center, dims, ns=None, dists=None)[source]

Calculate the polar coordinates for an image of given shape.

Center is assumed to be in real coordinates (if not just fake the dims). Distortions are corrected if ns and corresponding dists are given.

Parameters:
  • center (np.ndarray/tuple) – Center of polar coordinate system.
  • dims (tuple) – Tuple of dimensions.
  • ns (tuple) – List of distortion orders to correct for.
  • dists (np.ndarray) – Parameters for distortions.
Returns:

Polar coordinates (r, theta) as two np.ndarrays with same dimensions as original image.

Return type:

tuple

ncempy.algo.radial_profile.calc_radialprofile(img, rs, rMax, dr, rsigma, mask=None)[source]

Calculate the radial profile using Gaussian kernel density estimator.

It is suggested to use an rMax such that all directions are still in the image, otherwise outer areas will contribute differently and lead to different signal-to-noise. A value of dr corresponding to 1/10 px and rsigma corresponding to 1 px are good parameters.

Parameters:
  • img (np.ndarray) – Image to take radial profile of intensity from.
  • rs (np.ndarray) – Array containing the radial distance, same shape as img.
  • rMax (float) – Maximum radial distance used for the radial profile.
  • dr (float) – Step size for r-axis of radial distance.
  • rsigma (float) – Sigma for Gaussian used as kernel density estimator.
  • mask (np.ndarray) – Mask image, values at 0 are excluded from radial profile.
Returns:

Tuple of radial and intensity axes.

Return type:

tuple

ncempy.algo.radial_profile.correct_distortion(img, dims, center, ns, dists)[source]

Give corrected version of diffraction pattern with respect to distortions.

Parameters:
  • img (np.ndarray) – Image (2D).
  • dims (tuple) – Dimensions tuple.
  • center (np.ndarray or tuple) – Center to be used.
  • ns (tuple) – List of distortion orders.
  • dists (np.ndarray) – Distortion parameters.
Returns:

Corrected diffraction pattern.

Return type:

np.ndarray

ncempy.algo.radial_profile.fit_radialprofile(r, intens, funcs, init_guess, maxfev=None)[source]

Fit the radial profile.

Convenience wrapper for fitting.

Parameters:
  • r (np.ndarray) – r-axis of radial profile.
  • intens (np.ndarray) – Intensity-axis of radial profile.
  • funcs (tuple) – List of functions.
  • init_guess (np.ndarray) – Initial guess for parameters of functions in funcs.
  • maxfev (int) – maximum function calls (for scipy.optimize.leastsq)
Returns:

Optimized parameters.

Return type:

np.ndarray

ncempy.algo.radial_profile.residuals_fit(param, r, intens, funcs)[source]

Residual function to fit radial profile with flexibility.

The arguments for the single functions are automatically cut from the fit parameters.

Parameters:
  • param (np.ndarray) – Fit parameters.
  • r (np.ndarray) – r-axis.
  • intens (np.ndarray) – Intensity-axis.
  • funcs (tuple) – List of functions to include.
Returns:

Residuals.

Return type:

np.ndarray

ncempy.algo.radial_profile.run_singleImage(img, dims, settings, show=False)[source]

Evaluate a single ring diffraction pattern with given settings.

Parameters:
  • img (np.ndarray) – Image.
  • dims (tuple) – Corresponding dim vectors.
  • settings (dict) – Dict of settings necessary for the evaluation.
  • show (bool) – Set to directly show plots interactively.
Returns:

(np.ndarray)

Return type:

Optimized parameters of fitting the radial profile according to the settings.

ncempy.algo.rebin module

Rebin a 2D array. From stackoverflow: https://stackoverflow.com/questions/4624112/grouping-2d-numpy-array-in-average

This should be used carefully. There is not antialiasing applied which could produce odd effects for some data sets (such as high resolution data)

param im:2D array to reduce
type im:ndarray
param f:The factor to rebin by. Must be integer
type f:int
param funcType:The type of reduction. mean or sum are implemented.
type funcType:str
returns:The 2D array with the new size.
rtype:ndarray
ncempy.algo.gaussND module

A module with Gaussian functions for creating Gaussian distributions in 1D, 2D and 3D. There are also functions with the same name and _FIT which can be used to fit an intensity distribution in a ndarray with the corresponding Gaussian distribution using scipy’s curve_fit.

A 1D Lorentz and a 1D gaussLorentz function is also included for fitting zero loss peaks in EELS spectra.

These functions were written assuming you use np.meshgrid() with indexing = ‘xy’. If you use ‘ij’ indexing then be careful about how you pass in the x and y coordinates.

ncempy.algo.gaussND.MatrixQuaternionRot(vector, theta)[source]

Quaternion tp rotate a given theta angle around the given vector.

adapted from code by Yongsoo Yang, yongsoo.ysyang@gmail.com

Parameters:
  • vector (ndarray, 3-element) – A non-zero 3-element numpy array representing rotation axis
  • theta (float, degrees) – Rotation angle in degrees
Returns:

Authoryongsoo.ysyang@gmail.com

Return type:

Yongsoo Yang, Dept. of Physics and Astronomy, UCLA

ncempy.algo.gaussND.gauss1D(x, x0, sigma)[source]

Returns the value of a gaussian at a 2D set of points for the given standard deviation with maximum normalized to 1.

Parameters:
  • x (ndarray) – A vector of size (N,) of points for the x values
  • x0 (float) – The center of the Gaussian.
  • sigma (float) – Standard deviation of the Gaussian.
Returns:

g – A vector of size (N,) of the Gaussian distribution evaluated at the input x values.

Return type:

ndarray

Note

Calculate the half width at half maximum (HWHM) as >> HWHM = sqrt(2*log(2))*stDev ~ 1.18*stDev or if x goes from -1 to 1 >> 0.5*size(x)*stDev

ncempy.algo.gaussND.gauss2D(x, y, x0, y0, sigma_x, sigma_y)[source]

Calculates the value of a Gaussian at a 2D set of points for the given standard deviations with maximum normalized to 1. The Gaussian axes are assumed to be 90 degrees from each other.

Parameters:
  • y (x,) – A 2D array of size (M,N) of points for the x values. Use x, y = np.meshgrid(range(M),range(N),indexing=’xy’).
  • x0 (float) – The center of the Gaussian along the x direction.
  • y0 (float) – The center of the Gaussian along the y direction.
  • sigma_x (float) – Standard deviation of the Gaussian along x.
  • sigma_y (float) – Standard deviation of the Gaussian along y.
Returns:

g – A ndarray of size (N, M) of the Gaussian distribution evaluated at the input x values.

Return type:

ndarray

Note

The Gaussian is normalized such that the peak == 1. To normalize the integral divide by 2 * np.pi * sigma_x * sigma_y

ncempy.algo.gaussND.gauss2D_FIT(xy, x0, y0, sigma_x, sigma_y)[source]

Version of gauss2D used for fitting (1 x_data input (xy) and flattened output). Returns the value of a gaussian at a 2D set of points for the given standard deviation with maximum normalized to 1. The Gaussian axes are assumed to be 90 degrees from each other.

Parameters:
  • xy (tuple) – A (N,M) shaped tuple containing the vectors of points for the evaluation points.
  • x0 (float) – The x center of the Gaussian.
  • y0 (float) – The y center of the Gaussian.
  • sigma_x (float) – The standard deviation of the Gaussian along x.
  • sigma_y (float) – The standard deviation of the Gaussian along y.
Returns:

g2_norm – The Gaussian distribution with maximum value normalized to 1. The 2D ndarray is reshaped into a (N*M,) array for use in fitting functions in numpy and scipy.

Return type:

ndarray

ncempy.algo.gaussND.gauss2D_poly_FIT(xy, x0, y0, A, B, C)[source]

Returns the flattened values of an elliptical gaussian at a 2D set of points for the given polynomial pre-factors with maximum normalized to 1. The Gaussian axes are assumed to be 90 degrees from each other.

The matrix looks like [[A,B],[B,C]]. See https://en.wikipedia.org/wiki/Gaussian_function

Parameters:
  • xy (tuple) – Evaluation points along x and y of shape (N,M)
  • x0 (float) – Center of the maximum along x.
  • y0 (float) – Center of the maximum along y.
  • A (float) – The A pre-factor for the polynomial expansion in the exponent.
  • B (float) – The B pre-factor for the polynomial expansion in the exponent.
  • C (float) – The C pre-factor for the polynomial expansion in the exponent.
Returns:

g2_norm – A (N, M) sized 2D ndarray with a Gaussaian distribution rotated.

Return type:

ndarray

ncempy.algo.gaussND.gauss2D_theta(x, y, x0, y0, sigma_x, sigma_y, theta)[source]

Returns the value of a generalized Gaussian at a 2D set of points for the given standard deviation with maximum normalized to 1. The Gaussian axes (assumed to be 90 degrees from each other) can be oriented at different angles to the output array axes using theta.

Parameters:
  • x (ndarray or list) – Evaluation points along x of size (N,)
  • y (ndarray) – Evaluation points along y of size (M,)
  • x0 (float) – Center of the maximum along x.
  • y0 (float) – Center of the maximum along y.
  • sigma_x (float) – Standard deviation along x.
  • sigma_y (float) – Standard deviation along y.
  • theta (float) – The rotation of the Gaussian principle axes from the array horizontal and vertical. In radians.
Returns:

g2_norm – A (N, M) sized 2D ndarray with a Gaussian distribution rotated.

Return type:

ndarray

ncempy.algo.gaussND.gauss2D_theta_FIT(xy, x0, y0, sigma_x, sigma_y, theta)[source]

Version of gauss2D_theta used for fitting (1 x_data input and flattened output). Returns the value of a gaussian at a 2D set of points for the given standard deviation with maximum normalized to 1. The Gaussian axes can be oriented at different angles (theta) in radians.

Parameters:
  • xy (tuple) – Evaluation points along x and y of shape (N,M)
  • x0 (float) – Center of the maximum along x.
  • y0 (float) – Center of the maximum along y.
  • sigma_x (float) – Standard deviation along x.
  • sigma_y (float) – Standard deviation along y.
  • theta (float) – The rotation of the Gaussian principle axes from the array horizontal and vertical. In radians.
Returns:

g2_norm – A (N, M) sized 2D ndarray with a Gaussian distribution rotated.

Return type:

ndarray

ncempy.algo.gaussND.gauss3D(x, y, z, x0, y0, z0, sigma_x, sigma_y, sigma_z)[source]

Returns the value of a gaussian at a 2D set of points for the given standard deviations with maximum normalized to 1. The Gaussian axes are assumed to be 90 degrees from each other.

Note

Be careful about the indexing used in meshgrid and the order in which you pass the x, y, z variables in.

Parameters:
  • y, z (x,) – 2D arrays of points (from meshgrid)
  • y0, z0 (x0,) – The x, y, z centers of the Gaussian
  • sigma_y, sigma_z (sigma_x,) – The standard deviations of the Gaussian.
Returns:

g3_norm – A 3D ndarray

Return type:

ndarray

ncempy.algo.gaussND.gauss3DGEN_FIT((x, y, z), x0, y0, z0, sigma_x, sigma_y, sigma_z, Angle1, Angle2, Angle3, BG, Height)[source]

Returns the value of a gaussian at a 3D set of points for the given sub-pixel positions with standard deviations, 3D Eular rotation angles, constant Background value, and Gaussian peak height.

adapted from code by Yongsoo Yang, yongsoo.ysyang@gmail.com

Note

This is a work in progress. Needs more testing.

Parameters:
  • xyz (tuple of 3 np.ndarray) – 3D arrays of points (from meshgrid) combined in a tuple
  • y0, z0 (x0,) – The x, y, z centers of the Gaussian
  • sigma_x,sigma_y,sigma_z (float) – standard deviations along x,y,z direction before 3D angular rotation
  • Angle2, Angle3 (Angle1,) – Tait-Bryan angles in ZYX convention for 3D rotation in degrees
  • BG (float) – Background
  • Height (float) – The peak height of the Gaussian function
Returns:

The 3D Gaussian

Return type:

ndarray, 3D

ncempy.algo.gaussND.gauss3D_FIT((x, y, z), x0, y0, z0, sigma_x, sigma_y, sigma_z)[source]

Returns the value of a gaussian at a 2D set of points for the given standard deviations with maximum normalized to 1. The Gaussian axes are assumed to be 90 degrees from each other.

xyz - x0, y0, z0 = the x, y, z centers of the Gaussian sigma_x, sigma_y, sigma_z = The std. deviations of the Gaussian.

Note

Be careful about the indexing used in meshgrid and the order in which you pass the x, y, z variables in.

Parameters:
  • xyz (tuple of ndarrays) – A tuple containing the 3D arrays of points (from meshgrid)
  • y0, z0 (x0,) – The x, y, z centers of the Gaussian
  • sigma_y, sigma_z (sigma_x,) – The standard deviations of the Gaussian.
Returns:

g3_norm – A flattened array for fitting.

Return type:

ndarray

ncempy.algo.gaussND.gauss3D_poly(x, y, z, x0, y0, z0, A, B, C, D, E, F)[source]

gauss3Dpoly_FIT((x,y,z),x0,y0,z0,A,B,C,D,E,F) Returns the value of a gaussian at a 2D set of points for the given standard deviations with maximum normalized to 1. The Gaussian axes are not locked to be 90 degrees.

Parameters:
  • y, z (x,) – 3D arrays of points (from meshgrid)
  • y0, z0 (x0,) – The x, y, z centers of the Gaussian
  • B, C, D, E, F (A,) – The polynomial values for the fit
Returns:

The 3D Gaussian

Return type:

ndarray, 3D

ncempy.algo.gaussND.gauss3D_poly_FIT(xyz, x0, y0, z0, A, B, C, D, E, F)[source]

gauss3Dpoly_FIT((x,y,z),x0,y0,z0,A,B,C,D,E,F) Returns the value of a gaussian at a 2D set of points for the given standard deviations with maximum normalized to 1. The Guassian axes are not locked to be 90 degrees.

Parameters:
  • xyz (tuple of ndarrays) – 3D arrays of points (from meshgrid) combined in a tuple
  • y0, z0 (x0,) – The x, y, z centers of the Gaussian
  • B, C, D, E, F (A,) – The polynomial values for the fit
Returns:

The 3D Gaussian

Return type:

ndarray, 3D

ncempy.algo.gaussND.gaussLorentz1D(x, x0, w)[source]

A Gaussian-Lorentzian function in one dimension.

Parameters:
  • x (ndarray or list) – A 1D vector of points for the dependent variable
  • x0 (float) – The center of the peak maximum
  • w (float) – The parameter that modifies the width of the distribution.
Returns:

lg – A vector of size (N,) of the Lorentzian Gaussian distribution evaluated at the input x values.

Return type:

ndarray

ncempy.algo.gaussND.lorentz1D(x, x0, w)[source]

Returns the probability density function of the Lorentzian (aka Cauchy) function with maximum at 1.

Parameters:
  • x (ndarray or list) – A 1D vector of points for the dependent variable
  • x0 (float) – The center of the peak maximum
  • w (float) – The parameter that modifies the width of the distribution.
Returns:

l – A vector of size (N,) of the Lorentzian distribution evaluated at the input x values.

Return type:

ndarray

ncempy.algo.peak_find module

Module to find local maxima in 2D and 3D data such as images and density maps.

Many functions are based on work by Colin Ophus, cophus@lbl.gov and his excellent RealSpaceLattice01.m code.

author: Peter Ercius, percius@lbl.gov

ncempy.algo.peak_find.applyLatticeLimit(lattice, bounds)[source]

Remove lattice points outside the data bounds. For 2D and 3D data.

Parameters:
  • lattice (ndarray; (N, 2) or (N, 3)) – From lattice2D
  • bounds (tuple,) – Minimum and maximum for axes 0 and 1 (min0, max0, min1, max1) or axes 0, 1 and 2 (min0, max0, min1, max1, min2, max2)
Returns:

Same as lattice input except only containing points within the bounds specified. M <= N

Return type:

ndarray; (M, 2) or (M, 3)

ncempy.algo.peak_find.calculate_unit_cell(image, lattice, u, v, unit_cell_size)[source]

IN PROGRESS Calculate a unit cell using the input lattice and lattice parameters. Any unit cell at the edge of the image is not used as part of the unit cell will be invalid.

Parameters:
  • image (ndarray) – The image containing the unit cell in a regular pattern (i.e. a lattice)
  • lattice (ndarray) – The starting position of each unit cell in the image.
  • v (u,) – The u and v vectors for the fitted lattice
  • unit_cell_size (int or tuple) – If int then this is the number of pixels in the unit cell along the u vector. The size of the v vector will be automatically calculated to make the unit cell pixels square. If tuple then the number of pixels along each u and v dimension are used directly.
Returns:

An array of size unit_cell_size which contains the average value of each position in the unit cell for each peak position and u,v, lattice point.

Return type:

ndarray

ncempy.algo.peak_find.doubleRoll(image, vec)[source]

Roll an 2D ndarray (e.g. an image) along each dimension.

Parameters:
  • image (np.ndarray) – The 2D array to roll.
  • vec (iterable) – The number of time to apply the roll to the rows and columns where the tuple has the structure (row, col). Values in tuple must be int.
Returns:

The rolled image.

Return type:

np.ndarray

ncempy.algo.peak_find.enforceMinDist(positions, intensities, minDistance)[source]

Remove peaks that violate a minimum distance requirement. If two peaks are closer together than the minDistance then the one with the highest intensity is chosen.

Works for 2D and 3D data sets.

Parameters:
  • positions (np.ndarray) – An array of shape Nx3 or Nx2 containing the positions of the peaks to test
  • intensities (np.ndarray) – An array of shape Nx1 with the peak intensities.
  • minDistance (int) – The minimum distance in pixels
Returns:

validPeaks – An ndarray of Mx3 or Mx2 (M <= N) with invalid peaks removed.

Return type:

np.ndarray

ncempy.algo.peak_find.fit_peaks_gauss2d(image, peaks, cutOut, init, bounds, remove_edge_peaks=True)[source]

Fit peaks to a 2D Gaussian. The Gaussian function is gaussND.gauss2D()

Todo: Write tests. This function is copied from test_peakFind.ipynb Todo: Add example test_peakFind.ipynb from my jupyter notebook folder

Parameters:
  • image (np.ndarray) – The 2D image with the intensities to fit to
  • peaks (np.ndarray) – The peaks in a ndarray of shape (N,2) where N is the number of peaks. Should match output of peakFind.peakFind2D.
  • cutOut (int) – The size (+/-) of the region around each peak to fit.
  • init (tuple) – A tuple of initial sigma values in x and y. (sigma_x, sigma_y)
  • bounds (tuple) – A tuple of 2 tuples indicating the upper and lower bounds for the fitting. See scipy.optimize.curve_fit for more information. Ordered as ( (center_x_low, center_y_low, sigma_x_low, sigma_y_low), (center_x_high, center_y_high, sigma_x_high, sigma_y_high))
  • remove_edge_peaks (bool, default = True) – Peak positions at the edge of the image are set to NaN and are removed. If you want to have those positions returned as nan set this to False.

Note

This function uses np.meshgrid and indexing=’ij’ internally.

Returns:
Returns a tuple of 3 arrays:
[0] optimized peak positions as a (M, 2) ndarray. M <= N. Peaks on the edge of the image are removed. [1] peak intensity value interpolated at the optimized peak position. [2] the fitting values for each peak
Return type:tuple
ncempy.algo.peak_find.fit_peaks_gauss3d(volume, peaks, cutOut, init, bounds, remove_edge_peaks=True)[source]

Fit peaks to a 2D Gaussian. The Gaussian function is gaussND.gauss2D()

Todo: Write tests. This function is copied from test_peakFind_3d.ipynb

Parameters:
  • volume (np.ndarray) – The 3D image with the intensities to fit to.
  • peaks (np.ndarray) – The peaks in a ndarray of shape (N, 3) where N is the number of peaks. Should match output of peakFind.peakFind3D.
  • cutOut (int) – The size (+/-) of the region around each peak to fit.
  • init (tuple) – A 3-tuple of initial sigma values in x and y. (sigma_x, sigma_y, sigma_z)
  • bounds (tuple) – A 2-tuple of 6-tuples indicating the upper and lower bounds for the fitting. See scipy.optimize.curve_fit for more information. Ordered as ( (center_x_low, center_y_low, center_z_low, sigma_x_low, sigma_y_low, sigma_z_low), (center_x_high, center_y_high, center_z_high, sigma_x_high, sigma_y_high, sigma_z_high))
  • remove_edge_peaks (bool, default = True) – Peak positions at the edge of the image are set to NaN and are removed. If you want to have those positions returned as nan set this to False.

Note

This function uses np.meshgrid and indexing=’ij’ internally.

Returns:
Returns a tuple of 3 arrays:
[0] optimized peak positions as a (M, 3) ndarray. M <= N. Peaks on the edge of the image are removed unless otherwise specified. Then the edge peaks are returned as NAN values. [1] peak intensity value interpolated at the optimized peak position. [2] the local fitting values for each peak.
Return type:tuple
ncempy.algo.peak_find.generateLatticeFromRefinement(origin, u, v, ab, fraction=(1, 1))[source]

Generate lattice positions from refined output of refineLattice2D. This can be used to generate all positions when fraction is used in refineLattice2D

Parameters:
  • origin (tuple) – Origin. 2-tuple
  • u (tuple) – Initial u vector. 2-tuple
  • v (tuple) – Initial v vector. 2-tuple
  • fraction (tuple) – The fraction used in refineLattice2D. 2-tuple
  • ab (np.ndarray) – The set of positions in terms of u and v lattice vectors. (number, position) (M, 2)
  • fraction – The fractional coordinates for the unit cell. 2-tuple
Returns:

The lattice site positions of the given lattice vectors of shape (M, 2)

Return type:

np.ndarray

ncempy.algo.peak_find.lattice2D(u, v, a, b, origin, num_points)[source]

A modified version of peakFind.lattice2D_norm which can use non normalized u,v,vectors

Parameters:
  • v (u,) – 2 element vectors defining the lattice directions
  • b (a,) – values to multiply each vector by (if u,v,w are not normalized then set these to 1)
  • origin (tuple) – The origin in the format (x0,y0)
  • num_points (tuple) – 2 element tuple of the number of repeats of the lattice along each direction
Returns:

The set of points in the lattice. Size (num_points[0] * num_points[1], 2)

Return type:

ndarray

ncempy.algo.peak_find.lattice2D_2(u, v, a, b, xy0, numPoints)[source]

A modified version of lattice2D which uses non-normalized u, v vectors

ncempy.algo.peak_find.lattice2D_norm(u, v, a, b, origin, num_points)[source]

Returns a set of points in a lattice according to the u, v unit vectors (vectors are normalized internally) and lengths a,b centered at origin. The lattice has num_points along each u,v vector.

Parameters:
  • v (u,) – 2 element tuples defining the lattice directions as vectors.
  • b (a,) – values to multiply each vector by (if u,v are not normalized then set these to 1)
  • origin (2-tuple) – The origin in the format (x0, y0)
  • num_points (2-tuple) – 2 element tuple of the number of repeats of the lattice along each direction
Returns:

The set of points in the lattice. Size (num_points[0] * num_points[1], 2)

Return type:

ndarray

ncempy.algo.peak_find.lattice3D(u, v, w, a, b, c, origin, num_points)[source]

Returns a set of points in a lattice according to the u, v, w vectors and lengths a,b centered at origin. The lattice has num_points along each u,v,w vector.

Parameters:
  • v, w (u,) – 3 element vectors defining the lattice directions
  • b, c (a,) – Values to multiply each vector by (if u,v,w are not normalized then set these to 1)
  • origin (tuple) – The origin in the format (x0,y0,z0)
  • num_points (tuple) – 3 element tuple of the number of repeats of the lattice along each direction
Returns:

xyz – A (N,3) shaped set of coordinates

Return type:

ndarray

ncempy.algo.peak_find.lattice3D_2(u, v, w, a, b, c, xyz0, numPoints)[source]

A modified version of lattice3D which uses non-normalized u,v vectors

ncempy.algo.peak_find.latticeDisplacements(peaks, u, v, origin)[source]

Find the displacements of the experimental peaks to the fitted lattice parameters

Parameters:
  • peaks (ndarray) – The experimental peaks with shape (num_peaks, 2). For example, the output of fit_peaks_gauss2d.
  • v (u,) – The u and v vectors for the fitted lattice
  • origin (tuple) – The origin of the fitted lattice for u, v vectors
Returns:

The displacement of each peak from the expected lattice position.

Return type:

ndarray

ncempy.algo.peak_find.match_lattice_peaks(peaks, u, v, origin)[source]

Find the matching lattice points to the experimental peaks. This is useful to generate all fitted lattice points for a set of experimental peaks. This can then be used directly to calculate displacements for example.

Parameters:
  • peaks (ndarray) – The experimental peaks as a (num_peaks, 2) ndarray
  • v (u,) – The u and v vectors for the fitted lattice
  • origin (tuple) – The origin of the fitted lattice for u, v vectors
Returns:

An array of shape (num_peaks, 2) where each coordinate is the closest lattice point for each input peak.

Return type:

ndarray

ncempy.algo.peak_find.peak3View(fg, vol, peak_positions)[source]

Plot 3 orthogonal slices and the corresponding peaks in each slice.

Not fully test. Unsure whether the slices and peaks are exactly the same.

todo: move this to viz

Parameters:
  • fg (matplotlib figure) – The figure to use. This allows the user to reuse figures
  • vol (ndarray (3D)) – The volume to slice along each orthogonal direction. Set up for [X,Z,Y] orientation
  • peak_positions (ndarray, list, tuple (1D)) – This selects the point shown in all three slices
Returns:

Return type:

None

ncempy.algo.peak_find.peakFind2D(image, threshold)[source]

Find peaks in a 2D image based on the roll technique. A threshold is first applied to remove low intensity noisy peaks.

Parameters:
  • image (ndarray) – An 2D ndarray of image intensities with peaks to find.
  • threshold (float) – Threshold [0,1) the image to remove noisy peaks based on max intensity in the image
Returns:

positions – A Nx2 array with the index locations of the n number of peaks

Return type:

array

ncempy.algo.peak_find.peakFind3D(vol, threshold)[source]

Find peaks in a 3D volume based on the roll technique. A threshold is first applied to remove low intensity noisy peaks.

Parameters:
  • vol (ndarray) – An 3D ndarray of image intensities with peaks to find.
  • threshold (float) – Threshold [0,1) the image to remove noisy peaks based on max intensity in the image
Returns:

positions – A Nx3 array with the index locations of the n number of peaks

Return type:

ndarray

ncempy.algo.peak_find.peakPlot3D(X, Y, Z, mkr, myAxes3D)[source]

Plot a set of peaks in a 3D plot using matplotlib. See the example below for how to set up a figure as input to this function using matplotlib Axes3D.

todo: move this to viz

Parameters:
  • X (np.ndarray) – The x positions as a 1D array.
  • Y (np.ndarray) – The y positions as a 1D array.
  • Z (np.ndarray) – The Z positions as a 1D array.
  • mkr (string) – String indicating the marks type and color. Passed directly to pyplot.plt command (ex. ‘go’ for green circles)
  • myAxes3D (mpl_toolkits.mplot3d.Axes3d) – An Axes3D object. (see example below)
Returns:

Return type:

None

Example

This function requires the input of a Axes3D to plot a set of peaks. See below how to import and set up a figure for use with this function. >> import matplotlib.pyplot as plt >> import mpl_toolkits.mplot3d >> from ncempy.algo import peak_find >> fg1 = plt.figure() >> ax1 = mpl_toolkits.mplot3d.Axes3D(fg1) >> peak_find.peakPlot3D(peakList[:,2], peakList[:,1], peakList[:,0], ‘go’, ax1)

ncempy.algo.peak_find.peaksToImage(peakList, imShape, gaussSigma, gaussSize, indexing='ij')[source]

Place 2D Gaussian at a set of peak positions.

Parameters:
  • peakList (np.ndarray) – Array of peak positions of size (numPeaks, 2).
  • imShape (tuple) – 2 element tuple of the shape of the image containing the peaks
  • gaussSigma (tuple) – 2 element tuple with that sigma parameters passed to the gaussND.gauss2D() function (sigX,sigY)
  • gaussSize (tuple) – 2 element tuple describing the 2D size of the simulated gaussian box. Must be odd.
  • indexing (str, default = 'ij') – ‘ij’ or ‘xy’ indexing to pass to np.meshgrid.
Returns:

A image containing Gaussian peaks at each position indicated in peakList parameter.

Return type:

np.ndarray

ncempy.algo.peak_find.peaksToVolume(peakList, volShape, gaussSigma, gaussSize, indexing='ij')[source]

Place 3D Gaussian at a set of peak positions.

Parameters:
  • peakList (np.ndarray) – Array of peak positions of size (numPeaks, 3).
  • volShape (tuple) – 3 element tuple of the shape of the volume containing the peaks
  • gaussSigma (tuple) – 3 element tuple with that sigma parameters passed to the gaussND.gauss3D() function (sigX,sigY,sigZ)
  • gaussSize (tuple) – 3 element tuple describing the 3D size of the simulated gaussian box. Must be odd.
  • indexing (str) – String to pass to np.meshgrid to indicate which indexing is used. Either ‘ij’ or ‘xy’ are possible. Default is ‘ij’ for this module.
Returns:

A volume containing gaussian peaks at each position indicated in peakList parameter.

Return type:

np.ndarray

ncempy.algo.peak_find.refineLattice2D(or0, u0, v0, pos, fraction=(1, 1), max_iter=30, refine_locally=True, verbose=False, num_unit_cells=(1, 1))[source]

Refine lattice based on measurements and initial guess at lattice vectors. This code is designed to work only with square lattices. Hexagonal and more complex lattice might not work.

Parameters:
  • or0 (tuple) – Initial guess of origin. 2-tuple
  • u0 (tuple) – Initial u vector. Use num_unit_cells to average over several lengths of the lattice vector in this direction. 2-tuple
  • v0 (tuple) – Initial v vector. Use num_unit_cells to average over several lengths of the lattice vector in this direction. 2-tuple
  • pos (np.ndarray) – The set of positions to fit the lattice to. (number, position) (M, 2)
  • refine_locally (bool) – Refine locally near the origin before using all positions Locally is considered 2 time the larger vector of u0 or v0.
  • fraction (tuple) – Site fraction to take into account peak positions inside the unit cell. FCC imaged along [100] for example would need site_fraction = (2, 2).
  • max_iter (int) – The maximum number of iterations to run to refine. This usually converges in a few iterations. Default is 30.
  • num_unit_cells (tuple) – The number of unit cells the initial guess of u0 and v0 extend over. 2-tuple (num_u0, num_v0) with default
  • verbose (bool) – Print out helpful information.
Returns:

Tuple of of 4 np.ndarray (origin, u, v, ab) where ab is the site positions in fractions of u and v.

Return type:

tuple

ncempy.algo.peak_find.refineLattice3D(or0, u0, v0, w0, pos, fraction=(1, 1, 1), max_iter=30, refine_locally=True)[source]

Refine lattice based on measurements and initial guess at lattice vectors

Warning

Tested code but not production yet. This is a work un progress.

Parameters:
  • or0 (3-element tuple) – Origin.
  • u0 (tuple) – Initial u vector.
  • v0 (tuple) – Initial v vector
  • w0 (tuple) – Initial w vector
  • pos (array (number, position) (M,3)) – The set of positions to fit the lattice to.
  • refine_locally (bool) – Refine locally near the origin before using all positions Locally is considered 2 time the larger vector of u0 or v0.
  • fraction (tuple) – Site fraction to take into account peak positions inside the unit cell. FCC imaged along [100] for example would need site_fraction = (2, 2).
  • max_iter (int) – Maximum number of iterations.
Returns:

Tuple of refined origin, u, v optimized values as np.ndarray and the site positions in fractions of u and v.

Return type:

tuple of 4 np.ndarray (origin, u, v, ab)

ncempy.algo.peak_find.remove_xrays(imageOriginal, threshold, size_median_filter=(3, 3))[source]

Find x-ray pixels in TEM images and replace the xray intensity with a local median value.

Parameters:
  • imageOriginal (ndarray) – The image as a 2D ndarray.
  • threshold (float) – Value above which to consider a pixel intensity an x-ray.
  • size_median_filter (list or tuple) – The footprint of the median filter over which to search for a median value.
Returns:

The image with x-ray values replaced with median local value.

Return type:

ndarray

ncempy.algo.peak_find.tripleRoll(vol, vec)[source]

Roll a 3D ndarray (e.g. an volume) along each dimension.

Parameters:
  • vol (np.ndarray) – The 3D array to roll
  • vec (iterable) – The number of time to apply the roll the tuple has the structure (axis0, axis1, axis2). Values in tuple must be int.
Returns:

The rolled volume.

Return type:

ndarray

ncempy.algo.peak_find.writeXYZ(filename, XYZ, element, comment)[source]

Write out a set of XYZ coordinates that can be read by various crystal viewing software such as Vesta.

todo: Move this to io

Parameters:
  • filename (str) – The name of the file to write to.
  • XYZ (np.ndarray) – Atom coordinates in a 2D ndarray with axes [atom number, position]. Should be shape [N, 3]
  • element (tuple) – Element symbol of each coordinate as a string (e.g. ‘C’). Should be length N.
  • comment (str) – A comment to add to the file.

Example

# Write out the peak positions as Carbon atoms: >> writeXYZ(‘name.xyz’, positions, [‘C’,] * len(positions), ‘writeXYZ example carbon atoms’)

ncempy.command_line package
Submodules
ncempy.command_line.ncem2png module

This scripts converts SER, DM3, and DM4 files into PNG

ncempy.command_line.ncem2png.dm_to_png(source_file, dest_file, fixed_dimensions=None)[source]

Saves the DM3 or DM4 source_file as PNG dest_file. If the data has three of four dimensions. The image taken is from the middle image in those dimensions.

ncempy.command_line.ncem2png.extract_dimension(img, fixed_dimensions=None)[source]
ncempy.command_line.ncem2png.main()[source]
ncempy.command_line.ncem2png.ser_to_png(source_file, dest_file)[source]

Saves the SER source_file as PNG dest_file.

Module contents
ncempy.edstomo package
Submodules
ncempy.edstomo.CharacteristicEmission module
ncempy.edstomo.CharacteristicEmission.ElamPrint(s)
ncempy.edstomo.CharacteristicEmission.GetElamFluorescenceLines(ElementName, E=None, I=None)[source]
ncempy.edstomo.CharacteristicEmission.GetFluorescenceLineEnergy(ElementName, Series='K', Line=None)[source]
ncempy.edstomo.CharacteristicEmission.GetWeightedSum(Line, Series)[source]
ncempy.edstomo.DoGenfire module
ncempy.edstomo.bruker module
ncempy.edstomo.bruker.ExtractRawSignalsFromBrukerSequence(InputDirectory=None, OutputEMD=None)[source]

Read in a set of Bruker bcf files containing EDS acquisitions and write an EMD file with the data.

Parameters:
  • InputDirectory (str) – Name of directory containing the bcf files.
  • OutputEMD (str) – Name of the output file (may include a path).
Returns:

None.

ncempy.edstomo.bruker.GetEnergyDimension(Dim)[source]

Return an energy value in eV.

ncempy.edstomo.bruker.GetSpatialDimension(Dim)[source]

Return a spatial value in meters.

ncempy.edstomo.bruker.GetTiltsFromBrukerSequence(Directory=None)[source]

Find the set of Bruker bcf files and figure out their tilts based on filenames.

There should be a set of bcf files in a directory with names: … -10.bcf, -5.bcf, 0.bcf, 5.bcf … The number represents the stage tilt of that stack. The angles do not matter but the names must be coded as above. They will be automatically read and sorted.

Parameters:Directory (str) – Name of directory containing the bcf files.
Returns:A list of tilt angles.
Return type:(list of floats)
ncempy.edstomo.postprocess module
ncempy.edstomo.preprocess module
Module contents
ncempy.eval package
Submodules
ncempy.eval.line_profile module
ncempy.eval.line_profile.line_profile(im0, p0, p1, num_points, width=0, step=0.5)[source]

Use interpolation to measure a line scan across a 2D image. map_coordinates uses a different convention from numpy indexing. So the x and y positions need to be flipped.

Parameters:
  • im0 (ndarray) – The 2D array of the image
  • p0 (tuple) – The start of the line scan (x0,y0) and (row,col)
  • p1 (tuple) – The end of the line scan (x1,y1)
  • num_points (int) – The number of points in the line scan
  • width (int) – The width of the line scan to average over. Must be an integer.
  • step (float) – The step size laterally for the interpolation. Must be < 1
Returns:

A tuple containing the line profile and the positions where the interpolation was performed

Return type:

tuple

Note

As an example: if width = 1 and step = 0.5 then there will be 3 line scans averaged.

Example

>> line, (xx, yy) = line_profile(image,(0,100),(175,100),50,step=0.5,width=5) >> fg, ax = plt.subplots(1,2) >> ax[0].plot(line,’*-‘) >> ax[1].imshow(image) >> ax[1].scatter(xx, yy)

ncempy.eval.ring_diff module

Module to process ring diffraction patterns.

ncempy.eval.ring_diff.cur_eva_vers = 'ring_diffraction_evaluation_vers0.1'

Used to identify evaluation hdf5 groups.

Type:str
ncempy.eval.ring_diff.cur_set_vers = 'ring_diffraction_setting_vers0.1'

Used to identify settings hdf5 groups.

Type:str
ncempy.eval.ring_diff.dummy_settings = {'back_init': (1.0, 1.0, 1.0), 'back_xs': (1.0, 2.0, 3.0), 'back_xswidth': 0.1, 'fit_funcs': ('voigt',), 'fit_init': (1.0, 1.0, 1.0, 1.0), 'fit_maxfev': 10, 'fit_rrange': (1.0, 2.0), 'lmax_cinit': (1, 1), 'lmax_r': 1, 'lmax_range': (1.0, 2.0), 'lmax_thresh': 1, 'mask': array([[1., 1., 1., 1., 1.], [1., 1., 1., 1., 1.], [1., 1., 1., 1., 1.], [1., 1., 1., 1., 1.], [1., 1., 1., 1., 1.]]), 'ns': (1,), 'plt_imgminmax': (0.0, 1.0), 'rad_dr': 0.1, 'rad_rmax': 1, 'rad_sigma': 0.1}

Dummy settings with all parameters set.

Type:dict
ncempy.eval.ring_diff.get_settings(parent)[source]

Get settings for radial profile evaluation.

Parameters:parent (h5py.Group) – Input group.
Returns:Settings read from parent.
Return type:dict
ncempy.eval.ring_diff.min_dummy_settings = {'back_init': (1.0, 1.0, 1.0), 'back_xs': (1.0, 2.0, 3.0), 'back_xswidth': 0.1, 'fit_funcs': ('voigt',), 'fit_init': (1.0, 1.0, 1.0, 1.0), 'fit_maxfev': None, 'fit_rrange': (1.0, 2.0), 'lmax_cinit': (1, 1), 'lmax_r': 1, 'lmax_range': (1.0, 2.0), 'lmax_thresh': 1, 'mask': None, 'ns': (1,), 'plt_imgminmax': None, 'rad_dr': None, 'rad_rmax': None, 'rad_sigma': None}

Dummy settings with all parameters set but all optional ones as Nones.

Type:dict
ncempy.eval.ring_diff.put_settings(parent, settings)[source]

Put settings for radial profile evaluation.

Creates a subgroup in parent holding the settings as attributes.

Parameters:
  • parent (h5py.Group) – Group to hold settings subgroup.
  • setting (dict) – Settings to write.
Returns:

Handle to settings group.

Return type:

h5py.Group

ncempy.eval.ring_diff.put_sglgroup(parent, label, data_grp)[source]

Adds the evaluation into parent.

Remember that the resource of the external link must not be already opened elsewhere to access data.

Parameters:
  • parent (h5py.Group) – Hdf5 group to add this evaluation group to.
  • label (str) – Label for the evaluation group.
  • data_grp (h5py.Group) – Emdtype group where to find the data.
Returns:

Handle to group.

Return type:

h5py.Group

ncempy.eval.ring_diff.run_all(parent, outfile, overwrite=False, verbose=False, showplots=False)[source]

Run on a set-up emd file to do evaluations and save results.

All evaluations within parent are run.

Parameters:
  • parent (h5py.Group) – Handle to parent.
  • outfile (ncempy.io.emd.fileEMD) – Emdfile for output.
  • overwrite (bool) – Set to overwrite existing results in outfile.
  • verbose (bool) – Set to get verbose output during run.
  • showplots (bool) – Set to directly show plots interactively.
ncempy.eval.ring_diff.run_sglgroup(group, outfile, overwrite=False, verbose=False, showplots=False)[source]

Run evaluation on a single group.

Parameters:
  • group (h5py.Group) – Handle to evaluation group to execute.
  • outfile (ncempy.io.emd.fileEMD) – Emdfile for output.
  • overwrite (bool) – Set to overwrite existing results in outfile.
  • verbose (bool) – Set to get verbose output during run.
  • showplots (bool) – Set to directly show plots interactively.
ncempy.eval.stack_align module
ncempy.eval.stack_align.stack_align(stack, align_type='static', **kwargs)[source]

Align a series of images by correlation using multicorr. All images are aligned to the start image or as dynamically starting with the start image and then each successive image. keyword arguments are passed to multicorr. Shifting is done in Fourier space which is very accurate, but wraps edges.

Parameters:
  • stack (ndarray, 3D,) – Stack of images to align. Shape [num, Y, X]
  • align_type (str) – static or dynamic alignment. Static aligns all images to the first image. Dynamic aligns each image to the previous image starting with the first image
Returns:

A tuple containing the aligned images as a 3D ndarray of shape [num, Y, X] and shifts as a 2D ndarray of shape [num, 2]

Return type:

tuple, aligned stack, shifts

ncempy.io package
Submodules
ncempy.io.read module
param filename:The path and name of the file to attempt to load. This chooses the Reader() function based on the file suffix.
type filename:str or pathlib.Path
param dsetNum:The data set number to load if there are multiple data sets in a file. This is implemented for EMD and DM files only.
type dsetNum:int
ncempy.io.emd module

This module provides an interface to the Berkeley EMD file format.

See https://emdatasets.com/ for more details.

  • It fully supports v0.2 for reading and writing

Note

General users:
Use the simplified emd.emdReader() function to load the data and meta data as a python dictionary.
Advanced users and developers:
Access the file internals through the emd.fileDEMD() class.
ncempy.io.emd.defaultDims(data, pixel_size=None)[source]

A helper function that can generate a properly setup dim tuple with default values to allow quick writing of EMD files without the need to create these dim vectors manually.

Parameters:
  • data (ndarray) – The data that will be written to the EMD file. This is used to get the number of dims and their shape
  • pixel_size (tuple, optional) – A tuple of pixel sizes. Must have same length as the number of dimensions of data.
Returns:

dims – A properly formatted tuple of dim vectors used as input to emd.emdFile.put_emdgroup()

Return type:

list

ncempy.io.emd.emdReader(filename, dsetNum=0)[source]

A simple helper function to read in the data and metadata in a structured format similar to the other ncempy readers.

Note

Note fully implemented yet. Work in progress.

Parameters:
  • filename (str or pathlib.Path) – The path to the file as a string.
  • dsetNum (int) – The index of the data set to load.
Returns:

Data and metadata as a dictionary similar to other ncempy readers.

Return type:

dict

Example

Simply load all data and metadata from a data set in an EMD Velox file
>> import ncempy.io as nio >> emd0 = nio.emd.emdReader(‘filename.emd’, dsetNum = 0)
ncempy.io.emd.emdWriter(filename, data, pixel_size=None, overwrite=False)[source]

Simple method to write data to a file formatted as an EMD v0.2. The only possible metadata to write is the pixel size for each dimension. Use the emd.fileEMD() class for more complex operations. The file must not already exist.

Parameters:
  • filename (str or pathlib.Path) – The path and file name to write the data to.
  • data (ndarray) – The data as an ndarray to write to the file.
  • pixel_size (tuple) – A tuple with the same length as the number of dims in data. len(pixel_size) == data.ndim
  • overwrite (boolean) – If file exists, overwrite it.
class ncempy.io.emd.fileEMD(filename, readonly=True)[source]

Class to represent Berkeley EMD files.

Implemented for spec 0.2 using the recommended layout for metadata.

file_name

The name of the file

Type:str
file_path

A pathlib.Path object for the open file

Type:pathlib.Path
file_hdl

The h5py file handle which provides direct access to the underlying h5py file structure.

Type:h5py.File
version

A 2-tuple of the major and minor EMD version as ints.

Type:tuple
data

A link to the data group in the EMD file.

Type:h5py.Group
microscope

A link to the microscope EMD group with metadata.

Type:h5py.Group
sample

A link to the sample EMD group with metadata.

Type:h5py.Group
user

A link to the user EMD group with metadata.

Type:h5py.Group
comments

A link to the comments in this EMD file. A time stamp is included with the comment.

Type:h5py.Group
list_emds

A list of the EMD groups in the file.

Type:list

Examples

Open an Berkeley EMD file using the advanced interface. See the emdReader function for a more convenient way to load the data. We show how to list the available data sets, load a 3D data set and plot the first image. >> from ncempy.io import emd >> with emd.fileEMD(‘filename.emd’) as emd1: >> [print(dataGroup.name) for dataGroup in emd1.list_emds] # print all available EMD datasets >> data1, dims1 = emd1.get_emdgroup(0) # load the first full data array and dimension information

find_emdgroups(parent)[source]

Find all emd_data_type groups within the group parent and return a list of references to their HDF5 groups.

Parameters:parent (h5py.Group) – Handle to the parent group.
Returns:A list of h5py.Group handles to children groups being emd_data_type groups.
Return type:list
get_emddims(group)[source]

Get the emdtype dimensions saved in in group.

Parameters:group (h5py.Group) – Reference to the emdtype HDF5 group.
Returns:List of dimension vectors plus labels and units.
Return type:tuple
get_emdgroup(group)[source]

Get the emd data saved in the requested group.

Note

The memmap keyword has been removed. Please use get_memmap().

Parameters:group (h5py._hl.group.Group or int) – Reference to the HDF5 group to load. If int is used then the item corresponding to self.list_emds is loaded
Returns:
None or tuple containing:
: np.ndarray
The data of the emdtype group.
: list
List of [0] dimension vectors, [1] labels and [2] units.
Return type:tuple or None
get_memmap(group)[source]

Opens a new fileEMD object and returns the requested EMD group. This prevents the file from being closed. The original HDF5 file is left open. The file pointer is lost, but it will be closed once the EMD group is deleted or removed. See also get_emdgroup().

Parameters:group (h5py._hl.group.Group or int) – Reference to the HDF5 group to load. If int is used then the item corresponding to self.list_emds is loaded
put_comment(msg, timestamp=None)[source]

Create a comment in the EMD file.

If timestamp already exists, the msg is appended to existing comment.

Parameters:
  • msg (str) – String of the message to save.
  • timestamp (str/None) – Timestamp used as the key, defaults to the current UTC time.
put_emdgroup(label, data, dims, parent=None, overwrite=False, **kwargs)[source]

Put an emdtype dataset into the EMD file.

Parameters:
  • label (str) – Label for the emdtype group containing the dataset.
  • data (np.ndarray) – Numpy array containing the data.
  • dims (tuple) – Tuple containing the necessary dims as ((vec, name, units), (vec, name, units), …)
  • parent (h5py.Group or None) – Parent for the emdtype group, if None it will be written to /data.
  • overwrite (bool) – Set to force overwriting entry in EMD file.
  • **kwargs (various) – Keyword arguments to be passed to h5py.create_dataset(), e.g. for compression.
Returns:

Group referencing this emdtype dataset or None if failed.

Return type:

h5py.Group or None

write_dim(label, dim, parent)[source]

Auxiliary function to write a dim dataset to parent.

Input is not checked for sanity, so handle exceptions in call.

Parameters:
  • label (str) – Label for dataset, usually dim1, dim2, dimN.
  • dim (tuple) – Tuple containing (data, name, units).
  • parent (h5py.Group) – HDF5 handle to parent group.
Returns:

HDF5 dataset handle referencing this dim.

Return type:

h5py.Group

ncempy.io.dm module

A module to load data and meta data from DM3 and DM4 files into python as written by Gatan’s Digital Micrograph program.

Note

General users:
Use the simplified dm.dmReader() function to load the data and meta data as a python dictionary.
Advanced users and developers:
Access the file internals through the dm.fileDM() class.
On Memory mode:
The fileDM support and “on memory” mode that pre-loads the file data in memory and read operations during header parsing are performed against memory. This can significantly improve performance when the file resides in a parallel file system (PFS) because latency of seek operations PFSs is very high.
ncempy.io.dm.dmReader(filename, dSetNum=0, verbose=False, on_memory=True)[source]

A simple function to parse the file and read the requested dataset. Most users will want to use this function to simplify reading data directly into memory and to retriece the spatial axes (i.e. energy axis).

Parameters:
  • filename (str) – The filename to open and read into memory
  • dSetNum (int) – The number of the data set to read. Almost always should be = 0. Default = 0
  • verbose (bool) – Allow extra printing to see file internals. Default = False
  • on_memory (bool, default True) – Whether to use the on_memory option of fileDM. Usually provides much faster data reading.

Notes

Use the coords key for spatial axes (i.e. energy loss for spectra).

Returns:A dictionary of keys where the data is in the ‘data’ key. Other metadata is contained in other named keys such as ‘pixelSize’. ‘coords’ contains the coordinate axes with the proper origin and scale applied (i.e. the energy loss axis for EELS data)
Return type:dict

Example

Load data from a single image dm3 file and display it >> from ncempy.io import dm >> im0 = dm.dmReader(‘filename.dm3’) >> plt.imshow(im0[‘data’]) #show the single image from the data file

class ncempy.io.dm.fileDM(filename, verbose=False, on_memory=True)[source]

Opens the file and reads in the header. Data is loaded using the getDataset method.

xSize, ySize, zSize

The shape of the data for each data set in the file. Each value is the same as the shape attribute for an ndarray.

Type:list
zSize2

The shape of the 4th dimension for a 4D file (i.e. 4D-STEM data)

Type:list
dataSize

The total number of bytes in each dataset. Similar to numpy’s nbytes attribute for an ndarray.

Type:list
dataOffset

The starting byte number for each dataset. This can provide fast access directly to the data if needed by seeking to this byte number.

Type:list
dataShape

The total number of dimensions in eahc dataset. Similar to numpy’s ndims attribute for an ndarray.

Type:list
file_name

The name of the file

Type:str
file_path

A pathlib.Path object for the open file

Type:pathlib.Path
fid

The file handle.

Type:file
numObjects

The number of datasets in the file. Often (but not always) the file contains a thumbnail and the raw data. The thumbnail will always be the first object. See thumbnail attribute.

Type:list
thumbnail

Tells whether the first object or dataset in the file is a thumbnail. If true then this object is always skipped in methods and it is assumed that the user wants to skip this. Can retrive this thumbnail using a built-in method if desired.

Type:bool
scale

The real size of the pixel. Real and reciprical space are supported.

Type:list
scaleUnit

The unit name as a string for each dimension of each dataset.

Type:list
origin

The origin for the real or reciprocal space scaling for each dimension. Be careful, this value is actually meant to be scaled by the scale before being used. See ncempy.io.dmReader() for proper handling of this especially for spectroscopy data.

Type:list
allTags

Contains all tags in the DM file as key value pairs.

Type:dictionary

Examples

Read data from a file containing a single image into memory >> from ncempy.io import dm >> with dm.fileDM(‘filename.dm4’) as dmFile1: >> dataSet = dmFile1.getDataset(0)

Example of reading a full multi-image DM3 file into memory: >> with dm.fileDM(‘imageSeries.dm3’)as dmFile2: >> series = dmFile2.getDataset(0)

fromfile(*args, **kwargs)[source]

Reads data from a file or memory map. Calls np.fromfile and np.frombuffer depending on the on_memory mode of the fileDM.

Note

This is essentially a passthrough function to Numpy’s frombuffer and fromfile depending on the class variable on_memory.

Parameters:
  • *args – dtype and count are required
  • **kwargs
Returns:

Data read from the file as a 1d ndarray.

Return type:

ndarray

getDataset(index)[source]

Retrieve a dataset from the DM file.

Notes

Most DM3 and DM4 files contain a small “thumbnail” as the first dataset written as RGB data. This function ignores that dataset if it exists. To retrieve the thumbnail use the getThumbnail() function.

The pixelOrigin returned is not actually the start of the coordinates. The start of the energy axis for EELS (for example) will be pixelSize * pixelOrigin. dmReader() returns the correct coordinates. The correct origin is: pixelSize * pixelOrigin and be careful about the sign as it seems some datasets might use -pixelOrigin in the previous equation.

Parameters:index (int) – The number of the data set to retrieve ignoring the thumbnail. If a thumbnail exists then index = 0 actually corresponds to the second data set in a DM file.
Returns:A dictionary of the data and meta data. The data is associated with the ‘data’ key in the dictionary.
Return type:dict
getMemmap(index)[source]

Return a numpy memmap object (read-only) for the dataset requested. This is very useful for very large datasets to avoid loading the entire data set into memory. No meta data is returned.

Parameters:index (int) – The number of the dataset in the DM file.
Returns:A read-only numpy memmap object with access to the data. The file will continue to be open as long as the memmap is open. Delete the memmap to close the file.
Return type:numpy.memmap
getMetadata(index)[source]

Get the useful metadata in the file. This parses the allTags dictionary and retrieves only the useful information about hte experimental parameters. This is a (useful) subset of the information contains in the allTags attribute.

Note: some DM files contain extra information called the Tecnai Microscope Info. This is added to the metadata dictionary as a string.

Parameters:index (int) – The number of the dataset to get the metadata from.
getSlice(index, sliceZ, sliceZ2=0)[source]

Retrieve a slice of a dataset from the DM file. The data set will have a shape according to 3D = [sliceZ,Y,X] or 4D: [sliceZ2,sliceZ,Y,X]

Note: Most DM3 and DM4 files contain a small “thumbnail” as the first dataset written as RGB data. This function ignores that dataset if it exists. To retrieve the thumbnail use the getThumbnail() function.

Warning: DM4 files with 4D data sets are written as [X,Y,Z1,Z2]. This code currently gets the [X,Y] slice. Getting the [Z1,Z2] slice is not yet implemented. Use the getMemmap() function to retrieve arbitrary slices of large data sets.

Parameters:
  • index (int) – The number of the dataset in the DM file.
  • sliceZ (int) – The slice to get along the first dimension (C-ordering) for 3D datasets or 4D datasets.
  • sliceZ2 (int) – For 4D dataset
Returns:

A dictionary containing meta data and the data.

Return type:

dict

getThumbnail()[source]

Read the thumbnail saved as the first dataset in the DM file as an RGB array. This is not fully tested. Be careful using this.

Returns:Numpy array of size [3,Y,X] which is an RGB thumbnail.
Return type:ndarray
parseHeader()[source]

Parse the header by reading the root tag group. This ensures the file pointer is in the correct place.

seek(fid, offset, from_what=0)[source]

Positions the reading head for fid. fid can be a file or memory map. Follows the same convention as file.seek

Parameters:
  • fid (file id) – File or memory map.
  • offset (int) – Number of bytes to move the head forward (positive value) or backwards (negative value).
  • from_what (int) – Reference point to use in the head movement. 0: for beginning of the file (default behavior), 1: from the current head position, and 2: from the end of the file.
tell()[source]

Return the current position in the file. Switches mode based on on_memory mode.

Returns:pos – The current position in the file.
Return type:int
writeTags(new_folder_path_for_tags=None)[source]

Write out all tags as human readable text to a text file in the same directory (or a user definable directory) and with a the same name as the DM file.

Parameters:new_folder_path_for_tags (str or pathlib.Path, Optional) – Allow user to define a different path than the directory of the current file.
ncempy.io.emdVelox module

Provides an interface to Velox EMD datasets. Not to be confused with Berkeley EMD data sets (see emd.py) instead.

The reader for EMD Berkeley and Velox files will be combined in the near future once they are fully tested separately.

Currently limited to only images. This file can not load spectra.

Note

General users:
Use the simplified emdVelox.emdVeloxReader() function to load the data and meta data as a python dictionary.
Advanced users and developers:
Access the file internals through the emd.fileEMDVelox() class.
ncempy.io.emdVelox.emdVeloxReader(filename, dsetNum=0)[source]

A simple helper function to read in the data and metadata in a structured format similar to the other ncempy readers.

Note

Not fully implemented yet. Work in progress. Important metadata is missing, but you can get the data.

Parameters:
  • filename (str or pathlib.Path) – The path to the file.
  • dsetNum (int, default = 0) – The index of the data set to load.
Returns:

Data and metadata as a dictionary similar to other ncempy readers.

Return type:

dict

Example

Load all data and metadata from a data set in an EMD file
>> import ncempy.io as nio >> emd0 = nio.emdVelox.emdVeloxReader(‘filename.emd’, dsetNum = 0)
class ncempy.io.emdVelox.fileEMDVelox(filename)[source]

Class to represent Velox EMD files. It uses the h5py caching functionality to increase the default cache size from 1MB to 10MB. This significantly improves file reading for EMDVelox files which are written with Fortran- style ordering and an inefficient choice of chunking.

list_data

A list containing each h5py data group that can be loaded.

Type:list
_file_hdl

The File handle from h5py.File.

Type:h5py.File
metaDataJSON

The full metadata for the most recently loaded data set. Note that you have to load a data set for this to be populated or run parseMetaData(num).

Type:dict
file_name

The name of the file

Type:str
file_path

A pathlib.Path object for the open file

Type:pathlib.Path

Examples

Open an EMD Velox file containing 1 image. >> import ncempy.io as nio >> with nio.emdVelox.fileEMDVelox(‘1435 1.2 Mx STEM HAADF-DF4-DF2-BF.emd’) as emd1: >> print(emd1) # print information about the file >> im0, metadata0 = emd1.get_dataset(0)

get_dataset(group, memmap=False)[source]

Get the data from a group and the associated metadata.

Parameters:
  • group (HDF5 dataset or int) – The link to the HDF5 dataset in the file or an integer for the number of the dataset. The list of datasets is held in the list_data attribute populated on class init.
  • memmap (bool, default = False) – If False (default), then a numpy ndarray is returned. If True the HDF5 data set object is returned and data is loaded from disk as needed.
Returns:

A tuple containing the data as a ndarray or a HDF5 dataset object. The second argument is a python dict of metadata.

Return type:

tuple (ndarray or HDF5 dataset, dict)

parseMetaData(group)[source]

Parse metadata in a data group. Determines the pixelSize and detector name. The EMDVelox data sets have extensive metadata stored as a JSON type string.

Parameters:group (h5py.Group or int) – The h5py group to load the metadata from which is easily retrived from the list_data attribute. If input is an int then the group corresponding to list_data attribute is used. The string metadata is loaded and parsed by the json module into a dictionary.
Returns:md – The JSON information in the file returned as a python dictionary.
Return type:dict
ncempy.io.mrc module

A module to read MRC files in python and numpy. Written according to MRC specification at http://bio3d.colorado.edu/imod/betaDoc/mrc_format.txt Also works with FEI MRC files which include a special header block with experimental information.

Note

General users:
Use the simplified mrc.mrcReader() function to load the data and meta data as a python dictionary.
Advanced users and developers:
Access the file internals through the mrc.fileMRC() class.
ncempy.io.mrc.appendData(filename, data)[source]

Append a binary set of data to the end of a MRC file. This should only be used in conjunction with writeHeader() above.

Parameters:
  • filename (str) – Name of the MRC file with pre-initiated header and some data already written.
  • data (ndarray) – Data to append to the file.
ncempy.io.mrc.emd2mrc(filename, dsetPath)[source]

Convert EMD data set into MRC data set. The final data type is float32 for convenience.

Parameters:
  • filename (str) – The name of the EMD file
  • dsetPath (str) – The HDF5 path to the top group holding the data. ex. ‘/data/raw/’
class ncempy.io.mrc.fileMRC(filename, verbose=False)[source]

Read in the data in MRC format and other useful information like metadata. Follows the specification published at http://bio3d.colorado.edu/imod/betaDoc/mrc_format.txt

file_name

The name of the file

Type:str
file_path

A pathlib.Path object for the open file

Type:pathlib.Path
fid

The file handle to the opened MRC file.

Type:file
mrcType

The internal MRC data type.

Type:int
dataType

The numpy dtype corresponding to the mrcType.

Type:np.dtype
dataSize

The number of pixels along each dimension. Corresponds to the shape attribute of a np.ndarray

Type:np.ndarray
gridSize

The size of the grid. Usually the same as dataSize

Type:np.ndarray
volumeSize

The size of the volume along each direction in Angstroms.

Type:np.ndarray
voxelSize

The size of the voxel along each direction in Angstroms.

Type:np.ndarray
cellAngles

The angles of the cell. Ignored in most cases including ncempy.

Type:np.ndarray
axisOrientations

Mapping the orientations of the data to real space directions X, Y, Z. Ignored by ncempy

Type:np.ndarray
minMaxMean

The minimum, maximum and mean value of the data to avoid computing this every time.

Type:np.ndarray
extra

Extra mbinary metadata if it exists.

Type:np.ndarray
FEIinfo

A dictionary of metadata used by FEI (Thermo Fischer) microsocpes for important metadata. This metadata overwrites the voxelsize attribute if it exists.

Type:dict
dataOffset

The integer offset in bytes to the start of the raw data.

Type:int
dataOut

Will hold the data and metadata to output to the user after getDataset() call.

Type:dict
v

More output for debugging. False by default

Type:bool

Examples

Read in all data and metadata into memory. >> import ncempy.io as nio >> mrc0 = nio.mrc.mrcReader(‘file.mrc’)

Low level operations to get 1 slice of the 3D data >> import ncempy.io as nio >> with nio.mrc.fileMRC(‘file.mrc’) as f1: >> single_slice = f1.getSlice(0)

getDataset()[source]

Read in the full data block and reshape to an ndarray with C-style ordering.

getMemmap()[source]

Return a numpy memmap object (read-only) for the dataset. This is very useful for very large datasets to avoid loading the entire data set into memory. No meta data is returned.

Returns:A read-only numpy memmap object with access to the data on disk.
Return type:numpy.core.memmap
getSlice(num)[source]

Read in a slice of an MRC file. Useful for parsing through a large file without reading the entire data set into memory.

Parameters:num (int) – Get the requested image.
Returns:out – A 2D slice or a 3D set of slices along the first index
Return type:ndarray
Raises:IndexError – If num > the number of slices.
parseHeader()[source]

Read the header information which includes data type, data size, data shape, and metadata.

Note

This header uses Fortran-style ordering. Numpy uses C-style ordering. The header is read in and then some attributes are reversed [::-1] at the end for output to the user to match C-ordering in numpy.

Todo

Implement special dtype to read the entire header at once. Read everything at once using special dtype. ~5x faster than multiple np.fromfile() reads:

Untested but works in theory headerDtype = np.dtype([(‘head1’,’10int32’),(‘head2’,’6float32’),(‘axisOrientations’,’3int32’), (‘minMaxMean’,’3int32’),(‘extra’,’32int32’)]) head = np.fromfile(self.fid,dtype=headerDtype,count=1)

ncempy.io.mrc.mrc2emd(file_name)[source]

Write an MRC file as an HDF5 file in EMD format with same file name and .emd ending. Header information is retained as attributes.

Parameters:file_name (str) – The name of the file to convert from MRC to EMD format.
Returns:out – 1 if successful.
Return type:int

Todo

Update this to use ncempy.emd class

ncempy.io.mrc.mrc2raw(file_name)[source]

Convert an MRC type file to raw binary. The entire data is read into memory and then written to a new raw file in the same location with the data type and shape (C-ordering) written into the filename. No other meta data is retained.

Parameters:file_name (str or pathlib.Path) – The file name to convert.
ncempy.io.mrc.mrcReader(file_name)[source]

A simple function to read open a MRC, parse the header, and read the full data set.

Parameters:file_name (str or pathlib.Path) – The path to the file to load.
Returns:out – A dictionary containing the data and interesting metadata. The data is attached to the ‘data’ key.
Return type:dict

Example

Read in all data from disk into memory. This assumes the dataset is 3 dimensional: >> from ncempy.io.mrc import mrcReader >> import matplotlib.pyplot as plt >> mrc1 = mrcReader(‘filename.mrc’) >> plt.imshow(mrc1[‘data’][0, :, :]) # show the first image in the data set

ncempy.io.mrc.mrcWriter(filename, data, pixelSize, forceWrite=False)[source]

Write out a MRC type file according to the specification at http://bio3d.colorado.edu/imod/doc/mrc_format.txt

Parameters:
  • filename (str or pathlib.Path) – The name or Path of the file to write out to.
  • data (ndarray) – The array data to write to disk.
  • pixelSize (tuple) – The size of the pixel along each direction (in Angstroms) as a 3 element vector (sizeZ,sizeY,sizeX).
  • forceWrite (bool) – This will write the data as a C-contiguous array. It is not suggested to use this option and it will be removed in future versions.
ncempy.io.mrc.writeHeader(filename, shape, dtype, pixelSize)[source]

Write out a MRC type file header according to the specification at http://bio3d.colorado.edu/imod/doc/mrc_format.txt. This is useful for initializing an MRC file and then writing to it manually or see appendData() function below.

Parameters:
  • filename (str) – The name of the EMD file
  • shape (tuple) – The shape of the data to write
  • dtype (numpy.dtype) – The dtype to write out the data as. Only some numpy dtypes are supported byt his format. It is suggested to use np.float32 in most cases for maximum compatibility.
  • pixelSize (tuple) – The size of the pixel along each direction (in Angstroms) as a 3 element vector (sizeX,sizeY,sizeZ). sizeZ could be the angular step for a tilt series
ncempy.io.ser module

This module provides an interface to the SER file format written by FEI and Therm0 Fischer’s program TIA.

It reads STEM and TEM images and other datasets.

It is based on information provided by Dr Chris Boothroyd and work by Peter Ercius’ original serReader code in Matlab.

Note

General users:
Use the simplified ser.serReader() function to load the data and meta data as a python dictionary.
Advanced users and developers:
Access the file internals through the ser.fileSER() class.
exception ncempy.io.ser.NotSERError[source]

Exception if a file is not in SER file format.

class ncempy.io.ser.fileSER(filename, verbose=False)[source]

Class to represent SER files (read only).

_file_hdl

The open file as a raw stream.

Type:file
_emi

A dictionary of metadata from the EMI file accompanying the SER file.

Type:dict
file_name

The name of the file

Type:str
file_path

A pathlib.Path object for the open file

Type:pathlib.Path
head

Header information for the SER file as a dictionary. Provides direct access to data offsets and other internal file information.

Type:dict

Note

For most users, we suggest using the ser.serReader() function to load the full data set into memory. Otherwise, this class provides low level access to the SER file data and metadata and internals.

Examples

Read data from a single image into memory using the low level API. >> import matplotlib.pyplot as plt >> import ncempy.io as nio >> with nio.ser.fileSER(‘filename.ser’) as ser1: >> data, metadata = ser1.getDataset(0)

SER files are internally structured such that each image in a series is a different data set. Thus, time series data should be read as the following: >> with ser.fileSER(‘filename_1.ser’) as ser1: >> image0, metadata0 = ser1.getDataset(0) >> image1, metadata1 = ser1.getDataset(1)

getDataset(index, verbose=False)[source]

Retrieve data and meta data for one image or spectra from the file.

Parameters:
  • index (int) – Index of dataset.
  • verbose (bool, optional) – True to get extensive output while reading the file.
Returns:

dataset – Tuple contains data as np.ndarray and metadata (pixel size, etc.) as a dict.

Return type:

tuple, 2 elements in form (data metadata)

readHeader(verbose=False)[source]

Read and return the SER files header.

Parameters:verbose (bool) – True to get extensive output while reading the file.
Returns:The header of the SER file as dict.
Return type:dict
writeEMD(filename)[source]

Write SER data to an EMD file.

Parameters:filename (str or pathlib.Path) – Name of the EMD file.
ncempy.io.ser.read_emi(filename)[source]

Read the meta data from an emi file.

Parameters:filename (str or pathlib.Path) – Path to the emi file.
Returns:Dictionary of experimental metadata stored in the EMI file.
Return type:dict
ncempy.io.ser.serReader(filename)[source]
Simple function to parse the file and read all datasets. This is a one function implementation to load all data
in a ser file.
Parameters:filename (str) – The filename of the SER file containing the data.
Returns:dataOut – A dictionary containing the data and meta data. The data is accessed using the ‘data’ key and is a 1, 2, 3, or 4 dimensional numpy ndarray.
Return type:dict

Examples

Load a single image data set and show the image: >> import ncempy.io as nio >> ser1 = nio.ser.serReader(‘filename_1.ser’) >> plt.imshow(ser1[‘data’]) # show the single image from the data file

ncempy.viz package

A set of visualization functions based on matplotlib.pyplot which are generally useful for S/TEM data visualization.

ncempy.viz.im_and_fft(im, d=1.0, fft=None)[source]

Show the image and its fft side by side. Uses imfft to show the fft.

Parameters:
  • im (np.ndarray) – The image to show in both real and FFT space
  • d (float) – The pixel spacing
  • fft (np.ndarray, optional) – The FFT to display. If not provided then np.fft.fft2 is used.
Returns:

The matplotlib.pyplot figure

Return type:

plt.figure

ncempy.viz.im_calibrated(im, d)[source]

Plot an image calibrated using the pixel size d. The centers of the pixels will be the the center of each measurement. So, if you plot positions in real coordinates the points will be plotted in the center of the pixel.

Parameters:
  • im (np.ndarray) – The image to show using imshow
  • d (float) – The pixel size in both directions. The pixel size must be isotropic.
Returns:

The figure containing the plot

Return type:

pyplot.figure

ncempy.viz.imfft(im, d=1.0, ax=None)[source]

Show a 2D FFT as a diffractogram with log scaling applied and zero frequency fftshifted tp the center. A new figure is created or an axis can be specified.

The diffracotgram is calculated from the original intensities (I) as

\[1 + 0.001 * I ^2\]
Parameters:
  • im (np.ndarray) – The 2D fft of the diffraction pattern to display as a diffractogram
  • d (float, optional, default = 1.0) – The real space pixel size of the image used to get the FFT
  • ax (pyplot axis, optional) – An axis to plot into.
Returns:

The AxesImage that contains the image displayed.

Return type:

matplotlib.image.AxesImage

Example

This example shows how to display a 2D ndarray (image) as a diffractogram. The image has a real space pixel size of 0.1 nanometer.

>> imageFFT = np.fft.fft2(im) >> ncempy.viz.imfft(imageFFT, d = 0.1)

ncempy.viz.imrfft(im, d=1.0, ax=None)[source]

Show a 2D rFFT (real FFT) as a diffractogram with log scaling applied and fftshift-ed along axis 0. See imfft for full details.

Parameters:
  • im (ndarray) – The 2D fft of the diffraction pattern to display as a diffractogram
  • d (float, optional, default = 1.0) – The real space pixel size of the image used to get the FFT
  • ax (pyplot axis, optional) – An axis to plot into.
Returns:

The AxesImage that contains the image displayed.

Return type:

matplotlib.image.AxesImage

ncempy.viz.imsd(im, vmin=-2, vmax=2, **kwargs)[source]

Show an array as an image with intensities compared to the standard deviation of the data. Other keyword args are passed to pyplot.imshow(). vmin and vmax are by default set to -2 and 2 respectively which are usually good values to set for S/TEM data.

Parameters:
  • im (np.ndarray) – The image to show.
  • vmax (vmin,) – The vmin and vmax values to pass to imshow.
Returns:

The handle to the created figure

Return type:

matplotlib.pyplot.Figure

ncempy.viz.plot_distpolar(points, dims, dists, ns, show=False)[source]

Plot the results of distortion fitting in polar coordinates.

Parameters:
  • points (np.ndarray) – Points in polar coords.
  • dims (tuple) – Dimensions, necessary to have unit information.
  • dists (np.ndarray) – Results of dist fitting, length according to ns.
  • ns (list) – List of used orders.
  • show (bool) – Set to directly show the plot in interactive mode.
Returns:

Image of the plot.

Return type:

np.ndarray

ncempy.viz.plot_fit(r, intens, dims, funcs, param, show=False)[source]

Plot the fit results to the radial profile.

Parameters:
  • r (np.ndarray) – r-axis of radial profile.
  • intens (np.ndarray) – Intensity-axis of radial profile.
  • dims (tuple) – Dimensions of original image to read out units.
  • funcs (tuple) – List of functions.
  • param (np.ndarray) – Parameters for functions in funcs.
  • show (bool) – Set to directly show plot interactively.
Returns:

Image of the plot.

Return type:

np.ndarray

ncempy.viz.plot_points(img, points, vminmax=(0, 1), dims=None, invert=False, show=False)[source]

Plot the detected points on the input image for checking.

Parameters:
  • img (np.ndarray) – Image.
  • points (np.ndarray) – Array containing the points.
  • vminmax (tuple) – Tuple of two values for relative lower and upper cut off to display image.
  • dims (tuple) – Tuple of dims to plot in dimensions.
  • invert (bool) – Set to invert the image.
  • show (bool) – Set to directly show the plot interactively.
Returns:

Image of the plot.

Return type:

np.ndarray

ncempy.viz.plot_radialprofile(r, intens, dims, show=False)[source]

Plot radial profile.

Parameters:
  • r (np.ndarray) – r-axis of radial profile.
  • intens (np.ndarray) – Intensity-axis of radial profile.
  • dims (tuple) – Dimensions of original image to read out units.
  • show (bool) – Set to directly show plot interactively.
Returns:

Image of the plot.

Return type:

np.ndarray

ncempy.viz.plot_ringpolar(points, dims, show=False)[source]

Plot points in polar coordinate system.

Parameters:
  • points (np.ndarray) – Positions in polar coords.
  • dims (tuple) – Dimension information to plot labels.
  • show (bool) – Set to directly show plot in interactive mode.
Returns:

Image of the plot.

Return type:

numpy.ndarray

class ncempy.viz.stack_view(stack, **kwargs)[source]

Bases: object

Class to allow a volume to be scrubbed through with a matplotlib slider widget. The first axis of the volume is the slicing axis. Other keyword args are passed directly to imshow upon the figure creation.

Parameters:
  • stack (numpy.ndarray, 3D stack) – The stack of to show as images
  • **kwargs – Passed directly to pyplot.imshow()

Examples

Simple file reading

For simple file loading use the convenient read function.

Here is an example of loading a Digital Micrograph file. It returns the data and the metadata in an easily human readable form.

>>> import ncempy.io as nio
>>> dmData = nio.read('path/to/file/data.dm3')
>>> print(dmData['data'].shape) # the shape of the data
>>> print(dmData['pixelSize']) # print the pixel size

For developers

Developers will want access to the entire internal structure of the file. Then you need to instantiate the class, parse the header and then access its contents

Here is how to open a DM file and have access to the internal file parameters

>>> import ncempy.io as nio
>>> with nio.dm.fileDM('/path/to/file/data.dm3') as dm0:
>>>     dmData = dm0.getDataset(0) #the full first data set
>>>     dmSlice = dm0.getSlice(0, 1) #the second image of the first data set; equal to dmData[1,:,:] from the line above
>>>     print(dm0.metaData) # all of the interesting metadata for this data set (pixel size, accelerating voltage, etc.)
>>>     print(dm0.allTags) # all of the tags in the file including tags specific to the Digital Micrpograph software program

License

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see https://www.gnu.org/licenses/.