bcdi.postprocessing: Post-processing of data after phase retrieval

facet_recognition

This module provides tools for facet segmentation on nanocrystals, and for plotting the stereographic projection of a diffraction peak or an object. For more details, see: Carnis et al. Small 17, 2007702 (2021) https://doi.org/10.1002/smll.202007702

API Reference

Functions related to facet recognition of nanocrystals.

bcdi.postprocessing.facet_recognition.calc_stereoproj_facet(projection_axis, vectors, radius_mean, stereo_center)

Calculate the coordinates of normals in the stereographic projection.

The calculation depends on the reference axis. See: Nanoscale 10, 4833 (2018).

Parameters:
  • projection_axis – the projection is performed on q plane perpendicular to that axis (0, 1 or 2)

  • vectors – array of vectors to be projected (nb_vectors rows x 3 columns)

  • radius_mean – q radius from which the projection will be done

  • stereo_center – offset of the projection plane along the reflection axis, in the same unit as radius_mean. If stereo_center = 0, the projection plane will be the equator.

Returns:

the coordinates of the stereographic projection for the projection from the South pole(1st and 2nd columns) and from the North pole (3rd and 4th columns) projection, rescaled from radius_mean to 90 degrees

bcdi.postprocessing.facet_recognition.detect_edges(faces)

Find indices of vertices defining non-shared edges.

Parameters:

faces – ndarray of m*3 faces

Returns:

1D list of indices of vertices defining non-shared edges (near hole…)

bcdi.postprocessing.facet_recognition.distance_threshold(fit, indices, plane_shape, max_distance=0.9)

Filter out pixels depending on their distance to a fit plane.

Parameters:
  • fit – coefficients of the plane (a, b, c, d) such that a*x + b*y + c*z + d = 0

  • indices – tuple or array of plane indices, x being the 1st tuple element or array row, y the 2nd tuple element or array row and z the third tuple element or array row

  • plane_shape – shape of the initial plane array

  • max_distance – max distance allowed from the fit plane in pixels

Returns:

the updated plane, a stop flag

bcdi.postprocessing.facet_recognition.equirectangular_proj(normals, intensity, cmap=matplotlib.colors.ListedColormap, bw_method=0.03, min_distance=10, background_threshold=-0.35, debugging=False)

Detect facets in an object.

It uses an equirectangular projection of normals to mesh triangles and watershed segmentation.

Parameters:
  • normals – normals array

  • intensity – intensity array

  • cmap – colormap used for plotting

  • bw_method – bw_method of gaussian_kde

  • min_distance – min_distance of corner_peaks()

  • background_threshold – threshold for background determination (depth of the KDE)

  • debugging – if True, show plots for debugging

Returns:

ndarray of labelled regions

bcdi.postprocessing.facet_recognition.find_facet(refplane_indices, surf_indices, original_shape, step_shift, plane_label, plane_coeffs, min_points, debugging=False)

Shift a fit plane along its normal until it reaches the surface of a faceted object.

Parameters:
  • refplane_indices – a tuple of 3 arrays (1D, length N) describing the coordinates of the plane voxels, x values being the 1st tuple element, y values the 2nd tuple element and z values the 3rd tuple element (output of np.nonzero)

  • surf_indices – a tuple of 3 arrays (1D, length N) describing the coordinates of the surface voxels, x values being the 1st tuple element, y values the 2nd tuple element and z values the 3rd tuple element (output of np.nonzero)

  • original_shape – the shape of the full dataset (amplitude object, eventually upsampled)

  • step_shift – the amplitude of the shift to be applied to the plane along its normal

  • plane_label – the label of the plane, used in comments

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

  • min_points – threshold, minimum number of points that should coincide between the fit plane and the object surface

  • debugging – True to see debugging plots

Returns:

the shift that needs to be applied to the fit plane in order to best match with the object surface

bcdi.postprocessing.facet_recognition.find_neighbours(vertices, faces)

Get the list of neighbouring vertices for each vertex.

Parameters:
  • vertices – ndarray of n*3 vertices

  • faces – ndarray of m*3 faces

Returns:

list of lists of indices

bcdi.postprocessing.facet_recognition.fit_plane(plane, label, debugging=False)

Fit a plane to labelled indices using the equation a*x+ b*y + c*z + d = 0.

Parameters:
  • plane – 3D binary array, where the voxels belonging to the plane are set to 1 and others are set to 0.

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

  • debugging – show plots for debugging

Returns:

fit parameters (a, b, c, d), plane indices after filtering, errors associated, a stop flag

bcdi.postprocessing.facet_recognition.grow_facet(fit, plane, label, support, max_distance=0.9, debugging=True)

Find voxels of the object which belong to a facet.

It uses the facet plane equation and the distance to the plane to find such voxels.

