Skip to content

Feature extraction

The qim3d library provides a set of methods for feature extraction on volumetric data.

General usage

  • All features assume a single connected object in the input.
  • All feature functions accept either a 3D volume or a mesh as input.
  • If a volume is provided, it is typically binarized (using the provided threshold or Otsu's method by default) and converted to a mesh before feature extraction.
  • A mask can be provided to restrict feature extraction to a specific region of interest in the volume.
  • If a mesh is provided, the threshold and/or mask arguments are ignored.

Efficient feature extraction

Before extracting multiple features, convert your input volume to a mesh using qim3d.mesh.from_volume for best performance. This avoids repeated volume-to-mesh conversions under the hood during feature extraction.

qim3d.features

qim3d.features.area

area(obj, mask=None, threshold='otsu')

Compute the surface area of an object from a volume or mesh using the Pygel3D library.

Parameters:

Name Type Description Default
obj ndarray or Manifold

Input np.ndarray volume or a mesh object of type pygel3d.hmesh.Manifold.

required
mask ndarray or None

Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

None
threshold (float, str)

Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

'otsu'

Returns:

Name Type Description
area float

The surface area of the object.

Raises:

Type Description
ValueError

If the mask shape does not match the volume shape.

Example

Compute area from a np.ndarray volume:

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume()

# Compute the surface area of the object
area = qim3d.features.area(synthetic_object)
print(area)
58535.06

Example

Compute area from a pygel3d.hmesh.Manifold mesh:

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume()

# Convert into a mesh
mesh = qim3d.mesh.from_volume(synthetic_object)

# Compute the surface area of the object
area = qim3d.features.area(mesh)
print(area)
58535.06

Source code in qim3d/features/_common_features_methods.py
def area(
    obj: np.ndarray | hmesh.Manifold,
    mask: np.ndarray | None = None,
    threshold: float | str = 'otsu',
) -> float:
    """
    Compute the surface area of an object from a volume or mesh using the Pygel3D library.

    Args:
        obj (np.ndarray or hmesh.Manifold): Input `np.ndarray` volume or a mesh object of type `pygel3d.hmesh.Manifold`.
        mask (numpy.ndarray or None): Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.
        threshold (float, str): Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

    Returns:
        area (float): The surface area of the object.

    Raises:
        ValueError: If the mask shape does not match the volume shape.

    Example:
        Compute area from a `np.ndarray` volume:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume()

        # Compute the surface area of the object
        area = qim3d.features.area(synthetic_object)
        print(area)
        ```
        58535.06

    Example:
        Compute area from a `pygel3d.hmesh.Manifold` mesh:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume()

        # Convert into a mesh
        mesh = qim3d.mesh.from_volume(synthetic_object)

        # Compute the surface area of the object
        area = qim3d.features.area(mesh)
        print(area)
        ```
        58535.06

    """
    # Prepare object
    mesh = prepare_obj(obj, threshold=threshold, mask=mask, return_mesh=True)

    # Compute area
    area = hmesh.area(mesh)

    return area

qim3d.features.volume

volume(obj, mask=None, threshold='otsu')

Compute the volume of an object from a volume or mesh using the Pygel3D library.

Parameters:

Name Type Description Default
obj ndarray or Manifold

Input np.ndarray volume or a mesh object of type pygel3d.hmesh.Manifold.

required
mask ndarray or None

Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

None
threshold (float, str)

Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

'otsu'

Returns:

Name Type Description
volume float

The volume of the object.

Raises:

Type Description
ValueError

If the mask shape does not match the volume shape.

Example

import qim3d
import numpy as np

# Generate a synthetic object
synthetic_object = qim3d.generate.volume(noise_scale=0.01, final_shape=(100, 100, 100))

# Create a mask for the bottom right corner
mask = np.zeros_like(synthetic_object, dtype=bool)
mask[50:100, 50:100, 50:100] = True

# Compute the volume of the object within the region of interest defined by the mask
volume = qim3d.features.volume(synthetic_object, threshold=50, mask=mask)
print(volume)
48774.99

Source code in qim3d/features/_common_features_methods.py
def volume(
    obj: np.ndarray | hmesh.Manifold,
    mask: np.ndarray | None = None,
    threshold: float | str = 'otsu',
) -> float:
    """
    Compute the volume of an object from a volume or mesh using the Pygel3D library.

    Args:
        obj (np.ndarray or hmesh.Manifold): Input `np.ndarray` volume or a mesh object of type `pygel3d.hmesh.Manifold`.
        mask (numpy.ndarray or None): Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.
        threshold (float, str): Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

    Returns:
        volume (float): The volume of the object.

    Raises:
        ValueError: If the mask shape does not match the volume shape.

    Example:
        ```python
        import qim3d
        import numpy as np

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume(noise_scale=0.01, final_shape=(100, 100, 100))

        # Create a mask for the bottom right corner
        mask = np.zeros_like(synthetic_object, dtype=bool)
        mask[50:100, 50:100, 50:100] = True

        # Compute the volume of the object within the region of interest defined by the mask
        volume = qim3d.features.volume(synthetic_object, threshold=50, mask=mask)
        print(volume)
        ```
        48774.99

    """
    # Prepare object
    mesh = prepare_obj(obj, threshold=threshold, mask=mask, return_mesh=True)

    # Compute volume
    volume = hmesh.volume(mesh)

    return volume

qim3d.features.size

size(obj, mask=None, threshold='otsu')

Compute the size (maximum side length of the bounding box enclosing the object) of an object from a volume or mesh.

Parameters:

Name Type Description Default
obj ndarray or Manifold

Input np.ndarray volume or a mesh object of type pygel3d.hmesh.Manifold.

required
mask ndarray or None

Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

None
threshold (float, str)

Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

'otsu'

Returns:

Name Type Description
size float

The size of the object, defined as the maximum side length of the bounding box enclosing the object.

Raises:

Type Description
ValueError

If the mask shape does not match the volume shape.

Example

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume(
    final_shape=(100,30,30),
    noise_scale=0.008,
    shape="cylinder"
    )

# Compute size of the object
size = qim3d.features.size(synthetic_object)
print(f"Size: {size}")

# Visualize the synthetic object
qim3d.viz.volumetric(synthetic_object)
Size: 100.0

Source code in qim3d/features/_common_features_methods.py
def size(
    obj: np.ndarray | hmesh.Manifold,
    mask: np.ndarray | None = None,
    threshold: float | str = 'otsu',
) -> float:
    """
    Compute the size (maximum side length of the bounding box enclosing the object) of an object from a volume or mesh.

    Args:
        obj (np.ndarray or hmesh.Manifold): Input `np.ndarray` volume or a mesh object of type `pygel3d.hmesh.Manifold`.
        mask (numpy.ndarray or None): Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.
        threshold (float, str): Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

    Returns:
        size: The size of the object, defined as the maximum side length of the bounding box enclosing the object.

    Raises:
        ValueError: If the mask shape does not match the volume shape.

    Example:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume(
            final_shape=(100,30,30),
            noise_scale=0.008,
            shape="cylinder"
            )

        # Compute size of the object
        size = qim3d.features.size(synthetic_object)
        print(f"Size: {size}")

        # Visualize the synthetic object
        qim3d.viz.volumetric(synthetic_object)
        ```
        Size: 100.0
        <iframe src="https://platform.qim.dk/k3d/size_feature_example.html" width="100%" height="500" frameborder="0"></iframe>

    """
    # Prepare object
    mesh = prepare_obj(obj, threshold=threshold, mask=mask, return_mesh=True)

    # Min and max corners of the bounding box
    bbox = hmesh.bbox(mesh)
    mins, maxs = bbox

    # Maximum side length of the bounding box
    side_lengths = maxs - mins
    size = np.max(side_lengths)

    return size

qim3d.features.sphericity

sphericity(obj, mask=None, threshold='otsu')

Compute the sphericity of an object from a volume or mesh.

Parameters:

Name Type Description Default
obj ndarray or Manifold

Input np.ndarray volume or a mesh object of type pygel3d.hmesh.Manifold.

required
mask ndarray or None

Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

None
threshold (float, str)

Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

'otsu'

Returns:

Name Type Description
sphericity float

The sphericity of the object. Higher values indicate a more spherical shape.

Raises:

Type Description
ValueError

If the mask shape does not match the volume shape.

Example

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume(noise_scale=0.005)

# Compute the sphericity of the object
sphericity = qim3d.features.sphericity(synthetic_object)
print(f"Sphericity: {sphericity:.4f}")

# Visualize the synthetic object
qim3d.viz.volumetric(synthetic_object)
Sphericity: 0.9058

Example

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume(noise_scale=0.008)

# Manipulate the object
synthetic_object = qim3d.operations.stretch(synthetic_object, z_stretch=50)
synthetic_object = qim3d.operations.curve_warp(synthetic_object, x_amp=10, x_periods=4)

# Compute the sphericity of the object
sphericity = qim3d.features.sphericity(synthetic_object)
print(f"Sphericity: {sphericity:.4f}")

# Visualize the synthetic object
qim3d.viz.volumetric(synthetic_object)
Sphericity: 0.6876

Source code in qim3d/features/_common_features_methods.py
def sphericity(
    obj: np.ndarray | hmesh.Manifold,
    mask: np.ndarray | None = None,
    threshold: float | str = 'otsu',
) -> float:
    """
    Compute the sphericity of an object from a volume or mesh.

    Args:
        obj (np.ndarray or hmesh.Manifold): Input `np.ndarray` volume or a mesh object of type `pygel3d.hmesh.Manifold`.
        mask (numpy.ndarray or None): Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.
        threshold (float, str): Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

    Returns:
        sphericity (float): The sphericity of the object. Higher values indicate a more spherical shape.

    Raises:
        ValueError: If the mask shape does not match the volume shape.

    Example:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume(noise_scale=0.005)

        # Compute the sphericity of the object
        sphericity = qim3d.features.sphericity(synthetic_object)
        print(f"Sphericity: {sphericity:.4f}")

        # Visualize the synthetic object
        qim3d.viz.volumetric(synthetic_object)
        ```
        Sphericity: 0.9058
        <iframe src="https://platform.qim.dk/k3d/sphericity_feature_example_2.html" width="100%" height="500" frameborder="0"></iframe>

    Example:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume(noise_scale=0.008)

        # Manipulate the object
        synthetic_object = qim3d.operations.stretch(synthetic_object, z_stretch=50)
        synthetic_object = qim3d.operations.curve_warp(synthetic_object, x_amp=10, x_periods=4)

        # Compute the sphericity of the object
        sphericity = qim3d.features.sphericity(synthetic_object)
        print(f"Sphericity: {sphericity:.4f}")

        # Visualize the synthetic object
        qim3d.viz.volumetric(synthetic_object)
        ```
        Sphericity: 0.6876
        <iframe src="https://platform.qim.dk/k3d/sphericity_feature_example_1.html" width="100%" height="500" frameborder="0"></iframe>

    """
    # Prepare object
    mesh = prepare_obj(obj, threshold=threshold, mask=mask, return_mesh=True)

    # Compute surface area and volume
    area = qim3d.features.area(mesh)
    volume = qim3d.features.volume(mesh)

    if area == 0 or volume == 0:
        log.warning('Surface area or volume is zero, sphericity is undefined.')
        return np.nan

    # Compute sphericity
    sphericity = (np.pi ** (1 / 3) * (6 * volume) ** (2 / 3)) / area

    return sphericity

qim3d.features.roughness

roughness(obj, mask=None, threshold='otsu')

Compute the roughness (ratio between surface area and volume) of an object from a volume or mesh.

Parameters:

Name Type Description Default
obj ndarray or Manifold

Input np.ndarray volume or a mesh object of type pygel3d.hmesh.Manifold.

required
mask ndarray or None

Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

None
threshold (float, str)

Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

'otsu'

Returns:

Name Type Description
roughness float

The roughness of the object, defined as the ratio between surface area and volume. Higher roughness indicates a more complex surface.

Raises:

Type Description
ValueError

If the mask shape does not match the volume shape.

Example

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume(
    base_shape=(128,128,128),
    noise_scale=0.019,
    )

# Compute the roughness of the object
roughness = qim3d.features.roughness(synthetic_object)
print(f"Roughness: {roughness:.4f}")

# Visualize the synthetic object
qim3d.viz.volumetric(synthetic_object)
Roughness: 0.1005

Example

import qim3d

# Generate a synthetic object
synthetic_object = qim3d.generate.volume(
    base_shape=(128,128,128),
    noise_scale=0.08,
    decay_rate=18,
    gamma=0.9,
    )

# Compute the roughness of the object
roughness = qim3d.features.roughness(synthetic_object)
print(f"Roughness: {roughness:.4f}")

# Visualize the synthetic object
qim3d.viz.volumetric(synthetic_object)
Roughness: 0.2534

Source code in qim3d/features/_common_features_methods.py
def roughness(
    obj: np.ndarray | hmesh.Manifold,
    mask: np.ndarray | None = None,
    threshold: float | str = 'otsu',
) -> float:
    """
    Compute the roughness (ratio between surface area and volume) of an object from a volume or mesh.

    Args:
        obj (np.ndarray or hmesh.Manifold): Input `np.ndarray` volume or a mesh object of type `pygel3d.hmesh.Manifold`.
        mask (numpy.ndarray or None): Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.
        threshold (float, str): Threshold value for binarization of the input volume. If 'otsu', Otsu's method is used. Defaults to 'otsu'.

    Returns:
        roughness (float): The roughness of the object, defined as the ratio between surface area and volume. Higher roughness indicates a more complex surface.

    Raises:
        ValueError: If the mask shape does not match the volume shape.

    Example:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume(
            base_shape=(128,128,128),
            noise_scale=0.019,
            )

        # Compute the roughness of the object
        roughness = qim3d.features.roughness(synthetic_object)
        print(f"Roughness: {roughness:.4f}")

        # Visualize the synthetic object
        qim3d.viz.volumetric(synthetic_object)
        ```
        Roughness: 0.1005
        <iframe src="https://platform.qim.dk/k3d/roughness_feature_example_v1.html" width="100%" height="500" frameborder="0"></iframe>

    Example:
        ```python
        import qim3d

        # Generate a synthetic object
        synthetic_object = qim3d.generate.volume(
            base_shape=(128,128,128),
            noise_scale=0.08,
            decay_rate=18,
            gamma=0.9,
            )

        # Compute the roughness of the object
        roughness = qim3d.features.roughness(synthetic_object)
        print(f"Roughness: {roughness:.4f}")

        # Visualize the synthetic object
        qim3d.viz.volumetric(synthetic_object)
        ```
        Roughness: 0.2534
        <iframe src="https://platform.qim.dk/k3d/roughness_feature_example_v2.html" width="100%" height="500" frameborder="0"></iframe>

    """
    # Prepare object
    mesh = prepare_obj(obj, threshold=threshold, mask=mask, return_mesh=True)

    # Compute surface area and volume
    area = qim3d.features.area(mesh)
    volume = qim3d.features.volume(mesh)

    if area == 0 or volume == 0:
        log.warning('Surface area or volume is zero, roughness is undefined.')
        return np.nan

    # Compute roughness
    roughness = area / volume

    return roughness

qim3d.features.mean_std_intensity

mean_std_intensity(volume, mask=None)

Compute the mean and standard deviation of intensities in a volume.

Parameters:

Name Type Description Default
volume ndarray

Input np.ndarray volume.

required
mask ndarray or None

Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

None

Returns:

Name Type Description
tuple tuple[float, float]

Mean and standard deviation of intensities.

Raises:

Type Description
ValueError

If the mask shape does not match the volume shape.

Note
  • The background (intensities of 0) is excluded from the computation.
  • If a mask is provided, it will only compute the mean and standard deviation for that region of interest.
Example

import qim3d

# Load a sample object
shell_object = qim3d.examples.shell_225x128x128

# Compute mean and standard deviation of intensities in the object
mean_intensity, std_intensity = qim3d.features.mean_std_intensity(shell_object)
print(f"Mean intensity: {mean_intensity:.4f}")
print(f"Standard deviation of intensity: {std_intensity:.4f}")

# Visualize slices of the object
qim3d.viz.slices_grid(shell_object, color_bar=True, color_bar_style="large")
Mean intensity: 114.6734
Standard deviation of intensity: 45.8481 mean_std_intensity_feature

Source code in qim3d/features/_common_features_methods.py
def mean_std_intensity(
    volume: np.ndarray,
    mask: np.ndarray | None = None,
) -> tuple[float, float]:
    """
    Compute the mean and standard deviation of intensities in a volume.

    Args:
        volume (numpy.ndarray): Input `np.ndarray` volume.
        mask (numpy.ndarray or None): Boolean mask to apply for a region of interest in the volume. Must match the shape of the input volume. Defaults to None.

    Returns:
        tuple: Mean and standard deviation of intensities.

    Raises:
        ValueError: If the mask shape does not match the volume shape.

    Note:
        - The background (intensities of 0) is excluded from the computation.
        - If a mask is provided, it will only compute the mean and standard deviation for that region of interest.

    Example:
        ```python
        import qim3d

        # Load a sample object
        shell_object = qim3d.examples.shell_225x128x128

        # Compute mean and standard deviation of intensities in the object
        mean_intensity, std_intensity = qim3d.features.mean_std_intensity(shell_object)
        print(f"Mean intensity: {mean_intensity:.4f}")
        print(f"Standard deviation of intensity: {std_intensity:.4f}")

        # Visualize slices of the object
        qim3d.viz.slices_grid(shell_object, color_bar=True, color_bar_style="large")
        ```
        Mean intensity: 114.6734  
        Standard deviation of intensity: 45.8481
        ![mean_std_intensity_feature](../../assets/screenshots/mean_std_intensity_feature_example.png)

    """

    # Mask the volume (if provided)
    volume = prepare_obj(volume, threshold=None, mask=mask, return_mesh=False)

    # Get only the non-zero intensities (i.e., ignoring the background)
    intensities = volume[volume > 0]

    # Compute mean and standard deviation
    mean_intensity = np.mean(intensities)
    std_intensity = np.std(intensities)

    return mean_intensity, std_intensity