bcdi.utils: various utilities for DFT registration and subpixel shift, data loading and fitting, validating parameters

utilities

Various non-specific utilities for i/o, fitting curves, generic functions for cropping or padding a dataset…

API Reference

Functions related to data loading, encoding, fitting, data manipulation.

bcdi.utils.utilities.apply_logical_array(arrays: Real | numpy.ndarray | Tuple[Real | numpy.ndarray, ...], frames_logical: numpy.ndarray | None) Real | numpy.ndarray | Tuple[Real | numpy.ndarray, ...]

Apply a logical array to a sequence of arrays.

Assuming a 1D array, it will be cropped where frames_logical is 0.

Parameters:
  • arrays – a list or tuple of numbers or 1D arrays

  • frames_logical – array of length the number of measured frames. In case of cropping/padding the number of frames changes. A frame whose index is set to 1 means that it is used, 0 means not used, -1 means padded (added) frame

Returns:

an array (if a single array was provided) or a tuple of cropped arrays

bcdi.utils.utilities.bin_data(array, binning, debugging=False, **kwargs)

Rebin a 1D, 2D or 3D array.

If its dimensions are not a multiple of binning, the array will be cropped. Adapted from PyNX.

Parameters:
  • array – the array to resize

  • binning – the rebin factor - pixels will be summed by groups of binning (x binning (x binning)). This can also be a tuple/list of rebin values along each axis, e.g. binning=(4,1,2) for a 3D array

  • debugging – boolean, True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’logger’: an optional logger

Returns:

the binned array

bcdi.utils.utilities.bin_parameters(binning: int, nb_frames: int, params: List[Any], debugging: bool = True) List[Any]

Bin some parameters.

It selects parameter values taking into account an eventual binning of the data. The use case is to bin diffractometer motor positions for a dataset binned along the rocking curve axis.

Parameters:
  • binning – binning factor for the axis corresponding to the rocking curve

  • nb_frames – number of frames of the rocking curve dimension

  • params – list of parameters

  • debugging – set to True to have printed parameters

Returns:

list of parameters (same list length), taking into account binning

bcdi.utils.utilities.cast(val: float | ~typing.List | numpy.ndarray, target_type: type = <class 'float'>, **kwargs) float | List | numpy.ndarray

Cast val to a number or an array of numbers of the target type.

Parameters:
  • val – the value to be converted

  • target_type – the type to convert to

  • kwargs

    • ‘logger’: an optional logger

bcdi.utils.utilities.catch_error(exception)

Process exception in asynchronous multiprocessing.

Parameters:

exception – the arisen exception

bcdi.utils.utilities.convert_str_target(value: Any, target: str, conversion_table: Dict[str, Any] | None = None) Any

Convert strings from value to the desired target.

Parameters:
  • value – an object containing strings to be converted

  • target – the target string, which has to be present in the conversion table

  • conversion_table – a dictionary for the conversion

Returns:

the converted object

bcdi.utils.utilities.crop_pad(array, output_shape, pad_value=0, pad_start=None, crop_center=None, debugging=False, **kwargs)

Crop or pad the 3D object depending on output_shape.

Parameters:
  • array – 3D complex array to be padded

  • output_shape – desired output shape (3D)

  • pad_value – will pad using this value

  • pad_start – for padding, tuple of 3 positions in pixel where the original array should be placed. If None, padding is symmetric along the respective axis

  • crop_center – for cropping, [z, y, x] position in the original array (in pixels) of the center of the output array. If None, it will be set to the center of the original array

  • debugging (bool) – set to True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’logger’: an optional logger

Returns:

myobj cropped or padded with zeros

bcdi.utils.utilities.crop_pad_1d(array, output_length, pad_value=0, pad_start=None, crop_center=None, extrapolate=False)

Crop or pad the 2D object depending on output_shape.

Parameters:
  • array – 1D complex array to be padded

  • output_length – int desired output length

  • pad_value – will pad using this value

  • pad_start – for padding, position in pixel where the original array should be placed. If None, padding is symmetric

  • crop_center – for cropping, position in pixels in the original array of the center of the ourput array. If None, it will be set to the center of the original array

  • extrapolate – set to True to extrapolate using the current spacing (supposed constant)

Returns:

myobj cropped or padded

bcdi.utils.utilities.crop_pad_2d(array, output_shape, pad_value=0, pad_start=None, crop_center=None, debugging=False)

Crop or pad the 2D object depending on output_shape.

Parameters:
  • array – 2D complex array to be padded

  • output_shape – list of desired output shape [y, x]

  • pad_value – will pad using this value

  • pad_start – for padding, tuple of 2 positions in pixel where the original array should be placed. If None, padding is symmetric along the respective axis

  • crop_center – for cropping, [y, x] position in the original array (in pixels) of the center of the ourput array. If None, it will be set to the center of the original array

  • debugging (bool) – set to True to see plots

Returns:

myobj cropped or padded with zeros

bcdi.utils.utilities.decode_json(dct)

Define the parameter object_hook in json.load function, supporting various types.

Parameters:

dct – the input dictionary of strings

Returns:

a dictionary

bcdi.utils.utilities.dos2unix(input_file: str, savedir: str) None

Convert DOS linefeeds (crlf) to UNIX (lf).

Parameters:
  • input_file – the original filename (absolute path)

  • savedir – path of the directory where to save

bcdi.utils.utilities.find_crop_center(array_shape, crop_shape, pivot)

Find the position of the center of the cropping window.

It finds the closest voxel to pivot which allows to crop an array of array_shape to crop_shape.

Parameters:
  • array_shape (tuple) – initial shape of the array

  • crop_shape (tuple) – final shape of the array

  • pivot (tuple) – position on which the final region of interest dhould be centered (center of mass of the Bragg peak)

Returns:

the voxel position closest to pivot which allows cropping to the defined shape.

bcdi.utils.utilities.find_file(filename: str | None, default_folder: str, **kwargs) str

Locate a file.

The filename can be either the name of the file (including the extension) or the full path to the file.

Parameters:
  • filename – the name or full path to the file

  • default_folder – it will look for the file in that folder if filename is not the full path.

  • kwargs

    • ‘logger’: an optional logger

Returns:

str, the path to the file

bcdi.utils.utilities.find_nearest(reference_array, test_values, width=None)

