quantimpy.minkowski module¶
Compute the Minkowski functionals and functions
This module can compute both the Minkowski functionals and functions for 2D and 3D Numpy 1 arrays. These computations can handle both isotropic and anisotropic image resolutions.
Notes
More information about the used algorithm can be found in the book “Statistical analysis of microstructures in materials science” by Joachim Ohser and Frank Mücklich 2.
- quantimpy.minkowski.functionals(image, res, norm)¶
Compute the Minkowski functionals in 2D or 3D.
This function computes the Minkowski functionals for the Numpy array image. Both 2D and 3D arrays are supported. Optionally, the (anisotropic) resolution of the array can be provided using the Numpy array res. When a resolution array is provided it needs to be of the same dimension as the image array.
- Parameters
image (ndarray, bool) – Image can be either a 2D or 3D array of data type bool.
res (ndarray, {int, float}, optional) – By default the resolution is assumed to be 1 <unit of length>/pixel in all directions. If a resolution is provided it needs to be of the same dimension as the image array.
norm (bool, defaults to False) – When norm=True the functionals are normalized with the total area or volume of the image. Defaults to norm=False.
- Returns
out – In the case of a 2D image this function returns a tuple of the area, length, and the Euler characteristic. In the case of a 3D image this function returns a tuple of the volume, surface, curvature, and the Euler characteristic. The return data type is float.
- Return type
tuple, float
See also
Notes
The definition of the Minkowski functionals follows the convention in the physics literature 3.
Considering a 2D body, \(X\), with a smooth boundary, \(\delta X\), the following functionals are computed:
\[\begin{split}M_{0} (X) &= \int_{X} d s, \\ M_{1} (X) &= \frac{1}{2 \pi} \int_{\delta X} d c, \text{ and } \\ M_{2} (X) &= \frac{1}{2 \pi^{2}} \int_{\delta X} \left[\frac{1}{R} \right] d c,\end{split}\]where \(d s\) is a surface element and \(d c\) is a circumference element. \(R\) is the radius of the local curvature. This results in the following definitions for the surface area, \(S = M_{0} (X)\), circumference, \(C = 2 \pi M_{1} (X)\), and the 2D Euler characteristic, \(\chi (X) = \pi M_{2} (X)\).
Considering a 3D body, \(X\), with a smooth boundary surface, \(\delta X\), the following functionals are computed:
\[\begin{split}M_{0} (X) &= V = \int_{X} d v, \\ M_{1} (X) &= \frac{1}{8} \int_{\delta X} d s, \\ M_{2} (X) &= \frac{1}{2 \pi^{2}} \int_{\delta X} \frac{1}{2} \left[\frac{1}{R_{1}} + \frac{1}{R_{2}}\right] d s, \text{ and } \\ M_{3} (X) &= \frac{3}{(4 \pi)^{2}} \int_{\delta X} \left[\frac{1}{R_{1} R_{2}}\right] d s,\end{split}\]where \(d v\) is a volume element and \(d s\) is a surface element. \(R_{1}\) and \(R_{2}\) are the principal radii of curvature of surface element \(d s\). This results in the following definitions for the volume, \(V = M_{0} (X)\), surface area, \(S = 8 M_{1} (X)\), integral mean curvature, \(H = 2 \pi^{2} M_{2} (X)\), and the 3D Euler characteristic, \(\chi (X) = 4 \pi/3 M_{3} (X)\).
Examples
These examples use the scikit-image Python package 4 and the Matplotlib Python package 5. For a 2D image the Minkowski functionals can be computed using the following example:
import numpy as np import matplotlib.pyplot as plt from skimage.morphology import (disk) from quantimpy import minkowski as mk image = np.zeros([128,128],dtype=bool) image[16:113,16:113] = disk(48,dtype=bool) plt.gray() plt.imshow(image[:,:]) plt.show() minkowski = mk.functionals(image) print(minkowski) # Compute Minkowski functionals for image with anisotropic resolution res = np.array([2, 1]) minkowski = mk.functionals(image,res) print(minkowski)
For a 3D image the Minkowski functionals can be computed using the following example:
import numpy as np import matplotlib.pyplot as plt from skimage.morphology import (ball) from quantimpy import minkowski as mk image = np.zeros([128,128,128],dtype=bool) image[16:113,16:113,16:113] = ball(48,dtype=bool) plt.gray() plt.imshow(image[:,:,64]) plt.show() minkowski = mk.functionals(image) print(minkowski) # Compute Minkowski functionals for image with anisotropic resolution res = np.array([2, 1, 3]) minkowski = mk.functionals(image,res) print(minkowski)
- quantimpy.minkowski.functions_close(closing, res, norm)¶
Compute the Minkowski functions in 2D or 3D.
This function computes the Minkowski functionals as function of the grayscale values in the Numpy array closing. Both 2D and 3D arrays are supported. Optionally, the (anisotropic) resolution of the array can be provided using the Numpy array res. When a resolution array is provided it needs to be of the same dimension as the ‘closing’ array.
The algorithm iterates over all grayscale values present in the array, starting at the largest value (white). For every grayscale value the array is converted into a binary image where values larger than the grayscale value become one (white) and all other values become zero (black). For each of these binary images the minkowski functionals are computed according to the function
functionals()
.This function can be used in combination with the
morphology()
module to compute the Minkowski functions of different morphological distance maps.- Parameters
closing (ndarray, float) – Closing can be either a 2D or 3D array of data type float.
res (ndarray, {int, float}, optional) – By default the resolution is assumed to be 1 <unit of length>/pixel in all directions. If a resolution is provided it needs to be of the same dimension as the image array.
norm (bool, defaults to False) – When norm=True the functions are normalized with the total area or volume of the image. Defaults to norm=False.
- Returns
out – In the case of a 2D image this function returns a tuple of Numpy arrays consisting of the distance (assuming one grayscale value is used per unit of length), area, length, and the Euler characteristic. In the case of a 3D image this function returns a tuple of Numpy arrays consistenting of the distance, volume, surface, curvature, and the Euler characteristic. The return data type is float.
- Return type
tuple, ndarray, float
See also
Examples
These examples use the scikit-image Python package 4 and the Matplotlib Python package 5. For a 2D image the Minkowski functions can be computed using the following example:
import numpy as np import matplotlib.pyplot as plt from skimage.morphology import (disk) from quantimpy import morphology as mp from quantimpy import minkowski as mk image = np.zeros([128,128],dtype=bool) image[16:113,16:113] = disk(48,dtype=bool) erosion_map = mp.erode_map(image) plt.gray() plt.imshow(image[:,:]) plt.show() plt.gray() plt.imshow(erosion_map[:,:]) plt.show() dist, area, length, euler = mk.functions_close(erosion_map) plt.plot(dist,area) plt.show() plt.plot(dist,length) plt.show() plt.plot(dist,euler) plt.show()
For a 3D image the Minkowski functionals can be computed using the following example:
import numpy as np import matplotlib.pyplot as plt from skimage.morphology import (ball) from quantimpy import morphology as mp from quantimpy import minkowski as mk image = np.zeros([128,128,128],dtype=bool) image[16:113,16:113,16:113] = ball(48,dtype=bool) erosion_map = mp.erode_map(image) plt.gray() plt.imshow(image[:,:,64]) plt.show() plt.gray() plt.imshow(erosion_map[:,:,64]) plt.show() dist, volume, surface, curvature, euler = mk.functions_close(erosion_map) plt.plot(dist,volume) plt.show() plt.plot(dist,surface) plt.show() plt.plot(dist,curvature) plt.show() plt.plot(dist,euler) plt.show()
- quantimpy.minkowski.functions_open(opening, res, norm)¶
Compute the Minkowski functions in 2D or 3D.
This function computes the Minkowski functionals as function of the grayscale values in the Numpy array opening. Both 2D and 3D arrays are supported. Optionally, the (anisotropic) resolution of the array can be provided using the Numpy array res. When a resolution array is provided it needs to be of the same dimension as the ‘opening’ array.
The algorithm iterates over all grayscale values present in the array, starting at the smallest value (black). For every grayscale value the array is converted into a binary image where values larger than the grayscale value become one (white) and all other values become zero (black). For each of these binary images the minkowski functionals are computed according to the function
functionals()
.This function can be used in combination with the
morphology()
module to compute the Minkowski functions of different morphological distance maps.- Parameters
opening (ndarray, float) – Opening can be either a 2D or 3D array of data type float.
res (ndarray, {int, float}, optional) – By default the resolution is assumed to be 1 <unit of length>/pixel in all directions. If a resolution is provided it needs to be of the same dimension as the image array.
norm (bool, defaults to False) – When norm=True the functions are normalized with the total area or volume of the image. Defaults to norm=False.
- Returns
out – In the case of a 2D image this function returns a tuple of Numpy arrays consisting of the distance (assuming one grayscale value is used per unit of length), area, length, and the Euler characteristic. In the case of a 3D image this function returns a tuple of Numpy arrays consistenting of the distance, volume, surface, curvature, and the Euler characteristic. The return data type is float.
- Return type
tuple, ndarray, float
See also
Examples
These examples use the Skimage Python package 4 and the Matplotlib Python package 5. For a 2D image the Minkowski functions can be computed using the following example:
import numpy as np import matplotlib.pyplot as plt from skimage.morphology import (disk) from quantimpy import morphology as mp from quantimpy import minkowski as mk image = np.zeros([128,128],dtype=bool) image[16:113,16:113] = disk(48,dtype=bool) erosion_map = mp.erode_map(image) plt.gray() plt.imshow(image[:,:]) plt.show() plt.gray() plt.imshow(erosion_map[:,:]) plt.show() dist, area, length, euler = mk.functions_open(erosion_map) plt.plot(dist,area) plt.show() plt.plot(dist,length) plt.show() plt.plot(dist,euler) plt.show()
For a 3D image the Minkowski functionals can be computed using the following example:
import numpy as np import matplotlib.pyplot as plt from skimage.morphology import (ball) from quantimpy import morphology as mp from quantimpy import minkowski as mk image = np.zeros([128,128,128],dtype=bool) image[16:113,16:113,16:113] = ball(48,dtype=bool) erosion_map = mp.erode_map(image) plt.gray() plt.imshow(image[:,:,64]) plt.show() plt.gray() plt.imshow(erosion_map[:,:,64]) plt.show() dist, volume, surface, curvature, euler = mk.functions_open(erosion_map) plt.plot(dist,volume) plt.show() plt.plot(dist,surface) plt.show() plt.plot(dist,curvature) plt.show() plt.plot(dist,euler) plt.show()
References
- 1
Charles R. Harris, K. Jarrod Millman, Stéfan J. van der Walt et al., “Array programming with NumPy”, Nature, vol. 585, pp 357-362, 2020, doi:10.1038/s41586-020-2649-2
- 2
Joachim Ohser and Frank Mücklich, “Statistical analysis of microstructures in materials science”, Wiley and Sons, New York (2000) ISBN: 0471974862
- 3
Klaus R. Mecke, “Additivity, convexity, and beyond: applications of Minkowski Functionals in statistical physics” in “Statistical Physics and Spatial Statistics”, pp 111–184, Springer (2000) doi: 10.1007/3-540-45043-2_6
- 4(1,2,3)
Stéfan van der Walt, Johannes L. Schönberger, Juan Nunez-Iglesias, François Boulogne, Joshua D. Warner, Neil Yager, Emmanuelle Gouillart, Tony Yu and the scikit-image contributors. “scikit-image: Image processing in Python.” PeerJ 2:e453 (2014) doi: 10.7717/peerj.453
- 5(1,2,3)
John D. Hunter, “Matplotlib: A 2D Graphics Environment”, Computing in Science & Engineering, vol. 9, no. 3, pp. 90-95, 2007. doi:10.1109/MCSE.2007.55