Dense regions decomposition

Functions to detect dense and bright regions (with potential clustered spots), then use gaussian simulations to correct a misdetection in these regions.

Decompose and simulate dense regions

The main function to decompose dense regions is:

It is also possible to perform the main steps of this decomposition separately:

See an example of application here.

bigfish.detection.decompose_dense(image, spots, voxel_size, spot_radius, kernel_size=None, alpha=0.5, beta=1, gamma=5)

Detect dense and bright regions with potential clustered spots and simulate a more realistic number of spots in these regions.

  1. We estimate image background with a large gaussian filter. We then remove the background from the original image to denoise it.

  2. We build a reference spot by aggregating predetected spots.

  3. We fit gaussian parameters on the reference spots.

  4. We detect dense regions to decompose.

  5. We simulate as many gaussians as possible in the candidate regions.

Parameters:
imagenp.ndarray

Image with shape (z, y, x) or (y, x).

spotsnp.ndarray

Coordinate of the spots with shape (nb_spots, 3) or (nb_spots, 2) for 3-d or 2-d images respectively.

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

spot_radiusint, float, Tuple(int, float) or List(int, float)

Radius of the spot, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same radius is applied to every dimensions.

kernel_sizeint, float, Tuple(float, int), List(float, int) or None

Standard deviation used for the gaussian kernel (one for each dimension), in pixel. If it’s a scalar, the same standard deviation is applied to every dimensions. If None, we estimate the kernel size from ‘spot_radius’, ‘voxel_size’ and ‘gamma’

alphaint or float

Intensity percentile used to compute the reference spot, between 0 and 1. The higher, the brighter are the spots simulated in the dense regions. Consequently, a high intensity score reduces the number of spots added. Default is 0.5, meaning the reference spot considered is the median spot.

betaint or float

Multiplicative factor for the intensity threshold of a dense region. Default is 1. Threshold is computed with the formula:

\[\mbox{threshold} = \beta * \mbox{max(median spot)}\]

With \(\mbox{median spot}\) the median value of all detected spot signals.

gammaint or float

Multiplicative factor use to compute the gaussian kernel size:

\[\mbox{kernel size} = \frac{\gamma * \mbox{spot radius}}{\mbox{ voxel size}}\]

We perform a large gaussian filter with such scale to estimate image background and remove it from original image. A large gamma increases the scale of the gaussian filter and smooth the estimated background. To decompose very large bright areas, a larger gamma should be set.

Returns:
spotsnp.ndarray

Coordinate of the spots detected, with shape (nb_spots, 3) or (nb_spots, 2). One coordinate per dimension (zyx or yx coordinates).

dense_regionsnp.ndarray, np.int64

Array with shape (nb_regions, 7) or (nb_regions, 6). One coordinate per dimension for the region centroid (zyx or yx coordinates), the number of RNAs detected in the region, the area of the region, its average intensity value and its index.

reference_spotnp.ndarray

Reference spot in 3-d or 2-d.

Notes

If gamma = 0 and kernel_size = None, image is not denoised.

bigfish.detection.get_dense_region(image, spots, voxel_size, spot_radius, beta=1)

Detect and filter dense and bright regions.

A candidate region has at least 2 connected pixels above a specific threshold.

Parameters:
imagenp.ndarray

Image with shape (z, y, x) or (y, x).

spotsnp.ndarray

Coordinate of the spots with shape (nb_spots, 3) or (nb_spots, 2).

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

spot_radiusint, float, Tuple(int, float) or List(int, float)

Radius of the spot, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same radius is applied to every dimensions.

betaint or float

Multiplicative factor for the intensity threshold of a dense region. Default is 1. Threshold is computed with the formula:

\[\mbox{threshold} = \beta * \mbox{max(median spot)}\]

With \(\mbox{median spot}\) the median value of all detected spot signals.

Returns:
dense_regionsnp.ndarray

Array with selected skimage.measure._regionprops._RegionProperties objects.

