aicsshparam package#

Subpackages#

Submodules#

aicsshparam.shparam module#

aicsshparam.shparam.get_shcoeffs(image: array, lmax: int, sigma: float = 0, compute_lcc: bool = True, alignment_2d: bool = True, make_unique: bool = False)[source]#

Compute spherical harmonics coefficients that describe an object stored as an image.

Calculates the spherical harmonics coefficients that parametrize the shape formed by the foreground set of voxels in the input image. The input image does not need to be binary and all foreground voxels (background=0) are used in the computation. Foreground voxels must form a single connected component. If you are sure that this is the case for the input image, you can set compute_lcc to False to speed up the calculation. In addition, the shape is expected to be centered in the input image.

Parameters:
  • image (ndarray) – Input image. Expected to have shape ZYX.

  • lmax (int) – Order of the spherical harmonics parametrization. The higher the order the more shape details are represented.

Returns:

  • coeffs_dict (dict) – Dictionary with the spherical harmonics coefficients and the mean square error between input and its parametrization

  • grid_rec (ndarray) – Parametric grid representing sh parametrization

  • image_ (ndarray) – Input image after pre-processing (lcc calculation, smooth and binarization).

  • mesh (vtkPolyData) – Polydata representation of image_.

  • grid_down (ndarray) – Parametric grid representing input object.

  • transform (tuple of floats) – (xc, yc, zc, angle) if alignment_2d is True or (xc, yc, zc) if alignment_2d is False. (xc, yc, zc) are the coordinates of the shape centroid after alignment; angle is the angle used to align the image

Other Parameters:
  • sigma (float, optional) – The degree of smooth to be applied to the input image, default is 0 (no smooth)

  • compute_lcc (bool, optional) – Whether to compute the largest connected component before applying the spherical harmonic parametrization, default is True. Set compute_lcc to False in case you are sure the input image contains a single connected component. It is crucial that parametrization is calculated on a single connected component object.

  • alignment_2d (bool) – Whether the image should be aligned in 2d. Default is True.

  • make_unique (bool) – Set true to make sure the alignment rotation is unique.

Notes

Alignment mode ‘2d’ allows for keeping the z axis unchanged which might be important for some applications.

Examples

import numpy as np
from aicsshparam import shparam, shtools

img = np.ones((32,32,32), dtype=np.uint8)

(coeffs, grid_rec), (image_, mesh, grid, transform) =
    shparam.get_shcoeffs(image=img, lmax=2)
mse = shtools.get_reconstruction_error(grid, grid_rec)

print('Coefficients:', coeffs)

>>> Coefficients: {'shcoeffs_L0M0C': 18.31594310878251, 'shcoeffs_L0M1C': 0.0,
'shcoeffs_L0M2C': 0.0, 'shcoeffs_L1M0C': 0.020438775421611564, 'shcoeffs_L1M1C':
-0.0030960466571801513, 'shcoeffs_L1M2C': 0.0, 'shcoeffs_L2M0C':
-0.0185688727281408, 'shcoeffs_L2M1C': -2.9925077712704384e-05,
'shcoeffs_L2M2C': -0.009087503958673892, 'shcoeffs_L0M0S': 0.0,
'shcoeffs_L0M1S': 0.0, 'shcoeffs_L0M2S': 0.0, 'shcoeffs_L1M0S': 0.0,
'shcoeffs_L1M1S': 3.799611612562637e-05, 'shcoeffs_L1M2S': 0.0,
'shcoeffs_L2M0S': 0.0, 'shcoeffs_L2M1S': 3.672543904347801e-07,
'shcoeffs_L2M2S': 0.0002230857005948496}

print('Error:', mse)

>>> Error: 2.3738182456948795

aicsshparam.shtools module#

aicsshparam.shtools.align_image_2d(image: array, alignment_channel: int | None = None, make_unique: bool = False, compute_aligned_image: bool = True)[source]#

Align a multichannel 3D image based on the channel specified by alignment_channel. The expected shape of image is (C,Z,Y,X) or (Z,Y,X).

Parameters:
  • image (np.array) – Input array of shape (C,Z,Y,X) or (Z,Y,X).

  • alignment_channel (int) – Number of channel to be used as reference for alignemnt. The alignment will be propagated to all other channels.

  • make_unique (bool) – Set true to make sure the alignment rotation is unique.

  • compute_aligned_image (bool) – Set false to only compute and return the alignment angle