Find the indices where original_array is nearest to array_values.

Parameters:
  • reference_array – a 1D array where to look for the nearest values

  • test_values – a number or a 1D array of numbers to be tested

  • width – if not None, it will look for the nearest element within the range [x-width/2, x+width/2[

Returns:

index or indices from original_array nearest to values of length len(test_values). Returns (index - 1). If there is no nearest neighbour in the range defined by width.

bcdi.utils.utilities.fit3d_poly1(x_axis, a, b, c, d)

Calculate the 1st order polynomial function on points in a 3D grid.

Parameters:
  • x_axis – (3xN) tuple or array of 3D coordinates

  • a – offset

  • b – 1st order parameter for the 1st coordinate

  • c – 1st order parameter for the 2nd coordinate

  • d – 1st order parameter for the 3rd coordinate

Returns:

the 1st order polynomial calculated on x_axis

bcdi.utils.utilities.fit3d_poly2(x_axis, a, b, c, d, e, f, g)

Calculate the 2nd order polynomial function on points in a 3D grid.

Parameters:
  • x_axis – (3xN) tuple or array of 3D coordinates

  • a – offset

  • b – 1st order parameter for the 1st coordinate

  • c – 1st order parameter for the 2nd coordinate

  • d – 1st order parameter for the 3rd coordinate

  • e – 2nd order parameter for the 1st coordinate

  • f – 2nd order parameter for the 2nd coordinate

  • g – 2nd order parameter for the 3rd coordinate

Returns:

the 2nd order polynomial calculated on x_axis

bcdi.utils.utilities.fit3d_poly3(x_axis, a, b, c, d, e, f, g, h, i, j)

Calculate the 3rd order polynomial function on points in a 3D grid.

Parameters:
  • x_axis – (3xN) tuple or array of 3D coordinates

  • a – offset

  • b – 1st order parameter for the 1st coordinate

  • c – 1st order parameter for the 2nd coordinate

  • d – 1st order parameter for the 3rd coordinate

  • e – 2nd order parameter for the 1st coordinate

  • f – 2nd order parameter for the 2nd coordinate

  • g – 2nd order parameter for the 3rd coordinate

  • h – 3rd order parameter for the 1st coordinate

  • i – 3th order parameter for the 2nd coordinate

  • j – 3th order parameter for the 3rd coordinate

Returns:

the 3rd order polynomial calculated on x_axis

bcdi.utils.utilities.fit3d_poly4(x_axis, a, b, c, d, e, f, g, h, i, j, k, m, n)

Calculate the 4th order polynomial function on points in a 3D grid.

Parameters:
  • x_axis – (3xN) tuple or array of 3D coordinates

  • a – offset

  • b – 1st order parameter for the 1st coordinate

  • c – 1st order parameter for the 2nd coordinate

  • d – 1st order parameter for the 3rd coordinate

  • e – 2nd order parameter for the 1st coordinate

  • f – 2nd order parameter for the 2nd coordinate

  • g – 2nd order parameter for the 3rd coordinate

  • h – 3rd order parameter for the 1st coordinate

  • i – 3th order parameter for the 2nd coordinate

  • j – 3th order parameter for the 3rd coordinate

  • k – 4th order parameter for the 1st coordinate

  • m – 4th order parameter for the 2nd coordinate

  • n – 4th order parameter for the 3rd coordinate

Returns:

the 4th order polynomial calculated on x_axis

bcdi.utils.utilities.function_lmfit(params, x_axis, distribution, iterator=0)

Calculate distribution defined by lmfit Parameters.

Parameters:
  • params – a lmfit Parameters object

  • x_axis – where to calculate the function

  • distribution – the distribution to use

  • iterator – the index of the relevant parameters

Returns:

the gaussian function calculated at x_axis positions

bcdi.utils.utilities.gaussian(x_axis, amp, cen, sig)

Gaussian line shape.

Parameters:
  • x_axis – where to calculate the function

  • amp – the amplitude of the Gaussian

  • cen – the position of the center

  • sig – HWHM of the Gaussian

Returns:

the Gaussian line shape at x_axis

bcdi.utils.utilities.gaussian_window(window_shape, sigma=0.3, mu=0.0, voxel_size=None, debugging=False)

Create a 2D or 3D Gaussian window using scipy.stats.multivariate_normal.

Parameters:
  • window_shape – shape of the window

  • sigma – float, sigma of the distribution

  • mu – float, mean of the distribution

  • voxel_size – tuple, voxel size in each dimension corresponding to window_shape. If None, it will default to 1/window_shape[ii] for each dimension so that it is independent of the shape of the window

  • debugging – True to see plots

Returns:

the Gaussian window

bcdi.utils.utilities.generate_frames_logical(nb_images: int, frames_pattern: List[int] | None) numpy.ndarray

Generate a logical array allowing ti discrad frames in the dataset.

Parameters:
  • nb_images – the number of 2D images in te dataset.

  • frames_pattern – user-provided list which can be: - a binary list of length nb_images - a list of the indices of frames to be skipped

Returns:

a binary numpy array of length nb_images

bcdi.utils.utilities.higher_primes(number, maxprime=13, required_dividers=(4,))

Find the closest larger number that meets some condition.

Find the closest integer >=n (or list/array of integers), for which the largest prime divider is <=maxprime, and has to include some dividers. The default values for maxprime is the largest integer accepted by the clFFT library for OpenCL GPU FFT. Adapted from PyNX.

Parameters:
  • number – the integer number

  • maxprime – the largest prime factor acceptable

  • required_dividers – a list of required dividers for the returned integer.

Returns:

the integer (or list/array of integers) fulfilling the requirements

bcdi.utils.utilities.image_to_ndarray(filename, convert_grey=True)

Convert an image to a numpy array using pillow.

Matplotlib only supports the PNG format.

Parameters:
  • filename – absolute path of the image to open

  • convert_grey – if True and the number of layers is 3, it will be converted to a single layer of grey

Returns:

bcdi.utils.utilities.in_range(point, extent)

Check if a point in within a certain range.

It returns a boolean depending on whether point is in the indices range defined by extent or not.

Parameters:
  • point – tuple of two real numbers (2D case) (y, x) or three real numbers (3D case) (z, y, x) representing the voxel indices to be tested

  • extent – tuple of four integers (2D case) (y_start, y_stop, x_tart, x_stop) or six integers (3D case) (z_start, z_stop, y_start, y_stop, x_tart, x_stop) representing the range of valid indices

Returns:

True if point belongs to extent, False otherwise

bcdi.utils.utilities.line(x_array, a, b)

Return y values such that y = a*x + b.

Parameters:
  • x_array – a 1D numpy array of length N

  • a – coefficient for x values

  • b – constant offset

Returns:

an array of length N containing the y values

bcdi.utils.utilities.linecut(array, point, direction, direction_basis='voxel', voxel_size=1)

Calculate iteratively a linecut through an array without interpolation.

Parameters:
  • array – 2D or 3D numpy array from which the linecut will be extracted

  • point – tuple of three integral indices, point by which the linecut pass.

  • direction – list of 2 (for 2D) or 3 (for 3D) vector components, direction of the linecut in units of pixels

  • direction_basis – ‘orthonormal’ if the vector direction is expressed in an orthonormal basis. In that case it will be corrected for the different voxel sizes in each direction. ‘voxel’ if direction is expressed in the non-orthonormal basis defined by the voxel sizes in each direction.

  • voxel_size – real positive number or tuple of 2 (for 2D) or 3 (for 3D) real positive numbers representing the voxel size in each dimension.

Returns:

distances (1D array, distance along the linecut in the unit given by voxel_size) and cut (1D array, linecut through array in direction passing by point)

bcdi.utils.utilities.load_background(background_file)

Load a background file.

Parameters:

background_file – the path of the background file

Returns:

a 2D background

bcdi.utils.utilities.load_file(file_path: str, fieldname: str | None = None) Tuple[numpy.ndarray, str]

Load a file.

In case of .cxi or .h5 file, it will use a default path to the data. ‘fieldname’ is used only for .npz files.

Parameters:
  • file_path – the path of the reconstruction to load. Format supported: .npy .npz .cxi .h5

  • fieldname – the name of the field to be loaded

Returns:

the loaded data and the extension of the file

bcdi.utils.utilities.load_flatfield(flatfield_file)

Load a flatfield file.

Parameters:

flatfield_file – the path of the flatfield file

Returns:

a 2D flatfield

bcdi.utils.utilities.load_hotpixels(hotpixels_file)

Load a hotpixels file.

Parameters:

hotpixels_file – the path of the hotpixels file

Returns:

a 2D array of hotpixels (1 for hotpixel, 0 for normal pixel)

bcdi.utils.utilities.lorentzian(x_axis, amp, cen, sig)

Lorentzian line shape.

Parameters:
  • x_axis – where to calculate the function

  • amp – the amplitude of the Lorentzian

  • cen – the position of the center

  • sig – HWHM of the Lorentzian

Returns:

the Lorentzian line shape at x_axis

bcdi.utils.utilities.make_support(arrays: numpy.ndarray | Sequence[numpy.ndarray], support_threshold: float) numpy.ndarray | Sequence[numpy.ndarray]

Create a support for each provided array, using a threshold on its modulus.

Parameters:
  • arrays – a sequence of numpy ndarrays

  • support_threshold – a float in [0, 1], normalized threshold that will be applied to the modulus of each array

Returns:

a tuple of numpy ndarrays, supports corresponding to each input array

bcdi.utils.utilities.mean_filter(data, nb_neighbours, mask=None, target_val=0, extent=1, min_count=3, interpolate='mask_isolated', debugging=False, **kwargs)

Mask or apply a mean filter to data.

The procedure is applied only if the empty pixel is surrounded by nb_neighbours or more pixels with at least min_count intensity per pixel.

Parameters:
  • data – 2D or 3D array to be filtered

  • nb_neighbours – minimum number of non-zero neighboring pixels for median filtering

  • mask – mask array of the same shape as data

  • target_val – value where to interpolate, allowed values are int>=0 or np.nan

  • extent – in pixels, extent of the averaging window from the reference pixel (extent=1, 2, 3 … corresponds to window width=3, 5, 7 … )

  • min_count – minimum intensity (inclusive) in the neighboring pixels to interpolate, the pixel will be masked otherwise.

  • interpolate – based on ‘nb_neighbours’, if ‘mask_isolated’ will mask isolated pixels, if ‘interp_isolated’ will interpolate isolated pixels

  • debugging (bool) – set to True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

Returns:

updated data and mask, number of pixels treated

bcdi.utils.utilities.move_log(result: Tuple[Path, Path, Logger | None])

Move log files to the desired location, after processing a file.

It can be used as a standard function or as a callback for multiprocessing.

Parameters:

result – the output of process_scan, containing the 2d data, 2d mask, counter for each frame; the file index; and an optional logger

bcdi.utils.utilities.objective_lmfit(params, x_axis, data, distribution)

Calculate the total residual for fits to several data sets.

Data sets should be held in a 2-D array (1 row per dataset).

Parameters:
  • params – a lmfit Parameters object

  • x_axis – where to calculate the distribution

  • data – data to fit

  • distribution – distribution to use for fitting

Returns:

the residuals of the fit of data using the parameters

bcdi.utils.utilities.pad_from_roi(arrays: numpy.ndarray | Tuple[numpy.ndarray, ...], roi: List[int], binning: Tuple[int, int], pad_value: float | Tuple[float, ...] = 0.0, **kwargs) numpy.ndarray | Tuple[numpy.ndarray, ...]

Pad a 3D stack of frames provided a region of interest.

The stacking is assumed to be on the first axis.

Parameters:
  • arrays – a 3D array or a sequence of 3D arrays of the same shape

  • roi – the desired region of interest of the unbinned frame. For an array in arrays, the shape is (nz, ny, nx), and roi corresponds to [y0, y1, x0, x1]

  • binning – tuple of two integers (binning along Y, binning along X)

  • pad_value – number or tuple of nb_arrays numbers, will pad using this value

  • kwargs

    • ‘logger’: an optional logger

Returns:

an array (if a single array was provided) or a tuple of arrays interpolated on an orthogonal grid (same length as the number of input arrays)

bcdi.utils.utilities.plane(xy_array, a, b, c)

Return z values such that z = a*x + b*y + c.

Parameters:
  • xy_array – a (2xN) numpy array, x values being the first row and y values the second row

  • a – coefficient for x values

  • b – coefficient for y values

  • c – constant offset

Returns:

an array of length N containing the z values

bcdi.utils.utilities.plane_dist(indices, params)

Calculate the distance of an ensemble of voxels to a plane given by its parameters.

Parameters:
  • indices – a (3xN) numpy array, x values being the 1st row, y values the 2nd row and z values the 3rd row

  • params – a tuple of coefficient (a, b, c, d) such that ax+by+cz+d=0

Returns:

a array of shape (N,) containing the distance to the plane for each voxel

bcdi.utils.utilities.plane_fit(indices, label='', threshold=1, debugging=False)

Fit a plane to the voxels defined by indices.

Parameters:
  • indices – a (3xN) numpy array, x values being the 1st row, y values the 2nd row and z values the 3rd row

  • label – int, label of the plane used for the title in the debugging plot

  • threshold – the fit will be considered as good if the mean distance of the voxels to the plane is smaller than this value

  • debugging – True to see plots

Returns:

a tuple of coefficient (a, b, c, d) such that ax+by+cz+d=0, the matrix of covariant values

bcdi.utils.utilities.primes(number)

Return the prime decomposition of n as a list. Adapted from PyNX.

Parameters:

number – the integer to be decomposed

Returns:

the list of prime dividers of number

bcdi.utils.utilities.pseudovoigt(x_axis, amp, cen, sig, ratio)

Pseudo Voigt line shape.

Parameters:
  • x_axis – where to calculate the function

  • amp – amplitude of the Pseudo Voigt

  • cen – position of the center of the Pseudo Voigt

  • sig – FWHM of the Pseudo Voigt

  • ratio – ratio of the Gaussian line shape

Returns:

the Pseudo Voigt line shape at x_axis

bcdi.utils.utilities.ref_count(address)

Get the reference count using ctypes module.

Parameters:

address – integer, the memory adress id

Returns:

the number of references to the memory address

bcdi.utils.utilities.remove_avg_background(array, q_values, avg_background, avg_qvalues, method='normalize')

Subtract the average 1D background to the 3D array using q values.

Parameters:
  • array – the 3D array. It should be sparse for faster calculation.

  • q_values – tuple of three 1D arrays (qx, qz, qy), q values for the 3D dataset

  • avg_background – average background data

  • avg_qvalues – q values for the 1D average background data

  • method – ‘subtract’ or ‘normalize’

Returns:

the 3D background array

bcdi.utils.utilities.remove_nan(data, mask=None)

Remove nan values from data.

Optionally, it can update a mask (masked data = 1 in the mask, 0 otherwise).

Parameters:
  • data – numpy ndarray

  • mask – if provided, numpy ndarray of the same shape as the data

Returns:

the filtered data and (optionally) mask

bcdi.utils.utilities.rgb2gray(rgb)

Convert a three layered RGB image in gray.

Parameters:

rgb – the image in RGB

Returns:

the image conveted to gray

bcdi.utils.utilities.rotate_crystal(arrays, axis_to_align=None, reference_axis=None, voxel_size=None, fill_value=0, rotation_matrix=None, is_orthogonal=False, reciprocal_space=False, debugging=False, **kwargs)

Rotate arrays to align axis_to_align onto reference_axis.

The pivot of the rotation is in the center of the arrays. axis_to_align and reference_axis should be in the order X Y Z, where Z is downstream, Y vertical and X outboard (CXI convention).

Parameters:
  • arrays – tuple of 3D real arrays of the same shape.

  • axis_to_align – the axis to be aligned (e.g. vector q), expressed in an orthonormal frame x y z

  • reference_axis – will align axis_to_align onto this vector, expressed in an orthonormal frame x y z

  • voxel_size – tuple, voxel size of the 3D array in z, y, and x (CXI convention)

  • fill_value – tuple of numeric values used in the RegularGridInterpolator for points outside of the interpolation domain. The length of the tuple should be equal to the number of input arrays.

  • rotation_matrix – optional numpy ndarray of shape (3, 3), rotation matrix to apply to arrays. If it is provided, the parameters axis_to_align and reference_axis will be discarded.

  • is_orthogonal – set to True is the frame is orthogonal, False otherwise (detector frame) Used for plot labels.

  • reciprocal_space – True if the data is in reciprocal space, False otherwise. Used for plot labels.

  • debugging – tuple of booleans of the same length as the number of input arrays, True to see plots before and after rotation

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’title’: tuple of strings, titles for the debugging plots, same length as the number of arrays

    • ’scale’: tuple of strings (either ‘linear’ or ‘log’), scale for the debugging plots, same length as the number of arrays

    • width_z: size of the area to plot in z (axis 0), centered on the middle of the initial array

    • width_y: size of the area to plot in y (axis 1), centered on the middle of the initial array

    • width_x: size of the area to plot in x (axis 2), centered on the middle of the initial array

Returns:

a rotated array (if a single array was provided) or a tuple of rotated arrays (same length as the number of input arrays)

bcdi.utils.utilities.rotate_vector(vectors: numpy.ndarray | Tuple[numpy.ndarray, ...], axis_to_align: numpy.ndarray | None = None, reference_axis: numpy.ndarray | None = None, rotation_matrix: numpy.ndarray | None = None) numpy.ndarray | Tuple[numpy.ndarray, ...]

Rotate vectors.

It calculates the vector components (3D) in the basis where axis_to_align and reference_axis are aligned. axis_to_align and reference_axis should be in the order X Y Z, where Z is downstream, Y vertical and X outboard (CXI convention).

Parameters:
  • vectors – the vectors to be rotated, tuple of three components (values or 1D-arrays) expressed in an orthonormal frame x y z

  • axis_to_align – the axis of myobj (vector q), expressed in an orthonormal frame x y z

  • reference_axis – will align axis_to_align onto this vector, expressed in an orthonormal frame x y z

  • rotation_matrix – optional numpy ndarray of shape (3, 3), rotation matrix to apply to vectors. If it is provided, the parameters axis_to_align and reference_axis will be discarded.

Returns:

tuple of three ndarrays in CXI convention z y x, each of shape (vectors[0].size, vectors[1].size, vectors[2].size). If a single vector is provided, returns a 1D array of size 3.

bcdi.utils.utilities.rotation_matrix_3d(axis_to_align, reference_axis)

Calculate the rotation matrix which aligns axis_to_align onto reference_axis in 3D.

Parameters:
  • axis_to_align – the 3D vector to be aligned (e.g. vector q), expressed in an orthonormal frame x y z

  • reference_axis – will align axis_to_align onto this 3D vector, expressed in an orthonormal frame x y z

Returns:

the rotation matrix as a np.array of shape (3, 3)

bcdi.utils.utilities.skewed_gaussian(x_axis, amp, loc, sig, alpha)

Skewed Gaussian line shape.

Parameters:
  • x_axis – where to calculate the function

  • amp – the amplitude of the Gaussian

  • loc – the location parameter

  • sig – HWHM of the Gaussian

  • alpha – the shape parameter

Returns:

the skewed Gaussian line shape at x_axis

bcdi.utils.utilities.smaller_primes(number, maxprime=13, required_dividers=(4,))

Find the closest smaller number that meets some condition.

Find the closest integer <=n (or list/array of integers), for which the largest prime divider is <=maxprime, and has to include some dividers. The default values for maxprime is the largest integer accepted by the clFFT library for OpenCL GPU FFT. Adapted from PyNX.

Parameters:
  • number – the integer number

  • maxprime – the largest prime factor acceptable

  • required_dividers – a list of required dividers for the returned integer.

Returns:

the integer (or list/array of integers) fulfilling the requirements

bcdi.utils.utilities.sum_roi(array, roi, debugging=False)

Sum the array intensities in the defined region of interest.

Parameters:
  • array – 2D or 3D array. If ndim=3, the region of interest is applied sequentially to each 2D frame, the iteration being peformed over the first axis.

  • roi – [Vstart, Vstop, Hstart, Hstop] region of interest for the sum

  • debugging – True to see plots

Returns:

a number (if array.ndim=2) or a 1D array of length array.shape[0] (if array.ndim=3) of summed intensities

bcdi.utils.utilities.try_smaller_primes(number, maxprime=13, required_dividers=(4,))

Check if a number meets some condition.

Check if the largest prime divider is <=maxprime, and optionally includes some dividers. Adapted from PyNX.

Parameters:
  • number – the integer number for which the prime decomposition will be checked

  • maxprime – the maximum acceptable prime number. This defaults to the largest integer accepted by the clFFT library for OpenCL GPU FFT.

  • required_dividers – list of required dividers in the prime decomposition. If None, this check is skipped.

Returns:

True if the conditions are met.

bcdi.utils.utilities.unpack_array(array: float | List[float] | numpy.ndarray) float | numpy.ndarray

Unpack an array or Sequence of length 1 into a single element.

bcdi.utils.utilities.update_frames_logical(frames_logical: numpy.ndarray, logical_subset: numpy.ndarray) numpy.ndarray

Update frames_logical with a logical array of smaller length.

The number of non-zero elements of frames_logical should be equal to the length of the logical subset.

bcdi.utils.utilities.upsample(array: numpy.ndarray | List, factor: int = 2) numpy.ndarray

Upsample an array.

Parameters:
  • array – the numpy array to be upsampled

  • factor – int, factor for the upsampling

Returns:

the upsampled numpy array

bcdi.utils.utilities.wrap(obj, start_angle, range_angle)

Wrap obj between start_angle and (start_angle + range_angle).

Parameters:
  • obj – number or array to be wrapped

  • start_angle – start angle of the range

  • range_angle – range

Returns:

wrapped angle in [start_angle, start_angle+range[

format

Functions and class related to string formatting for JSON dump and repr.

API Reference

Functions related to formatting for string representation.

class bcdi.utils.format.CustomEncoder(*, skipkeys=False, ensure_ascii=True, check_circular=True, allow_nan=True, sort_keys=False, indent=None, separators=None, default=None)

Class to handle the serialization of np.ndarrays, sets.

default(obj)

Override the JSONEncoder.default method to support more types.

bcdi.utils.format.create_repr(obj: Any, cls: type) str

Generate the string representation of the object.

It uses the parameters given to __init__, except self, args and kwargs.

Parameters:
  • obj – the object for which the string representation should be generated

  • cls – the cls from which __init__ parameters should be extracted (e.g., base class in case of inheritance)

Returns:

the string representation

bcdi.utils.format.format_repr(value: Any | None, quote_mark: bool = True) str

Format strings for the __repr__ method.

Parameters:
  • value – string or None

  • quote_mark – True to put quote marks around strings

Returns:

a string

bcdi.utils.format.format_value(value: Any) Tuple[Any, bool]

Format the value for a proper representation.

bcdi.utils.format.ndarray_to_list(array: numpy.ndarray) List

Convert a numpy ndarray of any dimension to a nested list.

Parameters:

array – the array to be converted

image_registration

Functions related to image registration and the alignment of datasets in direct or reciprocal space.

API Reference

Functions related to the registration and alignement of two arrays.

bcdi.utils.image_registration.align_arrays(reference_array, shifted_array, shift_method='modulus', interpolation_method='subpixel', support_threshold=None, precision=1000, verbose=True, debugging=False, **kwargs)

Align two arrays using dft registration and subpixel shift.

The shift between arrays can be determined either using the modulus of the arrays or a support created from it using a threshold.

Parameters:
  • reference_array – 3D array, reference object

  • shifted_array – 3D array to be aligned, same shape as reference_obj

  • shift_method – ‘raw’, ‘modulus’, ‘support’ or ‘skip’. Object to use for the determination of the shift. If ‘raw’, it uses the raw, eventually complex array. if ‘modulus’, it uses the modulus of the array. If ‘support’, it uses a support created by threshold the modulus of the array.

  • interpolation_method – ‘subpixel’ for interpolating using subpixel shift, ‘rgi’ for interpolating using a RegularGridInterpolator, ‘roll’ to shift voxels by an integral number (the shifts are rounded to the nearest integer)

  • support_threshold – all points where the normalized modulus is larger than this value will be set to 1 in the support.

  • precision – precision for the DFT registration in 1/pixel

  • verbose – boolean, True to print comments

  • debugging – boolean, set to True to see plots

Returns:

  • the aligned array

  • the shift: tuple of floats

bcdi.utils.image_registration.align_diffpattern(reference_data, data, mask=None, shift_method='raw', interpolation_method='roll', verbose=True, debugging=False, **kwargs)

Align two diffraction patterns.

The alignement can be based either on the shift of the center of mass or on dft registration.

Parameters:
  • reference_data – the first 3D or 2D diffraction intensity array which will serve as a reference.

  • data – the 3D or 2D diffraction intensity array to align.

  • mask – the 3D or 2D mask corresponding to data

  • shift_method – ‘raw’, ‘modulus’, ‘support’ or ‘skip’. Object to use for the determination of the shift. If ‘raw’, it uses the raw, eventually complex array. if ‘modulus’, it uses the modulus of the array. If ‘support’, it uses a support created by threshold the modulus of the array.

  • interpolation_method – ‘rgi’ for RegularGridInterpolator or ‘subpixel’ for subpixel shift

  • verbose – boolean, True to print comments

  • debugging – boolean, set to True to see plots

Returns:

  • the shifted data

  • the shifted mask

  • if return_shift, returns a tuple containing the shifts

bcdi.utils.image_registration.average_arrays(avg_obj, ref_obj, obj, support_threshold=0.25, correlation_threshold=0.9, aligning_option='dft', space='reciprocal_space', debugging=False, **kwargs)

Average two reconstructions after aligning it.

This function can be used to average a series of arrays within a loop. Alignment is performed using either DFT registration or the shift of the center of mass of the array. Averaging is processed only if their Pearson cross-correlation after alignment is larger than the correlation threshold.

Parameters:
  • avg_obj – 3D array of complex numbers, current average

  • ref_obj – 3D array of complex numbers, used as a reference for the alignment

  • obj – 3D array of complex numbers, array to be aligned with the reference and to be added to avg_obj

  • support_threshold – normalized threshold for the definition of the support. It is applied on the modulus of the array

  • correlation_threshold – float in [0, 1], minimum correlation between two dataset to average them

  • aligning_option – ‘com’ for center of mass, ‘dft’ for dft registration and subpixel shift

  • space – ‘direct_space’ or ‘reciprocal_space’, in which space the average will be performed

  • debugging – boolean, set to True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’width_z’: size of the area to plot in z (axis 0), centered on the middle of the initial array

    • ’width_y’: size of the area to plot in y (axis 1), centered on the middle of the initial array

    • ’width_x’: size of the area to plot in x (axis 2), centered on the middle of the initial array

    • ’reciprocal_space’: True if the object is in reciprocal space, it is used only for defining labels in plots

    • ’is_orthogonal’: True if the data is in an orthonormal frame. Used for defining default plot labels.

Returns:

the average complex density

bcdi.utils.image_registration.calc_new_positions(old_positions: list, shift: Sequence[float]) numpy.ndarray

Transform old_positions depending on the shift.

Parameters:
  • old_positions – list of 1D arrays corresponding to the position of the voxels in a regular grid before transformation, in each dimension. For example, if the array to be interpolated is 3D, old_positions will be a list of three 1D arrays [array0, array1, array2], where array0 describes the voxel positions along axis 0, array1 along axis 1 and array2 along axis 2. It does not work if the grid is not regular (each coordinate would need to be described by a 3D array instead).

  • shift – a tuple of floats, corresponding to the shift in each dimension that need to be applied to array

Returns:

the shifted positions where to interpolate, for the RegularGridInterpolator

bcdi.utils.image_registration.dft_registration(buf1ft, buf2ft, ups_factor=100)

Efficient subpixel image registration by cross-correlation.

This code gives the same precision as the FFT upsampled cross correlation in a small fraction of the computation time and with reduced memory requirements. It obtains an initial estimate of the cross-correlation peak by an FFT and then refines the shift estimation by upsampling the DFT only in a small neighborhood of that estimate by means of a matrix-multiply DFT. With this procedure all the image points are used to compute the upsampled cross-correlation. Manuel Guizar - Dec 13, 2007

Portions of this code were taken from code written by Ann M. Kowalczyk and James R. Fienup. J.R. Fienup and A.M. Kowalczyk, “Phase retrieval for a complex-valued object by using a low-resolution image,” J. Opt. Soc. Am. A 7, 450-458 (1990).

Citation for this algorithm: Manuel Guizar-Sicairos, Samuel T. Thurman, and James R. Fienup, “Efficient subpixel image registration algorithms,” Opt. Lett. 33, 156-158 (2008).

Parameters:
  • buf1ft – Fourier transform of reference image, DC in (1,1) [DO NOT FFTSHIFT]

  • buf2ft – Fourier transform of image to register, DC in (1,1) [DO NOT FFTSHIFT]

  • ups_factor – upsampling factor (integer). Images will be registered to within 1/ups_factor of a pixel. For example ups_factor = 20 means the images will be registered within 1/20 of a pixel. (default = 1)

Returns:

  • output: [error,diff_phase,net_row_shift,net_col_shift]

  • error: translation invariant normalized RMS error between f and g

  • diff_phase: global phase difference between the two images (should be zero if images are non-negative).

  • row_shift, col_shift: pixel shifts between images

bcdi.utils.image_registration.dftups(array, output_row_nb, output_column_nb, ups_factor=1, row_offset=0, column_offset=0)

Upsampled DFT by matrix multiplies.

It can compute an upsampled DFT in just a small region. Receives DC in upper left corner, image center must be in (1,1). Manuel Guizar - Dec 13, 2007 Modified from dftups, by J.R. Fienup 7/31/06

This code is intended to provide the same result as if the following operations were performed:

  • Embed the array “in” in an array that is usfac times larger in each dimension. ifftshift to bring the center of the image to (1,1).

  • Take the FFT of the larger array

  • Extract an [output_row_nb, output_column_nb] region of the result. Starting with the [row_offset+1 column_offset+1] element.

It achieves this result by computing the DFT in the output array without the need to zero pad. Much faster and memory efficient than the zero-padded FFT approach if [nor noc] are much smaller than [nr*usfac nc*usfac]

Parameters:
  • array – 2D input array

  • output_row_nb – number of pixels in the output upsampled DFT, in units of upsampled pixels (default = size(in))

  • output_column_nb – number of pixels in the output upsampled DFT, in units of upsampled pixels (default = size(in))

  • ups_factor – upsampling factor (default ups_factor = 1)

  • row_offset – row offset, allow to shift the output array to a region of interest on the DFT (default = 0)

  • column_offset – column offset, allow to shift the output array to a region of interest on the DFT (default = 0)

bcdi.utils.image_registration.get_shift(reference_array: numpy.ndarray, shifted_array: numpy.ndarray, shift_method: str = 'modulus', precision: int = 1000, support_threshold: None | float = None, verbose: bool = True, **kwargs) Sequence[float]

Calculate the shift between two arrays.

The shift is calculated using dft registration. If a threshold for creating a support is not provided, DFT registration is performed on the arrays themselves. If a threshold is provided, shifts are calculated using the support created by thresholding the modulus of the arrays.

Parameters:
  • reference_array – numpy ndarray

  • shifted_array – numpy ndarray of the same shape as reference_array

  • shift_method – ‘raw’, ‘modulus’ or ‘support’. Object to use for the determination of the shift. If ‘raw’, it uses the raw, eventually complex array. if ‘modulus’, it uses the modulus of the array. If ‘support’, it uses a support created by threshold the modulus of the array.

  • precision – precision for the DFT registration in 1/pixel

  • support_threshold – optional normalized threshold in [0, 1]. If not None, it will be used to define a support. The center of mass will be calculated for that support instead of the modulus.

  • verbose – True to print comment

  • kwargs

    • ‘logger’: an optional logger

Returns:

list of shifts, of length equal to the number of dimensions of the arrays. These are shifts that need to be applied to shifted_array in order to align it with reference_array (no need to flip signs)

bcdi.utils.image_registration.getimageregistration(array1, array2, precision=10)

Calculate the registration (shift) between two arrays.

Parameters:
  • array1 – the reference array

  • array2 – the array to register

  • precision – subpixel precision of the registration. Images will be registered to within 1/precision of a pixel

Returns:

the list of shifts that need to be applied to array2 in order to align it with array1 (no need to flip signs)

bcdi.utils.image_registration.index_max(mydata)

Look for the data max and location.

bcdi.utils.image_registration.index_max1(mydata)

Look for the data maximum locations.

bcdi.utils.image_registration.interp_rgi_translation(array: numpy.ndarray, shift: Sequence[float]) numpy.ndarray

Interpolate the shifted array on new positions using a RegularGridInterpolator.

Parameters:
  • array – a numpy array expressed in a regular grid

  • shift – a tuple of floats, corresponding to the shift in each dimension that need to be applied to array

Returns:

the shifted array

bcdi.utils.image_registration.shift_array(array: numpy.ndarray, shift: Sequence[float], interpolation_method: str = 'subpixel') numpy.ndarray

Shift array using the defined method given the offsets.

Parameters:
  • array – a numpy ndarray

  • shift – tuple of floats, shifts of the array in each dimension

  • interpolation_method – ‘subpixel’ for interpolating using subpixel shift, ‘rgi’ for interpolating using a RegularGridInterpolator, ‘roll’ to shift voxels by an integral number (the shifts are rounded to the nearest integer) created by thresholding the modulus of the array.

Returns:

the shifted array

bcdi.utils.image_registration.subpixel_shift(array, z_shift, y_shift, x_shift=None)

Shift array by the shift values.

Adapted from the Matlab code of Jesse Clark.

Parameters:
  • array – array to be shifted

  • z_shift – shift in the first dimension

  • y_shift – shift in the second dimension

  • x_shift – shift in the third dimension

Returns:

the shifted array

io_helper

Class and decorators related to the safe opening of files.

API Reference

Module containing decorators and context manager classes for input-output.

class bcdi.utils.io_helper.ContextFile(filename: str, open_func: type | Callable, scan_number: int | None = None, mode: str = 'r', encoding: str = 'utf-8', longname: str | None = None, shortname: str | None = None, directory: str | None = None, **kwargs)

Convenience context manager to open files.

The supported opening callables are silx.io.specfile.Specfile, io.open, h5py.File and SIXS (nxsReady.Dataset, ReadNxs3.Dataset).

Parameters:

kwargs

  • ‘logger’: an optional logger

bcdi.utils.io_helper.safeload(func: Callable) Callable

Decorate a class method to safely opening files.

Parameters:

func – a class method accessing the file

bcdi.utils.io_helper.safeload_static(func: Callable) Callable

Decorate a class static method or a function to safely opening files.

Parameters:

func – a class static method accessing the file

parameters

Classes and functions related to the validation of configuration parameters.

API Reference

Validation of configuration parameters.

The validation is performed only on the expected parameters. Other parameters are simply discarded.

class bcdi.utils.parameters.CDIPreprocessingChecker(initial_params: Dict[str, Any], default_values: Dict[str, Any] | None = None, match_length_params: Tuple = (), required_params: Tuple = ())

Configure preprocessing-dependent parameters for the ‘CDI’ case.

class bcdi.utils.parameters.ConfigChecker(initial_params: Dict[str, Any], default_values: Dict[str, Any] | None = None, match_length_params: Tuple = (), required_params: Tuple = ())

Validate and configure parameters.

Parameters:
  • initial_params – the dictionary of parameters to validate and configure

  • default_values – an optional dictionary of default values for keys in initial_params

  • match_length_params – a tuple of keys from initial_params which should match a certain length (e.g. the number of scans)

  • required_params – a tuple of keys that have to be present in initial_params

check_config() Dict[str, Any]

Check if the provided config is consistent.

exception bcdi.utils.parameters.MissingKeyError

Custom Exception for a missing required key in the config dictionary.

exception bcdi.utils.parameters.ParameterError(key, value, allowed)

Custom Exception for a dictionary of parsed parameters.

Parameters:
  • key – a key of the dictionary

  • value – the corresponding value

  • allowed – allowed values for value

class bcdi.utils.parameters.PostprocessingChecker(initial_params: Dict[str, Any], default_values: Dict[str, Any] | None = None, match_length_params: Tuple = (), required_params: Tuple = ())

Configure postprocessing-dependent parameters.

check_config() Dict[str, Any]

Check if the provided config is consistent.

class bcdi.utils.parameters.PreprocessingChecker(initial_params: Dict[str, Any], default_values: Dict[str, Any] | None = None, match_length_params: Tuple = (), required_params: Tuple = ())

Configure preprocessing-dependent parameters.

bcdi.utils.parameters.valid_param(key: str, value: Any) Tuple[Any, bool]

Validate a key value pair corresponding to an input parameter.

It will raise an exception if the check fails.

Parameters:
  • key – name of the parameter

  • value – the value of the parameter

Returns:

a tuple (formatted_value, is_valid). is_valid is True if the key is valid, False otherwise.

parser

Class and functions related to loading of the config file and command line parameters.

API Reference

Parsing of command-line arguments and config files.

class bcdi.utils.parser.ConfigParser(file_path: str, command_line_args: Dict[str, Any] | None = None)

Base class for parsing arguments.

If provided, command line arguments will override the parameters from the config file. Some validation is also realized in the class.

Parameters:
  • file_path – path of the configuration file that contains the arguments, str.

  • command_line_args – an optional dictionary of command-line arguments

property command_line_args

Return an optional dictionary of command line parameters.

property file_path

Path of the configuration file.

static filter_dict(dic: Dict | None, filter_value: Any = None) Dict

Filter out key where the value is None.

Parameters:
  • dic – a dictionary

  • filter_value – value to be filtered out

Returns:

a dictionary with only keys where the value is not None

load_arguments() Dict

Parse the byte string, eventually override defaults and check parameters.

bcdi.utils.parser.add_cli_parameters(argument_parser: ArgumentParser) ArgumentParser

Add generic parameters to the argument parser.

Parameters:

argument_parser – an instance of argparse.ArgumentParser

Returns:

the updated instance

bcdi.utils.parser.str_to_list(string: str, item_type: Type) List

Convert a comma-separated string to a list.

Parameters:
  • string – the string to convert, e.g. “11,12,13” or “file1.spec,file2.spec”

  • item_type – type to which the elements should be cast

Returns:

a list of converted values

snippets_logging

API Reference

Classes and functions for the configuration of logging.

class bcdi.utils.snippets_logging.ColorLogFormatter(fmt=None, datefmt=None, style='%', validate=True, *, defaults=None)

Format colored logs.

format(record)

Format log records.

It uses a default prefix and suffix to terminal color codes that corresponds to the log level name.

class bcdi.utils.snippets_logging.LoggingColor

Define terminal color codes.

bcdi.utils.snippets_logging.configure_logging(path: str | None = None, console_level: int = 20, file_level: int = 10) None

Create handlers and configure logging.

Parameters:
  • path – path of the file where to write logs

  • console_level – valid level for logging to the console

  • file_level – valid level for logging to a file

text

This module provides the class used to define filename when saving data and figures during processing.

API Reference

Implementation of the class used for comments in plots and file names.

class bcdi.utils.text.Comment(text: str)

Comment used in logging message during the analysis.

validation

This module contains all validation functions that are used throughout the code, in order to test the input parameters of functions.

API Reference

Functions related to the validation of input parameters.

bcdi.utils.validation.is_float(string)

Return True is the string represents a number.

Parameters:

string – the string to be checked

Returns:

True of False

bcdi.utils.validation.valid_1d_array(array, length=None, min_length=None, allow_none=True, allowed_types=None, allowed_values=None, name=None)

Check if the array is 1D and satisfies the requirements.

Parameters:
  • array – the numpy array to be checked

  • length – int, required length

  • min_length – int, minimum length of the array

  • allow_none – bool, True if the array can be None

  • allowed_types – list or tuple of valid types

  • allowed_values – Sequence of allowed values for the array

  • name – name of the calling object appearing in exception messages

Returns:

bool, True if all checks passed

bcdi.utils.validation.valid_container(obj, container_types, length=None, min_length=None, max_length=None, item_types=None, name=None, min_included=None, min_excluded=None, max_included=None, max_excluded=None, allow_none=False)

Check that the input object as three elements fulfilling the defined requirements.

Parameters:
  • obj – the object to be tested

  • container_types – list of the allowed types for obj

  • length – int, required length

  • min_length – mininum length (inclusive)

  • max_length – maximum length (inclusive)

  • item_types – list of the allowed types for the object items

  • min_included – minimum allowed value (inclusive)

  • min_excluded – minimum allowed value (exclusive)

  • max_included – maximum allowed value (inclusive)

  • max_excluded – maximum allowed value (exclusive)

  • allow_none – True if the container items are allowed to be None

  • name – name of the calling object appearing in exception messages

Returns:

True if checks pass, raise some error otherwise

bcdi.utils.validation.valid_item(value, allowed_types, min_included=None, min_excluded=None, max_included=None, max_excluded=None, allow_none=False, name=None)

Check that the input object as three elements fulfilling the defined requirements.

Parameters:
  • value – the value to be tested

  • allowed_types – allowed types of the object values

  • min_included – minimum allowed value (inclusive)

  • min_excluded – minimum allowed value (exclusive)

  • max_included – maximum allowed value (inclusive)

  • max_excluded – maximum allowed value (exclusive)

  • allow_none – True if the container items are allowed to be None

  • name – name of the calling object appearing in exception messages

Returns:

True if checks pass, raise some error otherwise

bcdi.utils.validation.valid_kwargs(kwargs, allowed_kwargs, name='kwargs')

Check if the provided parameters belong to the set of allowed kwargs.

Parameters:
  • kwargs – dictionnary of kwargs to check

  • allowed_kwargs – set of allowed keys

  • name – name of the calling object appearing in exception messages

Returns:

True if checks pass, raise some error otherwise

bcdi.utils.validation.valid_ndarray(arrays: numpy.ndarray | Tuple[numpy.ndarray, ...], ndim: int | Tuple[int, ...] | None = None, shape: Tuple[int, ...] | None = None, fix_ndim: bool = True, fix_shape: bool = True, name: str | None = None) bool

Check that arrays have the same shape and the correct number of dimensions.

Parameters:
  • arrays – a sequence of numpy ndarrays

  • ndim – int, the number of dimensions to be compared with

  • shape – sequence of int, shape to be compared with

  • fix_ndim – bool, if True the shape of all arrays should be equal

  • fix_shape – bool, if True the shape of all arrays should be equal

  • name – name of the calling object appearing in exception messages

Returns:

True if checks pass, raise some error otherwise