spots_out_regionnp.ndarray

Coordinate of the spots detected out of dense regions, with shape (nb_spots, 3) or (nb_spots, 2). One coordinate per dimension (zyx or yx coordinates).

max_sizeint

Maximum size of the regions.

bigfish.detection.simulate_gaussian_mixture(image, candidate_regions, voxel_size, sigma, amplitude=100, background=0, precomputed_gaussian=None)

Simulate as many gaussians as possible in the candidate dense regions in order to get a more realistic number of spots.

Parameters:
imagenp.ndarray

Image with shape (z, y, x) or (y, x).

candidate_regionsnp.ndarray

Array with filtered skimage.measure._regionprops._RegionProperties.

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

sigmaint, float, Tuple(int, float) or List(int, float)

Standard deviation of the gaussian, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

amplitudefloat

Amplitude of the gaussian.

backgroundfloat

Background minimum value.

precomputed_gaussianTuple[np.ndarray]

Tuple with tables of precomputed values for the erf, with shape (nb_value, 2). One table per dimension.

Returns:
spots_in_regionsnp.ndarray, np.int64

Coordinate of the spots detected inside dense regions, with shape (nb_spots, 4) or (nb_spots, 3). One coordinate per dimension (zyx or yx coordinates) plus the index of the region.

regionsnp.ndarray, np.int64

Array with shape (nb_regions, 7) or (nb_regions, 6). One coordinate per dimension for the region centroid (zyx or yx coordinates), the number of RNAs detected in the region, the area of the region, its average intensity value and its index.


Modelize a reference spot

To simulate additional spots in the dense regions it is necessary to:

  1. Build a reference spot.

  2. Modelize this reference spot by fitting gaussian parameters.

  3. Simulate gaussian signal.

bigfish.detection.build_reference_spot(image, spots, voxel_size, spot_radius, alpha=0.5)

Build a median or mean spot in 3 or 2 dimensions as reference.

Reference spot is computed from a sample of uncropped detected spots. If such sample is not possible, an empty frame is returned.

Parameters:
imagenp.ndarray

Image with shape (z, y, x) or (y, x).

spotsnp.ndarray

Coordinate of the spots with shape (nb_spots, 3) for 3-d images or (nb_spots, 2) for 2-d images.

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

spot_radiusint, float, Tuple(int, float) or List(int, float)

Radius of the spot, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same radius is applied to every dimensions.

alphaint or float

Intensity score of the reference spot, between 0 and 1. If 0, reference spot approximates the spot with the lowest intensity. If 1, reference spot approximates the brightest spot. Default is 0.5.

Returns:
reference_spotnp.ndarray

Reference spot in 3-d or 2-d.

bigfish.detection.modelize_spot(reference_spot, voxel_size, spot_radius, return_coord=False)

Fit a gaussian function on the reference spot.

Parameters:
reference_spotnp.ndarray

A 3-d or 2-d image with detected spot and shape (z, y, x) or (y, x).

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

spot_radiusint, float, Tuple(int, float) or List(int, float)

Radius of the spot, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same radius is applied to every dimensions.

return_coordbool

Return gaussian coordinates.

Returns:
parameters_fittedTuple[float]
  • mu_zfloat (optional)

    Coordinate of the gaussian center along the z axis, in pixel.

  • mu_yfloat (optional)

    Coordinate of the gaussian center along the y axis, in pixel.

  • mu_xfloat (optional)

    Coordinate of the gaussian center along the x axis, in pixel.

  • sigma_zfloat

    Standard deviation of the gaussian along the z axis, in pixel. Available only for a 3-d modelization.

  • sigma_yxfloat

    Standard deviation of the gaussian in the yx plan, in pixel.

  • amplitudefloat

    Amplitude of the gaussian.

  • backgroundfloat

    Background minimum value of the image.

bigfish.detection.precompute_erf(ndim, voxel_size, sigma, max_grid=200)