Returns:

  • img_aligned (np.array) – Aligned image

  • angle (float) – Angle used for align the shape.

aicsshparam.shtools.apply_image_alignment_2d(image: array, angle: float)[source]#

Apply an existing set of alignment parameters to a multichannel 3D image. The expected shape of image is (C,Z,Y,X) or (Z,Y,X).

Parameters:
  • image (np.array) – Input array of shape (C,Z,Y,X) or (Z,Y,X).

  • angle (float) – 2D rotation angle in degrees

Returns:

img_aligned – Aligned image

Return type:

np.array

aicsshparam.shtools.convert_coeffs_dict_to_matrix(coeffs_dict, lmax=32)[source]#

Convert a dictionary of SH coefficients to a matrix of SH coefficients. The dictionary should have keys in the format “shcoeffs_L{L}M{M}{C}” where L and M are the degree and order of the coefficient and C is either “C” or “S” for cosine or sine coefficients, respectively. :Parameters: * coeffs_dict (dict) – Dictionary of SH coefficients

  • lmax (int) – Maximum degree of the SH coefficients

Returns:

coeffs – Matrix of SH coefficients

Return type:

np.array

aicsshparam.shtools.get_even_reconstruction_from_coeffs(coeffs: array, lrec: int = 0, npoints: int = 512)[source]#

Converts a set of spherical harmonic coefficients into a 3d mesh using the Fibonacci grid for generating a mesh with a more even distribution of points

Parameters:

coeffs (np.array) – Input array of spherical harmonic coefficients. These array has dimensions 2xLxM, where the first dimension is 0 for cosine-associated coefficients and 1 for sine-associated coefficients. Second and thrid dimensions represent the expansion parameters (l,m).

Returns:

mesh – Mesh that represents the input parametric grid.

Return type:

vtkPolyData

Other Parameters:
  • lrec (int, optional) – Only coefficients l<lrec will be used for creating the mesh, default is 0 meaning all coefficients available in the matrix coefficients will be used.

  • npoints (int optional) – Number of points in the initial spherical mesh

Notes

The mesh resolution is set by the size of the coefficients matrix and therefore not affected by lrec.

aicsshparam.shtools.get_even_reconstruction_from_grid(grid: array, npoints: int = 512, centroid: Tuple = (0, 0, 0))[source]#

Converts a parametric 2D grid of type (lon,lat,rad) into a 3d mesh. lon in [0,2pi], lat in [0,pi]. The method uses a spherical mesh with an even distribution of points. The even distribution is obtained via the Fibonacci grid rule.

Parameters:

grid (np.array) – Input grid where the element grid[i,j] represents the radial coordinate at longitude i*2pi/grid.shape[0] and latitude j*pi/grid.shape[1].

Returns:

mesh – Mesh that represents the input parametric grid.

Return type:

vtkPolyData

Other Parameters:
  • npoints (int) – Number of points in the initial spherical mesh

  • centroid (tuple of floats, optional) – x, y and z coordinates of the centroid where the mesh will be translated to, default is (0,0,0).

aicsshparam.shtools.get_mesh_from_image(image: array, sigma: float = 0, lcc: bool = True, translate_to_origin: bool = True)[source]#

Converts a numpy array into a vtkImageData and then into a 3d mesh using vtkContourFilter. The input is assumed to be binary and the isosurface value is set to 0.5.

Optionally the input can be pre-processed by i) extracting the largest connected component and ii) applying a gaussian smooth to it. In case smooth is used, the image is binarize using thershold 1/e.

A size threshold is applying to garantee that enough points will be used to compute the SH parametrization.

Also, points as the edge of the image are set to zero (background) to make sure the mesh forms a manifold.

Parameters:

image (np.array) – Input array where the mesh will be computed on

Returns:

  • mesh (vtkPolyData) – 3d mesh in VTK format

  • img_output (np.array) – Input image after pre-processing

  • centroid (np.array) – x, y, z coordinates of the mesh centroid

Other Parameters:
  • lcc (bool, optional) – Wheather or not to compute the mesh only on the largest connected component found in the input connected component, default is True.

  • sigma (float, optional) – The degree of smooth to be applied to the input image, default is 0 (no smooth).

  • translate_to_origin (bool, optional) – Wheather or not translate the mesh to the origin (0,0,0), default is True.

aicsshparam.shtools.get_reconstruction_error(grid_input: array, grid_rec: array)[source]#