Parameters:
  • fit – coefficients of the plane (a, b, c, d) such that a*x + b*y + c*z + d = 0

  • plane – 3D binary support of the plane, with shape of the full dataset

  • label – the label of the plane processed

  • support – 3D binary support of the reconstructed object, with shape of the full dataset

  • max_distance – in pixels, maximum allowed distance to the facet plane of a voxel

  • debugging – set to True to see plots

Returns:

the updated plane, a stop flag

bcdi.postprocessing.facet_recognition.offset_plane(indices, offset, plane_normal)

Shift plane indices by the offset value in order to scan perpendicular to the plane.

Parameters:
  • indices – tuple of 3 1D ndarrays (array shape = nb_points)

  • offset – offset to be applied to the indices (offset of the plane)

  • plane_normal – ndarray of 3 elements, normal to the plane

Returns:

offseted indices

bcdi.postprocessing.facet_recognition.remove_duplicates(vertices, faces, debugging=False)

Remove duplicates in a list of vertices and faces.

A face is a triangle made of three vertices.

Parameters:
  • vertices – a ndarray of vertices, shape (N, 3)

  • faces – a ndarray of vertex indices, shape (M, 3)

  • debugging – True to see which vertices are duplicated and how lists are modified

Returns:

the updated vertices and faces with duplicates removed in place

bcdi.postprocessing.facet_recognition.stereographic_proj(normals, intensity, max_angle, savedir, voxel_size, projection_axis, min_distance=10, background_south=-1000, background_north=-1000, save_txt=False, cmap=matplotlib.colors.ListedColormap, planes_south=None, planes_north=None, plot_planes=True, scale='linear', comment_fig='', debugging=False)

Detect facets in an object.

It uses a stereographic projection of normals to mesh triangles and watershed segmentation.

Parameters:
  • normals – array of normals to mesh triangles (nb_normals rows x 3 columns)

  • intensity – array of intensities (nb_normals rows x 1 column)

  • max_angle – maximum angle in degree of the stereographic projection (should be larger than 90)

  • savedir – directory for saving figures

  • voxel_size – tuple of three numbers corresponding to the real-space voxel size in each dimension

  • projection_axis – the projection is performed on a plane perpendicular to that axis (0, 1 or 2)

  • min_distance – min_distance of corner_peaks()

  • background_south – threshold for background determination in the projection from South

  • background_north – threshold for background determination in the projection from North

  • save_txt – if True, will save coordinates in a .txt file

  • cmap – colormap used for plotting pole figures

  • planes_south – dictionnary of crystallographic planes, e.g. {‘111’:angle_with_reflection}

  • planes_north – dictionnary of crystallographic planes, e.g. {‘111’:angle_with_reflection}

  • plot_planes – if True, will draw circles corresponding to crystallographic planes in the pole figure

  • scale – ‘linear’ or ‘log’, scale for the colorbar of the plot

  • comment_fig – string, comment for the filename when saving figures

  • debugging – show plots for debugging

Returns:

  • labels_south and labels_north as 2D arrays for each projection from South and North

  • a (Nx4) array: projected coordinates of normals from South (u column 0, v column 1) and North (u column2 , v column 3). The coordinates are in degrees, not indices.

  • the list of rows to remove

bcdi.postprocessing.facet_recognition.surface_gradient(points, support, width=2)

Calculate the support gradient at point.

Parameters:
  • points – tuple or list of tuples of 3 integers (z, y, x), position where to calculate the gradient vector

  • support – 3D numpy binary array, being 1 in the crystal and 0 outside

  • width – half-width of the window where the gradient will be calculated (the support gradient is nonzero on a single layer, it avoids missing it)

Returns:

a list of normalized vector(s) (array(s) of 3 numbers) oriented towards the exterior of the cristal

bcdi.postprocessing.facet_recognition.surface_indices(surface, plane_indices, margin=3)

Find surface indices potentially belonging to a plane.

It crops the surface around the plane with a certain margin, and find corresponding surface indices.

Parameters:
  • surface – the 3D surface binary array

  • plane_indices – tuple of 3 1D-arrays of plane indices

  • margin – margin to include aroung plane indices, in pixels

Returns:

3*1D arrays of surface indices

bcdi.postprocessing.facet_recognition.taubin_smooth(faces, vertices, cmap=matplotlib.colors.ListedColormap, iterations=10, lamda=0.33, mu=0.34, radius=0.1, debugging=False)

Perform Taubin’s smoothing of a mesh.

It performs a back and forward Laplacian smoothing “without shrinking” of a triangulated mesh, as described by Gabriel Taubin (ICCV ‘95)

Parameters:
  • faces – m*3 ndarray of m faces defined by 3 indices of vertices

  • vertices – n*3 ndarray of n vertices defined by 3 positions

  • cmap – colormap used for plotting

  • iterations – number of iterations for smoothing

  • lamda – smoothing variable 0 < lambda < mu < 1

  • mu – smoothing variable 0 < lambda < mu < 1

  • radius – radius around which the normals are integrated in the calculation of the density of normals

  • debugging – show plots for debugging