Precompute different values for the erf with a nanometer resolution.

Parameters:
ndimint

Number of dimensions to consider (2 or 3).

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

sigmaint, float, Tuple(int, float) or List(int, float)

Standard deviation of the gaussian, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

max_gridint

Maximum size of the grid on which we precompute the erf, in pixel.

Returns:
table_erfTuple[np.ndarray]

Tuple with tables of precomputed values for the erf, with shape (nb_value, 2). One table per dimension. First column is the coordinate along the table dimension. Second column is the precomputed erf value.

bigfish.detection.initialize_grid(image_spot, voxel_size, return_centroid=False)

Build a grid in nanometer to compute gaussian function values over a full volume or surface.

Parameters:
image_spotnp.ndarray

An image with detected spot and shape (z, y, x) or (y, x).

voxel_sizeint, float, Tuple(int, float) or List(int, float)

Size of a voxel, in nanometer. One value per spatial dimension (zyx or yx dimensions). If it’s a scalar, the same value is applied to every dimensions.

return_centroidbool

Compute centroid estimation of the grid.

Returns:
gridnp.ndarray, np.float32

A grid with the shape (3, z * y * x) or (2, y * x), in nanometer.

centroid_coordTuple[float]

Estimated centroid of the spot, in nanometer. One element per dimension.

bigfish.detection.gaussian_2d(grid, mu_y, mu_x, sigma_yx, voxel_size_yx, amplitude, background, precomputed=None)

Compute the gaussian function over the grid representing a surface S with shape (S_y, S_x).

Parameters:
gridnp.ndarray, np.float

Grid data to compute the gaussian function for different voxel within a surface S. In nanometer, with shape (2, S_y * S_x).

mu_yfloat

Estimated mean of the gaussian signal along y axis, in nanometer.

mu_xfloat

Estimated mean of the gaussian signal along x axis, in nanometer.

sigma_yxint or float

Standard deviation of the gaussian in the yx plan, in nanometer.

voxel_size_yxint or float

Size of a voxel in the yx plan, in nanometer.

amplitudefloat

Estimated pixel intensity of the gaussian signal.

backgroundfloat

Estimated pixel intensity of the background.

precomputedTuple[np.ndarray] or None

Tuple with tables of precomputed values for the erf, with shape (nb_value, 2). One table per dimension.

Returns:
valuesnp.ndarray, np.float

Value of each voxel within the surface S according to the 2-d gaussian parameters. Shape (S_y * S_x,).

bigfish.detection.gaussian_3d(grid, mu_z, mu_y, mu_x, sigma_z, sigma_yx, voxel_size_z, voxel_size_yx, amplitude, background, precomputed=None)

Compute the gaussian function over the grid representing a volume V with shape (V_z, V_y, V_x).

Parameters:
gridnp.ndarray, np.float

Grid data to compute the gaussian function for different voxel within a volume V. In nanometer, with shape (3, V_z * V_y * V_x).

mu_zfloat

Estimated mean of the gaussian signal along z axis, in nanometer.

mu_yfloat

Estimated mean of the gaussian signal along y axis, in nanometer.

mu_xfloat

Estimated mean of the gaussian signal along x axis, in nanometer.

sigma_zint or float

Standard deviation of the gaussian along the z axis, in nanometer.

sigma_yxint or float

Standard deviation of the gaussian in the yx plan, in nanometer.

voxel_size_zint or float

Size of a voxel along the z axis, in nanometer.

voxel_size_yxint or float

Size of a voxel in the yx plan, in nanometer.

amplitudefloat

Estimated pixel intensity of the gaussian signal.

backgroundfloat

Estimated pixel intensity of the background.

precomputedTuple[np.ndarray] or None

Tuple with tables of precomputed values for the erf, with shape (nb_value, 2). One table per dimension.

Returns:
valuesnp.ndarray, np.float

Value of each voxel within the volume V according to the 3-d gaussian parameters. Shape (V_z * V_y * V_x,).