imars3d.backend.diagnostics package

iMars3D’s diagnostics module.

Submodules

imars3d.backend.diagnostics.rotation module

iMars3D: rotation center finding module.

class imars3d.backend.diagnostics.rotation.find_rotation_center(*args, **params)[source]

Bases: ParameterizedFunction

Automatically find the rotation center from a given radiograph stack.

Parameters:
  • arrays (np.ndarray) – 3D array of images, the first dimension is the rotation angle omega

  • angles (np.ndarray) – array of angles in degrees or radians, which must match the order of arrays

  • in_degrees (bool = True) – whether angles are in degrees or radians, default is True (degrees)

  • atol_deg (Optional[float] = None) – tolerance for the search of 180 deg paris, default is None (“auto”)

  • num_pairs (int = 1) – Number of pairs to look for. Specifying -1 means as many pairs as possible. The pairs will be equally spaced if possible. If the number of pairs requested is more than half of what is available, it will take the first n-piars.

  • max_workers (int = 0) – number of cores to use for parallel median filtering, default is 0, which means using all available cores.

  • tqdm_class (panel.widgets.Tqdm) – Class to be used for rendering tqdm progress

Return type:

rotation center in pixels

num_pairs = 1

imars3d.backend.diagnostics.tilt module

iMars3D’s tilt correction module.

class imars3d.backend.diagnostics.tilt.apply_tilt_correction(*args, **params)[source]

Bases: ParameterizedFunction

Apply the tilt correction to the given array.

For a 2 deg tilted rotation axis, this function will rotate each image -2 deg so that the rotation axis is upright.

Parameters:
  • arrays (np.ndarray) – The array for tilt correction

  • tilt (float) – The rotation axis tilt angle in degrees

  • max_workers (int) – Number of cores to use for parallel median filtering, default is 0, which means using all available cores.

  • tqdm_class (panel.widgets.Tqdm) – Class to be used for rendering tqdm progress

Returns:

The tilt corrected array

Return type:

np.ndarray

imars3d.backend.diagnostics.tilt.calculate_dissimilarity(tilt: float, image0: ndarray, image1: ndarray) float[source]

Calculate the dissimilarity between two images with given tilt.

Calculate the p3-norm based dissimilarity between the two images at a given tilt angle (in degrees).

Parameters:
  • tilt – The tilt angle in degrees.

  • image0 – The first image for comparison, which is often the radiograph taken at omega (< 180 deg)

  • image1 – The second image for comparison, which is often the radiograph taken at omega + 180 deg

Return type:

The p3-norm based dissimilarity between the two images

Note

1. The rotation of an image is carried out by scikit-image, which is using bilinear interpolation (order=1) by default. This introduces artifacts: the image is slightly distorted at a non-zero tilt angle, therefore perfect matching pairs with a tilted angles greater than 2.0 might return a non-zero dissimilarity.

2. In case of non-centered rotation axis, the image cropped to ensure only both images have the object in the center of FOV, however this might fail if the object is partially out of the FOV in both images as the image registration does not work for two different partials of the same object.

imars3d.backend.diagnostics.tilt.calculate_shift(reference_image: ndarray, moving_image: ndarray) float[source]

Calculate relative shift between two images.

Calculate the amount of shift needed to move the moving image to the reference image. This method requires the rotating object to be fully in the field of view.

Parameters:
  • reference_image – The reference image, often the radiograph taken at omega (< 180 deg)

  • moving_image – The moving image, often the radiograph taken at omega + 180 deg

Return type:

The amount of shift in pixels

imars3d.backend.diagnostics.tilt.calculate_tilt(image0: ndarray, image180: ndarray, low_bound: float = -5.0, high_bound: float = 5.0) OptimizeResult[source]

Use optimization to find the in-plane tilt angle.

Parameters:
  • image0 – The first image for comparison, which is often the radiograph taken at omega (< 180 deg)

  • image180 – The second image for comparison, which is often the radiograph taken at omega + 180 deg

  • low_bound – The lower bound of the tilt angle search space

  • high_bound – The upper bound of the tilt angle search space

Return type:

The optimization results from scipy.optimize.minimize_scalar

imars3d.backend.diagnostics.tilt.find_180_deg_pairs_idx(angles: ndarray, atol: float = None, in_degrees: bool = True) Tuple[ndarray, ndarray][source]

Return the indices of the 180 degree pairs from given list of angles.

Parameters:
  • angles – The list of angles as a 1d array.

  • atol – The absolute tolerance for the 180 degree pairs.

  • in_degrees – Whether the angles are in degrees or radians, default is in degrees.

Return type:

The indices of the 180 degree pairs (low_range, high_range)

class imars3d.backend.diagnostics.tilt.tilt_correction(*args, **params)[source]

Bases: ParameterizedFunction

Detect and correct the rotation axis tilt from the given radiograph stack.

Parameters:
  • arrays (np.ndarray) – The radiograph stack fro tilt correction

  • rot_angles (np.ndarray) – The list of rotation angles in radians (follow tomopy convention)

  • low_bound (float) – The lower bound of the tilt angle search space

  • high_bound (float) – The upper bound of the tilt angle search space

  • cut_off_angle_deg (float) – The angle in degrees to cut off the rotation axis tilt correction, i.e. skip applying tilt correction for tilt angles that are too small.

  • max_workers – Number of cores to use for parallel median filtering, default is 0, which means using all available cores.

  • tqdm_class (panel.widgets.Tqdm) – Class to be used for rendering tqdm progress

Returns:

The tilt corrected array

Return type:

np.ndarray