Returns:

smoothened vertices (ndarray n*3), normals to triangle (ndarray m*3), weighted density of normals, updated faces, errors

bcdi.postprocessing.facet_recognition.update_logfile(support, strain_array, summary_file, allpoints_file, label=0, angle_plane=numpy.nan, plane_coeffs=(0, 0, 0, 0), plane_normal=(0, 0, 0))

Update log files use in the facet_strain.py script.

Parameters:
  • support – the 3D binary support defining voxels to be saved in the logfile

  • strain_array – the 3D strain array

  • summary_file – the handle for the file summarizing strain statistics per facet

  • allpoints_file – the handle for the file giving the strain and the label for each voxel

  • label – the label of the plane

  • angle_plane – the angle of the plane with the measurement direction

  • plane_coeffs – the fit coefficients (a,b,c,d) of the plane such that ax+by+cz+d=0

  • plane_normal – the normal to the plane

Returns:

nothing

bcdi.postprocessing.facet_recognition.upsample(array, upsampling_factor, voxelsizes=None, title='', debugging=False)

Upsample array using a factor of upsampling.

Parameters:
  • array – the real array to be upsampled

  • upsampling_factor – int, the upsampling factor

  • voxelsizes – list, the voxel sizes of array

  • title – title for the debugging plot

  • debugging – True to see plots

Returns:

the upsampled array

analysis

This module provides the classes used in postprocessing strain analysis. It is mainly an abstration layer making analysis steps more understandable.

API Reference

Implementation of postprocessing analysis classes.

class bcdi.postprocessing.analysis.Analysis(scan_index: int, parameters: Dict[str, Any], setup: Setup, **kwargs)

Base class for the post-processing analysis workflow.

average_reconstructions() None

Create an average reconstruction out of many, after proper alignment.

center_object_based_on_modulus(**kwargs)

Center the object based on one of the methods ‘max’, ‘com’ and ‘max_com’.

Parameters:

kwargs

  • ‘centering_method’: centering method name among ‘max’, ‘com’ and ‘max_com’

property get_interplanar_distance: float | None

Calculate the interplanar distance in nm.

get_optical_path() numpy.ndarray

Calculate the optical path through the crystal.

abstract get_q_bragg_crystal_frame() numpy.ndarray | None

Return the q vector of the Bragg peak in the crystal frame.

abstract interpolate_into_crystal_frame()

Interpolate the direct space object into the crystal frame.

The exact steps depend on which frame the data lies in.

property voxel_sizes

List of three positive numbers, voxel size in nm in each dimension.

class bcdi.postprocessing.analysis.DetectorFrameLinearization(scan_index: int, parameters: Dict[str, Any], setup: Setup, **kwargs)

Analysis workflow for interpolation based on the linearization matrix.

The data before interpolation is in the detector frame (axis 0 = rocking angle, axis 1 detector Y, axis 2 detector X). The data after interpolation is in the pseudo-crystal frame (only ‘ref_axis_q’is an axis of the crystal frame, there is an unknown inplane rotation around it).

get_q_bragg_crystal_frame() numpy.ndarray | None

Return the q vector of the Bragg peak in the crystal frame.

interpolate_into_crystal_frame() None

Interpolate the data in the pseudo crystal frame.

class bcdi.postprocessing.analysis.InterpolatedCrystal(modulus: numpy.ndarray, phase: numpy.ndarray, strain: numpy.ndarray, planar_distance: float, parameters: Dict[str, Any], voxel_sizes: List[float], **kwargs)

Process the strain, modulus and phase for visualization.

flatten_sample_circles(setup: Setup) None

Send all sample circles to zero degrees.

Arrays are rotated such that all circles of the sample stage are at their zero position.

get_bulk(method: str = 'threshold') numpy.ndarray

Calculate the support representing the crystal without the surface layer.

property norm_of_q: float

Calculate the norm of q in 1/A (planar distance in nm).

rescale_q() None

Multiply the normalized diffusion vector by its original norm.

rotate_crystal(axis_to_align: List[float], reference_axis: numpy.ndarray) None

Rotate the modulus, phase and strain along the reference axis.

rotate_vector_to_saving_frame(vector: numpy.ndarray, axis_to_align: numpy.ndarray, reference_axis: numpy.ndarray) None

Calculate the diffusion vector in the crystal frame.

class bcdi.postprocessing.analysis.OrthogonalFrame(scan_index: int, parameters: Dict[str, Any], setup: Setup, **kwargs)

Analysis workflow for data already in an orthogonal frame.

calculate_voxel_sizes() List[float]

Calculate the direct space voxel sizes based on loaded q values.

get_q_bragg_crystal_frame() numpy.ndarray | None

Return the q vector of the Bragg peak in the crystal frame.

interpolate_into_crystal_frame() None

Regrid and rotate the data if necessary.

property q_values: List[numpy.ndarray] | None

List of three 1D arrays [qx, qz, qy].