Compute mean square error between two parametric grids. When applied to the input parametric grid and its corresponsing reconstructed version, it gives an idea of the quality of the parametrization with low values indicate good parametrization.

Parameters:
  • grid_input (np.array) – Parametric grid

  • grid_rec (np.array) – Reconstructed parametric grid

Returns:

mse – Mean square error

Return type:

float

aicsshparam.shtools.get_reconstruction_from_coeffs(coeffs: array, lrec: int = 0)[source]#

Converts a set of spherical harmonic coefficients into a 3d mesh.

Parameters:

coeffs (np.array) – Input array of spherical harmonic coefficients. These array has dimensions 2xLxM, where the first dimension is 0 for cosine-associated coefficients and 1 for sine-associated coefficients. Second and thrid dimensions represent the expansion parameters (l,m).

Returns:

mesh – Mesh that represents the input parametric grid.

Return type:

vtkPolyData

Other Parameters:

lrec (int, optional) – Degree of the reconstruction. If lrec<l, then only coefficients l<lrec will be used for creating the mesh. If lrec>l, then the mesh will be oversampled. Default is 0 meaning all coefficients available in the matrix coefficients will be used.

Notes

The mesh resolution is set by the size of the coefficients matrix and therefore not affected by lrec.

aicsshparam.shtools.get_reconstruction_from_grid(grid: array, centroid: Tuple = (0, 0, 0))[source]#

Converts a parametric 2D grid of type (lon,lat,rad) into a 3d mesh. lon in [0,2pi], lat in [0,pi].

Parameters:

grid (np.array) – Input grid where the element grid[i,j] represents the radial coordinate at longitude i*2pi/grid.shape[0] and latitude j*pi/grid.shape[1].

Returns:

mesh – Mesh that represents the input parametric grid.

Return type:

vtkPolyData

Other Parameters:

centroid (tuple of floats, optional) – x, y and z coordinates of the centroid where the mesh will be translated to, default is (0,0,0).

aicsshparam.shtools.rotate_image_2d(image: array, angle: float, interpolation_order: int = 0)[source]#

Rotate multichannel image in 2D by a given angle. The expected shape of image is (C,Z,Y,X). The rotation will be done clock-wise around the center of the image.

Parameters:
  • angle (float) – Angle in degrees

  • interpolation_order (int) – Order of interpolation used during the image rotation

Returns:

img_rot – Rotated image

Return type:

np.array

aicsshparam.shtools.save_polydata(mesh: vtkPolyData, filename: str)[source]#

Saves a mesh as a vtkPolyData file.

Parameters:
  • mesh (vtkPolyData) – Input mesh

  • filename (str) – File path where the mesh will be saved

  • output_type (vtk or ply) – Format of output polydata file

aicsshparam.shtools.update_mesh_points(mesh: vtkPolyData, x_new: array, y_new: array, z_new: array)[source]#

Updates the xyz coordinates of points in the input mesh with new coordinates provided.

Parameters:
  • mesh (vtkPolyData) – Mesh in VTK format to be updated.

  • x_new, y_new and z_new (np.array) – Array containing the new coordinates.

Returns:

mesh_updated – Mesh with updated coordinates.

Return type:

vtkPolyData

Notes

This function also re-calculate the new normal vectors for the updated mesh.

aicsshparam.shtools.voxelize_mesh(imagedata: vtkImageData, shape: Tuple, mesh: vtkPolyData, origin: List)[source]#

Voxelize a triangle mesh into an image.

Parameters:
  • imagedata (vtkImageData) – Imagedata that will be uses as support for voxelization.

  • shape (tuple) – Shape that imagedata scalars will take after voxelization.

  • mesh (vtkPolyData) – Mesh to be voxelized

  • origin (List) – xyz specifying the lower left corner of the mesh.

Returns:

img – Binary array.

Return type:

np.array

aicsshparam.shtools.voxelize_meshes(meshes: List)[source]#

List of meshes to be voxelized into an image. Usually the input corresponds to the cell membrane and nuclear shell meshes.

Parameters:

meshes (List) – List of vtkPolydatas representing the meshes to be voxelized into an image.

Returns:

  • img (np.array) – 3D image where voxels with value i represent are those found in the interior of the i-th mesh in the input list. If a voxel is interior to one or more meshes form the input list, it will take the value of the right most mesh in the list.

  • origin – Origin of the meshes in the voxelized image.

Module contents#

Top-level package for aics-shparam.

aicsshparam.get_module_version()[source]#