feret#

Calculate feret distances for 2-D objects.

This code comes from a gist written by @VolkerH under BSD-3 License

https://gist.github.com/VolkerH/0d07d05d5cb189b56362e8ee41882abf

During testing it was discovered that sorting points prior to derivation of upper and lower convex hulls was problematic so this step was removed.

Module Contents#

Functions#

orientation(→ int)

Determine the orientation of three points as either clockwise, counter-clock-wise or colinear.

sort_coords(→ numpy.typing.NDArray)

Sort the coordinates.

hulls(→ tuple[list, list])

Graham scan to find upper and lower convex hulls of a set of 2-D points.

all_pairs(→ list[tuple[list, list]])

Given a list of 2-D points, finds all ways of sandwiching the points.

rotating_calipers(→ collections.abc.Generator)

Given a list of 2-D points, finds all ways of sandwiching the points between two parallel lines.

triangle_height(→ float)

Calculate the height of triangle formed by three points.

_min_feret_coord(→ numpy.typing.NDArray)

Calculate the coordinate opposite the apex that is prependicular to the base of the triangle.

_angle_between(→ float)

Calculate the angle between the apex and base of the triangle.

sort_clockwise(→ numpy.typing.NDArray)

Sort an array of coordinates in a clockwise order.

min_max_feret(→ dict[float, tuple[int, int], float, ...)

Given a list of 2-D points, returns the minimum and maximum feret diameters.

get_feret_from_mask(→ tuple[float, tuple[int, int], ...)

Calculate the minimum and maximum feret diameter of the foreground object of a binary mask.

get_feret_from_labelim(→ dict)

Calculate the minimum and maximum feret and coordinates of each connected component within a labelled image.

plot_feret(, plot_calipers, plot_triangle_heights, ...)

Plot upper and lower convex hulls with rotating calipers and optionally the feret distances.

Attributes#

LOGGER

feret.LOGGER#
feret.orientation(p: numpy.typing.NDArray, q: numpy.typing.NDArray, r: numpy.typing.NDArray) int#

Determine the orientation of three points as either clockwise, counter-clock-wise or colinear.

Parameters:
  • p (npt.NDArray) – First point (assumed to have a length of 2).

  • q (npt.NDArray) – Second point (assumed to have a length of 2).

  • r (npt.NDArray) – Third point (assumed to have a length of 2).

Returns:

Returns a positive value if p-q-r are clockwise, neg if counter-clock-wise, zero if colinear.

Return type:

int

feret.sort_coords(points: numpy.typing.NDArray, axis: int = 1) numpy.typing.NDArray#

Sort the coordinates.

Parameters:
  • points (npt.NDArray) – Array of coordinates.

  • axis (int) – Which axis to axis coordinates on 0 for row; 1 for columns (default).

Returns:

Array sorted by row then column.

Return type:

npt.NDArray

feret.hulls(points: numpy.typing.NDArray, axis: int = 1) tuple[list, list]#

Graham scan to find upper and lower convex hulls of a set of 2-D points.

Points should be sorted in asecnding order first.

Graham scan <https://en.wikipedia.org/wiki/Graham_scan>

Parameters:
  • points (npt.NDArray) – 2-D Array of points for the outline of an object.

  • axis (int) – Which axis to sort coordinates on 0 for row; 1 for columns (default).

Returns:

Tuple of two Numpy arrays of the original coordinates split into upper and lower hulls.

Return type:

tuple[list, list]

feret.all_pairs(points: numpy.typing.NDArray) list[tuple[list, list]]#

Given a list of 2-D points, finds all ways of sandwiching the points.

Calculates the upper and lower convex hulls and then finds all pairwise combinations between each set of points.

Parameters:

points (npt.NDArray) – Numpy array of coordinates defining the outline of an object.

Returns:

List of all pair-wise combinations of points between lower and upper hulls.

Return type:

List[tuple[int, int]]

feret.rotating_calipers(points: numpy.typing.NDArray, axis: int = 0) collections.abc.Generator#

Given a list of 2-D points, finds all ways of sandwiching the points between two parallel lines.

This yields the sequence of pairs of points touched by each pair of lines across all points around the hull of a polygon.

Rotating Calipers <https://en.wikipedia.org/wiki/Rotating_calipers>_

Parameters:
  • points (npt.NDArray) – Numpy array of coordinates defining the outline of an object.

  • axis (int) – Which axis to sort coordinates on, 0 for row (default); 1 for columns.

Returns:

Numpy array of pairs of points.

Return type:

Generator

feret.triangle_height(base1: numpy.typing.NDArray | list, base2: numpy.typing.NDArray | list, apex: numpy.typing.NDArray | list) float#

Calculate the height of triangle formed by three points.

Parameters:
  • base1 (int) – Coordinate of first base point of triangle.

  • base2 (int) – Coordinate of second base point of triangle.

  • apex (int) – Coordinate of the apex of the triangle.

Returns:

Height of the triangle.

Return type:

float

Examples

>>> min_feret([4, 0], [4, 3], [0,0])
    4.0
feret._min_feret_coord(base1: numpy.typing.NDArray, base2: numpy.typing.NDArray, apex: numpy.typing.NDArray, round_coord: bool = False) numpy.typing.NDArray#

Calculate the coordinate opposite the apex that is prependicular to the base of the triangle.

Code courtesy of @SylviaWhittle.

Parameters:
  • base1 (npt.NDArray) – Coordinates of one point on base of triangle, these are on the same side of the hull.

  • base2 (npt.NDArray) – Coordinates of second point on base of triangle, these are on the same side of the hull.

  • apex (npt.NDArray) – Coordinate of the apex of the triangle, this is on the opposite hull.

  • round_coord (bool) – Whether to round the point to the nearest NumPy index relative to the apex’s position (i.e. either floor or ceiling).

Returns:

Coordinates of the point perpendicular to the base line that is opposite the apex, this line is the minimum feret distance for acute triangles (but not scalene triangles).

Return type:

npt.NDArray

feret._angle_between(apex: numpy.typing.NDArray, b: numpy.typing.NDArray) float#

Calculate the angle between the apex and base of the triangle.

Parameters:
  • apex (npt.NDArray) – Difference between apex and base1 coordinates.

  • b (npt.NDArray) – Difference between base2 and base1 coordinates.

Returns:

The angle between the base and the apex.

Return type:

float

feret.sort_clockwise(coordinates: numpy.typing.NDArray) numpy.typing.NDArray#

Sort an array of coordinates in a clockwise order.

Parameters:

coordinates (npt.NDArray) – Unordered array of points. Typically a convex hull.

Returns:

Points ordered in a clockwise direction.

Return type:

npt.NDArray

Examples

>>> import numpy as np
>>> from topostats.measure import feret
>>> unordered = np.asarray([[0, 0], [5, 5], [0, 5], [5, 0]])
>>> feret.sort_clockwise(unordered)
array([[0, 0],

[0, 5], [5, 5], [5, 0]])

feret.min_max_feret(points: numpy.typing.NDArray, axis: int = 0, precision: int = 13) dict[float, tuple[int, int], float, tuple[int, int]]#

Given a list of 2-D points, returns the minimum and maximum feret diameters.

Feret diameter <https://en.wikipedia.org/wiki/Feret_diameter>

Parameters:
  • points (npt.NDArray) – A 2-D array of points for the outline of an object.

  • axis (int) – Which axis to sort coordinates on, 0 for row (default); 1 for columns.

  • precision (int) – Number of decimal places passed on to in_polygon() for rounding.

Returns:

Tuple of the minimum feret distance and its coordinates and the maximum feret distance and its coordinates.

Return type:

dictionary

feret.get_feret_from_mask(mask_im: numpy.typing.NDArray, axis: int = 0) tuple[float, tuple[int, int], float, tuple[int, int]]#

Calculate the minimum and maximum feret diameter of the foreground object of a binary mask.

The outline of the object is calculated and the pixel coordinates transformed to a list for calculation.

Parameters:
  • mask_im (npt.NDArray) – Binary Numpy array.

  • axis (int) – Which axis to sort coordinates on, 0 for row (default); 1 for columns.

Returns:

Returns a tuple of the minimum feret and its coordinates and the maximum feret and its coordinates.

Return type:

Tuple[float, Tuple[int, int], float, Tuple[int, int]]

feret.get_feret_from_labelim(label_image: numpy.typing.NDArray, labels: None | list | set = None, axis: int = 0) dict#

Calculate the minimum and maximum feret and coordinates of each connected component within a labelled image.

If labels is None, all labels > 0 will be analyzed.

Parameters:
  • label_image (npt.NDArray) – Numpy array with labelled connected components (integer).

  • labels (None | list) – A list of labelled objects for which to calculate.

  • axis (int) – Which axis to sort coordinates on, 0 for row (default); 1 for columns.

Returns:

Labels as keys and values are a tuple of the minimum and maximum feret distances and coordinates.

Return type:

dict

feret.plot_feret(points: numpy.typing.NDArray, axis: int = 0, plot_points: str | None = 'k', plot_hulls: tuple | None = ('g-', 'r-'), plot_calipers: str | None = 'y-', plot_triangle_heights: str | None = 'b:', plot_min_feret: str | None = 'm--', plot_max_feret: str | None = 'm--', filename: str | pathlib.Path | None = './feret.png', show: bool = False) tuple#

Plot upper and lower convex hulls with rotating calipers and optionally the feret distances.

Plot varying levels of details in constructing convex hulls and deriving the minimum and maximum feret.

For format strings see the Notes section of matplotlib.pyplot.plot <https://matplotlib.org/stable/api/_as_gen/matplotlib.pyplot.plot.html>.

Parameters:
  • points (npt.NDArray) – Points to be plotted which form the shape of interest.

  • axis (int) – Which axis to sort coordinates on, 0 for row (default); 1 for columns. (Should give the same results!).

  • plot_points (str | None) – Format string for plotting points. If ‘None’ points are not plotted.

  • plot_hulls (tuple | None) – Tuple of length 2 of format strings for plotting the convex hull, these should differe to allow distinction between hulls. If ‘None’ hulls are not plotted.

  • plot_calipers (str | None) – Format string for plotting calipers. If ‘None’ calipers are not plotted.

  • plot_triangle_heights (str | None) – Format string for plotting the triangle heights used in calulcating the minimum feret. These should cross the opposite edge perpendicularly. If ‘None’ triangle heights are not plotted.

  • plot_min_feret (str | None) – Format string for plotting the minimum feret. If ‘None’ the minimum feret is not plotted.

  • plot_max_feret (str | None) – Format string for plotting the maximum feret. If ‘None’ the maximum feret is not plotted.

  • filename (str | Path | None) – Location to save the image to.

  • show (bool) – Whether to display the image.

Returns:

Returns Matplotlib figure and axis objects.

Return type:

fig, ax

Examples

>>> from skimage import draw
>>> from topostats.measure import feret
>>> tiny_quadrilateral = np.asarray(
    [
        [0, 0, 0, 0, 0, 0],
        [0, 0, 1, 0, 0, 0],
        [0, 1, 0, 0, 1, 0],
        [0, 0, 0, 0, 0, 0],
        [0, 0, 0, 0, 0, 0],
        [0, 0, 1, 0, 0, 0],
        [0, 0, 0, 0, 0, 0],
    ],
    dtype=np.uint8)
>>> feret.plot_feret(np.argwhere(tiny_quadrilateral == 1))
>>> holo_ellipse_angled = np.zeros((8, 10), dtype=np.uint8)
    rr, cc = draw.ellipse_perimeter(4, 5, 1, 3, orientation=np.deg2rad(30))
    holo_ellipse_angled[rr, cc] = 1
>>> feret.plot_feret(np.argwhere(holo_ellipse_angled == 1), plot_heights = None)
>>> another_triangle = np.asarray([[5, 4], [2, 1], [8,2]])
>>> feret.plot_feret(another_triangle)