regrid(new_voxelsizes: List[float]) None

Regrid the data based on user-defined voxel sizes.

class bcdi.postprocessing.analysis.PhaseManipulator(data: numpy.ndarray, parameters: Dict[str, Any], original_shape: Tuple[int, ...], wavelength: float, save_directory: str | None = None, **kwargs)

Process the phase of the data.

add_ramp(sign: int = 1) None

Add a linear ramp to the phase.

apodize() None

Apply a filtering window to the phase.

average_phase() None

Apply an averaging window to the phase.

center_phase() None

Wrap the phase around its mean.

compensate_refraction(optical_path: numpy.ndarray) None

Compensate the phase shift due to refraction through the crystal medium.

extract_phase_modulus() Tuple[numpy.ndarray, numpy.ndarray]

Get the phase and the modulus out of the data.

remove_offset() None

Remove a phase offset to the phase.

remove_ramp() None

Remove the linear trend in the phase based on a modulus support.

bcdi.postprocessing.analysis.create_analysis(scan_index: int, parameters: Dict[str, Any], setup: Setup, **kwargs) Analysis

Create the correct analysis class depending on the parameters.

bcdi.postprocessing.analysis.define_analysis_type(data_frame: str) str

Define the correct analysis type depending on the parameters.

facet_analysis

This module provides tools to import and stores the output of the facet analyzer plugin from Paraview for further analysis. See Ultramicroscopy 122, 65-75 (2012) https://doi.org/10.1016/j.ultramic.2012.07.024

One can extract the strain component and the displacement on the facets, and retrieve the correct facet normals based on a user input (geometric transformation into the crystal frame). Only cubic lattices are supported.

API Reference

Postprocessing of the output from the facet analyzer plugin for Paraview.

class bcdi.postprocessing.facet_analysis.Facets(filename: str, pathdir: str = './', savedir: str | None = None, lattice: float = 3.912, **kwargs)

Import and stores data output of facet analyzer plugin for further analysis.

Extract the mean and standard deviation of the strain and displacement distribution on the facets. Retrieves the correct facet normals based on a user input (geometric transformation into the crystal frame). It requires as input a VTK file extracted from the FacetAnalyser plugin from ParaView. See: https://doi.org/10.1016/j.ultramic.2012.07.024

Original tutorial on how to open vtk files: http://forrestbao.blogspot.com/2011/12/reading-vtk-files-in-python-via-python.html

Expected directory structure:
  • vtk file should have been saved in in Sxxxx/postprocessing

  • the analysis output will be saved in Sxxxx/postprocessing/facet_analysis

Several plotting options are attributes of this class, feel free to change them (cmap, strain_range, disp_range_avg, disp_range, strain_range_avg, comment, title_fontsize, axes_fontsize, legend_fontsize, ticks_fontsize)

Parameters:
  • filename – str, name of the VTK file

  • pathdir – str, path to the VTK file

  • savedir – str, path where to save results. If None, they will be saved in pathdir/facets_analysis

  • lattice – float, atomic spacing of the material in angstroms (only cubic lattices are supported).

evolution_curves(ncol: int = 1) None

Plot strain and displacement evolution for each facet.

Parameters:

ncol – number of columns in the plot

extract_facet(facet_id: int, plot: bool = False, elev: int = 0, azim: int = 0, output: bool = True, save: bool = True) Dict | None

Extract data from one facet.

It extracts the facet direction [x, y, z], the strain component, the displacement and their means, and also plots it.

Parameters:
  • facet_id – id of facet in paraview

  • plot – True to see plots:

  • elev – elevation angle in the z plane (in degrees).

  • azim – azimuth angle in the (x, y) plane (in degrees).

  • output – True to return facet data

  • save – True to save plot

property filename

Path to the VTK file.

fixed_reference(hkl_reference: Tuple[float, float, float] = (1, 1, 1), plot: bool = True) None

Compute the interplanar angles between each normal and a fixed reference vector.

Parameters:
  • hkl_reference – tuple of three real numbers, reference crystallographic direction

  • plot – True to see plots

load_vtk() None

Load the VTK file.

In paraview, the facets have an index that starts at 1, the index 0 corresponds to the edges and corners of the facets.

plot_displacement(figsize: Tuple[float, float] = (12, 10), elev: int = 0, azim: int = 0, save: bool = True) None

Plot two views of the surface dispalcement of the nanocrystal.

The first one with the surface coloured by the mean displacement per facet. The second one with the surface coloured by the displacement per voxel.

Parameters:
  • figsize – figure size in inches (width, height)

  • elev – elevation angle in the z plane (in degrees).

  • azim – azimuth angle in the (x, y) plane (in degrees).

  • save – True to save the figures

plot_strain(figsize: Tuple[float, float] = (12, 10), elev: int = 0, azim: int = 0, save: bool = True) None

Plot two views of the surface strain of the nanocrystal.

The first one with the surface coloured by the mean strain per facet. The second one with the surface coloured by the strain per voxel.

Parameters:
  • figsize – figure size in inches (width, height)

  • elev – elevation angle in the z plane (in degrees).

  • azim – azimuth angle in the (x, y) plane (in degrees).

  • save – True to save the figures

rotate_particle()

Rotate the nanocrystal.

The rotation is so that the base of the normals to the facets is computed with the new rotation matrix.

save_data(path_to_data: str) None

Save the field data as a csv file.

Parameters:

path_to_data – path where to save the data

save_edges_corners_data() None

Extract the edges and corners data, i.e. the mean strain and displacement.

set_rotation_matrix(u0: numpy.ndarray, v0: numpy.ndarray, u: numpy.ndarray, v: numpy.ndarray) None

Define the rotation matrix.

u and v should be the vectors perpendicular to two facets. The rotation matrix is then used if the argument rotate_particle is set to True in the method load_vtk.

Parameters:
  • u0 – numpy.ndarray, shape (3,)

  • v0 – numpy.ndarray, shape (3,)

  • u – numpy.ndarray, shape (3,)

  • v – numpy.ndarray, shape (3,)

test_rotation_matrix(vec: numpy.ndarray) None

Computes value of a vector passed through the rotation matrix.

Parameters:

vec – numpy ndarray of shape (3,). e.g. np.array([-0.833238, -0.418199, -0.300809])

to_hdf5(path_to_data: str) None

Save the facets object as an hdf5 file.

Can be combined with the file saved by the gwaihir gui.

Parameters:

path_to_data – path where to save the data

view_particle(facet_id_range: Tuple[int, int], elev_axis: str, show_edges_corners: bool, elev: int = 0, azim: int = 0) None

Visualization of the nanocrystal.

x, y and z correspond to the frame used in paraview before saving the facet analyser plugin data.

Parameters:
  • elev – elevation angle in the z plane (in degrees).

  • azim – azimuth angle in the (x, y) plane (in degrees).

  • facet_id_range – tuple of two facets numbers, facets with numbers between these two values will be plotted (higher boundary is excluded)

  • elev_axis – “x”, “y” or “z”

  • show_edges_corners – set it to True to plot also edges and corners

postprocessing_runner

This module provides the function which manage the whole postprocessing.

API Reference

Main runner for BCDI data postprocessing, after phase retrieval.

bcdi.postprocessing.postprocessing_runner.initialize_parameters(parameters: Dict[str, Any]) Dict[str, Any]

Configure and validate the existing dictionary of parameters.

bcdi.postprocessing.postprocessing_runner.run(prm: Dict[str, Any], procedure: str = 'strain_computation') None

Run the postprocessing defined by the configuration parameters.

It assumes that the dictionary of parameters was validated via a ConfigChecker instance.

Parameters:
  • prm – the parsed parameters

  • procedure – “orthogonalization” to do only the interpolation, “strain_computation” to use the full workflow

postprocessing_utils

This module provides methods used for post-processing data after phase retrieval. For example (but not limited to): phase offset and ramp removal, centering, cropping, padding, aligning reconstructions, filtering…

API Reference

Functions related to data postprocessing after phase retrieval.

exception bcdi.postprocessing.postprocessing_utils.ModulusBelowThreshold

Custom exception when there are no voxels above a threshold in an array.

bcdi.postprocessing.postprocessing_utils.apodize(amp, phase, initial_shape, window_type, debugging=False, **kwargs)

Apodize the complex array based on the window of the same shape.

Parameters:
  • amp – 3D array, amplitude before apodization

  • phase – 3D array, phase before apodization

  • initial_shape – shape of the FFT used for phasing

  • window_type – window filtering function, ‘normal’ or ‘tukey’ or ‘blackman’

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

  • kwargs

    • ‘cmap’: str, name of the colormap

    • for the normal distribution: ‘sigma’ and ‘mu’ of the 3d multivariate normal distribution, tuples of 3 floats

    • for the Tuckey window: ‘alpha’ (shape parameter) of the 3d Tukey window, tuple of 3 floats

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

    • ’logger’: an optional logger

Returns:

filtered amplitude, phase of the same shape as myamp

bcdi.postprocessing.postprocessing_utils.blackman_window(shape, normalization=1)

Create a 3d Blackman window based on shape.

Parameters:
  • shape – tuple, shape of the 3d window

  • normalization – value of the integral of the backman window

Returns:

the 3d Blackman window

bcdi.postprocessing.postprocessing_utils.bragg_temperature(spacing, reflection, spacing_ref=None, temperature_ref=None, use_q=False, material=None)

Calculate the temperature from Bragg peak position.

Parameters:
  • spacing – q or planar distance, in inverse angstroms or angstroms

  • reflection – measured reflection, e.g. np.array([1, 1, 1])

  • spacing_ref – reference spacing at known temperature (include substrate-induced strain)

  • temperature_ref – in K, known temperature for the reference spacing

  • use_q (bool) – set to True to use q, False to use planar distance

  • material – at the moment only ‘Pt’

Returns:

calculated temperature

bcdi.postprocessing.postprocessing_utils.calc_coordination(support, kernel=numpy.ones, width_z=None, width_y=None, width_x=None, debugging=False)

Calculate the coordination number of voxels in a support (numbe of neighbours).

Parameters:
  • support – 3D support array

  • kernel – kernel used for convolution with the support

  • 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

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

Returns:

the coordination matrix

bcdi.postprocessing.postprocessing_utils.center_com(array: numpy.ndarray, debugging: bool = False, **kwargs) numpy.ndarray

Center array based on center_of_mass(abs(array)) using pixel shift.

Parameters:
  • array – 3D array to be centered using the center of mass of abs(array)

  • debugging – boolean, True to see plots

  • kwargs

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

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

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

Returns:

array centered by pixel shift

bcdi.postprocessing.postprocessing_utils.center_max(array: numpy.ndarray, debugging: bool = False, **kwargs) numpy.ndarray

Center array based on max(abs(array)) using pixel shift.

Parameters:
  • array – 3D array to be centered based on max(abs(array))

  • debugging – boolean, True to see plots

  • kwargs

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

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

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

Returns:

array centered by pixel shift

bcdi.postprocessing.postprocessing_utils.center_object(method: str, obj: numpy.ndarray) numpy.ndarray

Center an object based on its modulus.

Parameters:
  • method – method used for centering. ‘com’ (center of mass), ‘max’, ‘max_com’ (max then com), ‘skip’.

  • obj – a 3D numpy array

Returns:

the centered array of the same shape as obj

bcdi.postprocessing.postprocessing_utils.filter_3d(array, filter_name='gaussian_highpass', kernel_length=21, debugging=False, **kwargs)

Apply a filter to the array by convoluting with a filtering kernel.

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

  • filter_name – name of the filter, ‘gaussian’, ‘gaussian_highpass’

  • kernel_length – length in pixels of the filtering kernel

  • debugging – True to see a plot of the kernel

  • kwargs

    • ‘sigma’: sigma of the gaussian kernel

bcdi.postprocessing.postprocessing_utils.find_bulk(amp: numpy.ndarray, support_threshold: float, method: str = 'threshold', width_z: int | None = None, width_y: int | None = None, width_x: int | None = None, debugging: bool = False, **kwargs) numpy.ndarray

Isolate the inner part of the crystal from the non-physical surface.

Parameters:
  • amp – 3D array, reconstructed object amplitude

  • support_threshold – threshold for isosurface determination

  • method – ‘threshold’ or ‘defect’. If ‘defect’, removes layer by layer using the coordination number.

  • 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

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

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’logger’: an optional logger

Returns:

the support corresponding to the bulk

bcdi.postprocessing.postprocessing_utils.find_datarange(array: numpy.ndarray, plot_margin: int | List[int] = 10, amplitude_threshold: float = 0.1, keep_size: bool = False) List[int]

Find the range where data is larger than a threshold.

It finds the meaningful range of the array where it is larger than the threshold, in order to later crop the array to that shape and reduce the memory consumption in processing. The range can be larger than the initial data size, which then will need to be padded.

Parameters:
  • array – a non-empty numpy array

  • plot_margin – user-defined margin to add on each side of the thresholded array

  • amplitude_threshold – threshold used to define a support from the amplitude

  • keep_size – set to True in order to keep the dataset full size

Returns:

a list of the ranges (centered in the middle of the array) to use in each dimension.

bcdi.postprocessing.postprocessing_utils.flip_reconstruction(obj, debugging=False, **kwargs)

Calculate the conjugate object giving the same diffracted intensity as ‘obj’.

Parameters:
  • obj – 3D reconstructed complex object

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

  • kwargs

    • ‘cmap’: str, name of the colormap

Returns:

the flipped complex object

bcdi.postprocessing.postprocessing_utils.gaussian_kernel(ndim, kernel_length=21, sigma=3, debugging=False)

Generate 2D or 3D Gaussian kernels.

Parameters:
  • ndim – number of dimensions of the kernel, 2 or 3

  • kernel_length – length in pixels of the filtering kernel

  • sigma – sigma of the gaussian pdf

  • debugging – True to see plots

Returns:

a 2D or 3D Gaussian kernel

bcdi.postprocessing.postprocessing_utils.get_opticalpath(support: numpy.ndarray, direction: str, k: numpy.ndarray, voxel_size: float | Tuple[float, float, float] | None = None, debugging: bool = False, **kwargs) numpy.ndarray

Calculate the optical path for refraction/absorption corrections in the crystal.

‘k’ should be in the same basis (crystal or laboratory frame) as the data. For xrayutilities, the data is orthogonalized in crystal frame.

Parameters:
  • support – 3D array, support used for defining the object

  • direction – “in” or “out” , incident or diffracted wave

  • k – vector for the incident or diffracted wave depending on direction, expressed in an orthonormal frame (without taking in to account the different voxel size in each dimension)

  • voxel_size – tuple, actual voxel size in z, y, and x (CXI convention)

  • debugging – boolena, 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

Returns:

the optical path in nm, of the same shape as mysupport

bcdi.postprocessing.postprocessing_utils.get_strain(phase, planar_distance, voxel_size, reference_axis='y', extent_phase=6.283185307179586, method='default', debugging=False, **kwargs)

Calculate the 3D strain array.

Parameters:
  • phase – 3D phase array (do not forget the -1 sign if the phasing algorithm is python or matlab-based)

  • planar_distance – the planar distance of the material corresponding to the measured Bragg peak

  • voxel_size – float or tuple of three floats, the voxel size of the phase array in nm

  • reference_axis – the axis of the array along which q is aligned: ‘x’, ‘y’ or ‘z’ (CXI convention)

  • extent_phase – range for phase wrapping, specify it when the phase spans over more than 2*pi

  • method – ‘default’ or ‘defect’. If ‘defect’, will offset the phase in a loop and keep the smallest value for the strain (Felix Hofmann’s method 2019).

  • debugging – True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

Returns:

the strain 3D array

bcdi.postprocessing.postprocessing_utils.match_data_range_for_interpolation(old_shape: Tuple[int, int, int], old_voxelsizes: Tuple[int, int, int], new_voxelsizes: Tuple[int, int, int]) Tuple[int, int, int]

Find the shape of the array that keeps the data span constant with new voxel sizes.

Parameters:
  • old_shape – current shape of the 3D array

  • old_voxelsizes – current voxel size of the 3D array in each dimension

  • new_voxelsizes – voxel size of the 3D array in each dimension for the interpolation

Returns:

the shape that keep the data span constant with the new voxel sizes.

bcdi.postprocessing.postprocessing_utils.mean_filter(array, support, half_width=0, width_z=None, width_y=None, width_x=None, vmin=numpy.nan, vmax=numpy.nan, title='Object', debugging=False, **kwargs)

Apply a mean filter to an object defined by a support.

Only voxels belonging to the object are taken into account, taking care of the object’s surface.

Parameters:
  • array – 3D array to be averaged

  • support – support used for averaging

  • half_width – half_width of the 2D square averaging window, 0 means no averaging, 1 is one pixel away…

  • 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

  • vmin – real number, lower boundary for the colorbar of the plots

  • vmax – real number, higher boundary for the colorbar of the plots

  • title – str, title for the plots

  • debugging – bool, True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’logger’: an optional logger

Returns:

averaged array of the same shape as the input array

bcdi.postprocessing.postprocessing_utils.ortho_modes(array_stack, nb_mode=None, method='eig', verbose=False)

Decompose an object into a set of orthogonal modes.

It finds modes from a N+1 dimensional array or a list/tuple of N-dimensional arrays. The decomposition is such that the total intensity (i.e. (abs(m)**2).sum( )) is conserved. Adapted from PyNX.

param array_stack:

the stack of modes to orthogonalize along the first dimension.

param nb_mode:

the maximum number of modes to be returned. If None, all are returned. This is useful if nb_mode is used, and only a partial list of modes is returned.

param method:

either ‘eig’ to use eigenvalue decomposition or ‘svd’ to use singular value decomposition.

param verbose:

set it to True to have more printed comments

return:

an array (modes) with the same shape as given in input, but with orthogonal modes, i.e. (mo[i]*mo[j].conj()).sum()=0 for i!=j. The modes are sorted by decreasing norm. If nb_mode is not None, only modes up to nb_mode will be returned.

bcdi.postprocessing.postprocessing_utils.regrid(array, old_voxelsize, new_voxelsize)

Interpolate real space data on a grid with a different voxel size.

Parameters:
  • array – 3D array, the object to be interpolated

  • old_voxelsize – tuple, actual voxel size in z, y, and x (CXI convention)

  • new_voxelsize – tuple, desired voxel size for the interpolation in z, y, and x (CXI convention)

Returns:

obj interpolated using the new voxel sizes

bcdi.postprocessing.postprocessing_utils.remove_offset(array, support, offset_method='com', phase_offset=0, offset_origin=None, title='', debugging=False, **kwargs)

Remove the offset in a 3D array based on a 3D support.

Parameters:
  • array – a 3D array

  • support – A 3D support of the same shape as array, defining the object

  • offset_method – ‘com’ or ‘mean’. If ‘com’, the value of array at the center of mass of the support will be subtracted to the array. If ‘mean’, the mean value of array on the support will be subtracted to the array.

  • phase_offset – value to add to the array

  • offset_origin – If provided, the value of array at this voxel will be subtracted to the array.

  • title – string, used in plot title

  • debugging – True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’reciprocal_space’: True if the object is in reciprocal space

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

Returns:

the processed array

bcdi.postprocessing.postprocessing_utils.remove_ramp(amp, phase, initial_shape, width_z=None, width_y=None, width_x=None, amplitude_threshold=0.25, threshold_gradient=0.2, method='gradient', ups_factor=2, debugging=False, **kwargs)

Remove the linear trend in the ramp using its gradient and a threshold n 3D dataset.

Parameters:
  • amp – 3D array, amplitude of the object

  • phase – 3D array, phase of the object to be detrended

  • initial_shape – shape of the FFT used for phasing

  • 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

  • amplitude_threshold – threshold used to define the support of the object from the amplitude

  • threshold_gradient – higher threshold used to select valid voxels in the gradient array

  • method – ‘gradient’ or ‘upsampling’

  • ups_factor – upsampling factor (the original shape will be multiplied by this value)

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

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’logger’: an optional logger

Returns:

normalized amplitude, detrended phase, ramp along z, ramp along y, ramp along x

bcdi.postprocessing.postprocessing_utils.remove_ramp_2d(amp, phase, initial_shape, width_y=None, width_x=None, amplitude_threshold=0.25, threshold_gradient=0.2, method='gradient', ups_factor=2, debugging=False)

Remove the linear trend in the ramp using its gradient and a threshold.

This function can be used for a 2D dataset.

Parameters:
  • amp – 2D array, amplitude of the object

  • phase – 2D array, phase of the object to be detrended

  • initial_shape – shape of the FFT used for phasing

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

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

  • amplitude_threshold – threshold used to define the support of the object from the amplitude

  • threshold_gradient – higher threshold used to select valid voxels in the gradient array

  • method – ‘gradient’ or ‘upsampling’

  • ups_factor – upsampling factor (the original shape will be multiplied by this value)

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

Returns:

normalized amplitude, detrended phase, ramp along y, ramp along x

bcdi.postprocessing.postprocessing_utils.sort_reconstruction(file_path, data_range, amplitude_threshold, sort_method='variance/mean') List[int]

Sort out reconstructions based on the metric ‘sort_method’.

Parameters:
  • file_path – path of the reconstructions to sort out

  • data_range – data will be cropped or padded to this range

  • amplitude_threshold – threshold used to define a support from the amplitude

  • sort_method – method for sorting the reconstructions: ‘variance/mean’, ‘mean_amplitude’, ‘variance’ or ‘volume’

Returns:

a list of sorted indices in ‘file_path’, from the best object to the worst.

bcdi.postprocessing.postprocessing_utils.tukey_window(shape, alpha=numpy.array)

Create a 3d Tukey window based on shape and the shape parameter alpha.

Parameters:
  • shape – tuple, shape of the 3d window

  • alpha – shape parameter of the Tukey window, tuple or ndarray of 3 values

Returns:

the 3d Tukey window

bcdi.postprocessing.postprocessing_utils.unwrap(obj: numpy.ndarray, support_threshold: float, seed: int = 0, debugging: bool = True, **kwargs) Tuple[numpy.ndarray, float]

Unwrap the phase of a complex object.

It is based on skimage.restoration.unwrap_phase. A mask can be applied by thresholding the modulus of the object.

Parameters:
  • obj – array to be unwrapped

  • support_threshold – relative threshold used to define a support from abs(obj)

  • seed – int, random seed. Use always the same value if you want a deterministic behavior.

  • debugging – set to True to see plots

  • kwargs

    • ‘cmap’: str, name of the colormap

    • ’reciprocal_space’: True if the object is in reciprocal space

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

Returns:

unwrapped phase, unwrapped phase range

process_scan

This module provides the function which manage the postprocessing for a single scan.

API Reference

Workflow for BCDI data postprocessing of a single scan, after phase retrieval.

bcdi.postprocessing.process_scan.process_scan(scan_idx: int, prm: Dict[str, Any]) Tuple[Path, Path, Logger | None]

Run the postprocessing defined by the configuration parameters for a single scan.

This function is meant to be run as a process in multiprocessing, although it can also be used as a normal function for a single scan. It assumes that the dictionary of parameters was validated via a ConfigChecker instance.

Parameters:
  • scan_idx – index of the scan to be processed in prm[“scans”]

  • prm – the parsed parameters

raw_orthogonalization

This module provides the function for BCDI data orthogonalization of a single scan, after phase retrieval (skipping centering, phase offset and ramp removal).

API Reference

Workflow for BCDI data orthogonalization of a single scan, after phase retrieval.

bcdi.postprocessing.raw_orthogonalization.orthogonalize(scan_idx: int, prm: Dict[str, Any]) Tuple[Path, Path, Logger | None]

Run the orthogonalization defined by the configuration parameters for a single scan.

This function is meant to be run as a process in multiprocessing, although it can also be used as a normal function for a single scan. It assumes that the dictionary of parameters was validated via a ConfigChecker instance.

Parameters:
  • scan_idx – index of the scan to be processed in prm[“scans”]

  • prm – the parsed parameters