metrics

Extract Quantitative Information

This submodule contains functions for determining key metrics about an image. Typically these are applied to an image after applying a filter, but a few functions can be applied directly to the binary image.

porespy.metrics.chord_counts(im) Finds the length of each chord in the supplied image and returns a list of their individual sizes
porespy.metrics.chord_length_distribution(im) Determines the distribution of chord lengths in an image containing chords.
porespy.metrics.linear_density(im[, bins, …]) Determines the probability that a point lies within a certain distance of the opposite phase along a specified direction
porespy.metrics.mesh_surface_area([mesh, …]) Calculates the surface area of a meshed region
porespy.metrics.phase_fraction(im[, normed]) Calculates the number (or fraction) of each phase in an image
porespy.metrics.pore_size_distribution(im[, …]) Calculate a pore-size distribution based on the image produced by the porosimetry or local_thickness functions.
porespy.metrics.porosity(im) Calculates the porosity of an image assuming 1’s are void space and 0’s are solid phase.
porespy.metrics.porosity_profile(im, axis) Returns a porosity profile along the specified axis
porespy.metrics.props_to_image(regionprops, …) Creates an image with each region colored according the specified prop, as obtained by regionprops_3d.
porespy.metrics.props_to_DataFrame(regionprops) Returns a Pandas DataFrame containing all the scalar metrics for each region, such as volume, sphericity, and so on, calculated by regionprops_3D.
porespy.metrics.radial_density(im[, bins, …]) Computes radial density function by analyzing the histogram of voxel values in the distance transform.
porespy.metrics.region_interface_areas(…) Calculates the interfacial area between all pairs of adjecent regions
porespy.metrics.region_surface_areas(regions) Extracts the surface area of each region in a labeled image.
porespy.metrics.regionprops_3D(im) Calculates various metrics for each labeled region in a 3D image.
porespy.metrics.representative_elementary_volume(im) Calculates the porosity of the image as a function subdomain size.
porespy.metrics.two_point_correlation_bf(im) Calculates the two-point correlation function using brute-force (see Notes)
porespy.metrics.two_point_correlation_fft(im) Calculates the two-point correlation function using fourier transforms
porespy.metrics.chord_counts(im)[source]

Finds the length of each chord in the supplied image and returns a list of their individual sizes

Parameters:im (ND-array) – An image containing chords drawn in the void space.
Returns:result – A 1D array with one element for each chord, containing its length.
Return type:1D-array

Notes

The returned array can be passed to plt.hist to plot the histogram, or to sp.histogram to get the histogram data directly. Another useful function is sp.bincount which gives the number of chords of each length in a format suitable for plt.plot.

porespy.metrics.chord_length_distribution(im, bins=None, log=False, voxel_size=1, normalization='count')[source]

Determines the distribution of chord lengths in an image containing chords.

Parameters:
  • im (ND-image) –

    An image with chords drawn in the pore space, as produced by apply_chords or apply_chords_3d.

    im can be either boolean, in which case each chord will be identified using scipy.ndimage.label, or numerical values in which case it is assumed that chords have already been identifed and labeled. In both cases, the size of each chord will be computed as the number of voxels belonging to each labelled region.

  • bins (scalar or array_like) – If a scalar is given it is interpreted as the number of bins to use, and if an array is given they are used as the bins directly.
  • log (Boolean) – If true, the logarithm of the chord lengths will be used, which can make the data more clear.
  • normalization (string) –

    Indicates how to normalize the bin heights. Options are:

    ’count’ or ‘number’ - (default) This simply counts the number of chords in each bin in the normal sense of a histogram. This is the rigorous definition according to Torquato [1].

    ’length’ - This multiplies the number of chords in each bin by the chord length (i.e. bin size). The normalization scheme accounts for the fact that long chords are less frequent than shorert chords, thus giving a more balanced distribution.

  • voxel_size (scalar) – The size of a voxel side in preferred units. The default is 1, so the user can apply the scaling to the returned results after the fact.
Returns:

result – A tuple containing the following elements, which can be retrieved by attribute name:

L or logL - chord length, equivalent to bin_centers

pdf - probability density function

cdf - cumulative density function

relfreq - relative frequency chords in each bin. The sum of all bin heights is 1.0. For the cumulative relativce, use cdf which is already normalized to 1.

bin_centers - the center point of each bin

bin_edges - locations of bin divisions, including 1 more value than the number of bins

bin_widths - useful for passing to the width argument of matplotlib.pyplot.bar

Return type:

named_tuple

References

[1] Torquato, S. Random Heterogeneous Materials: Mircostructure and Macroscopic Properties. Springer, New York (2002) - See page 45 & 292

porespy.metrics.linear_density(im, bins=25, voxel_size=1, log=False)[source]

Determines the probability that a point lies within a certain distance of the opposite phase along a specified direction

This relates directly the radial density function defined by Torquato [1], but instead of reporting the probability of lying within a stated distance to the nearest solid in any direciton, it considers only linear distances along orthogonal directions.The benefit of this is that anisotropy can be detected in materials by performing the analysis in multiple orthogonal directions.

Parameters:
  • im (ND-array) – An image with each voxel containing the distance to the nearest solid along a linear path, as produced by distance_transform_lin.
  • bins (int or array_like) – The number of bins or a list of specific bins to use
  • voxel_size (scalar) – The side length of a voxel. This is used to scale the chord lengths into real units. Note this is applied after the binning, so bins, if supplied, should be in terms of voxels, not length units.
Returns:

result

Return type:

named_tuple

References

[1] Torquato, S. Random Heterogeneous Materials: Mircostructure and Macroscopic Properties. Springer, New York (2002)

porespy.metrics.mesh_surface_area(mesh=None, verts=None, faces=None)[source]

Calculates the surface area of a meshed region

Parameters:
  • mesh (tuple) – The tuple returned from the mesh_region function
  • verts (array) – An N-by-ND array containing the coordinates of each mesh vertex
  • faces (array) – An N-by-ND array indicating which elements in verts form a mesh element.
Returns:

surface_area – The surface area of the mesh, calculated by skimage.measure.mesh_surface_area

Return type:

float

Notes

This function simply calls scikit-image.measure.mesh_surface_area, but it allows for the passing of the mesh tuple returned by the mesh_region function, entirely for convenience.

porespy.metrics.phase_fraction(im, normed=True)[source]

Calculates the number (or fraction) of each phase in an image

Parameters:
  • im (ND-array) – An ND-array containing integer values
  • normed (Boolean) – If True (default) the returned values are normalized by the total number of voxels in image, otherwise the voxel count of each phase is returned.
Returns:

result – A array of length max(im) with each element containing the number of voxels found with the corresponding label.

Return type:

1D-array

See also

porosity()

porespy.metrics.pore_size_distribution(im, bins=10, log=True, voxel_size=1)[source]

Calculate a pore-size distribution based on the image produced by the porosimetry or local_thickness functions.

Parameters:
  • im (ND-array) – The array of containing the sizes of the largest sphere that overlaps each voxel. Obtained from either porosimetry or local_thickness.
  • bins (scalar or array_like) – Either an array of bin sizes to use, or the number of bins that should be automatically generated that span the data range.
  • log (boolean) – If True (default) the size data is converted to log (base-10) values before processing. This can help
  • voxel_size (scalar) – The size of a voxel side in preferred units. The default is 1, so the user can apply the scaling to the returned results after the fact.
Returns:

result – A named-tuple containing several values:

R or logR - radius, equivalent to bin_centers

pdf - probability density function

cdf - cumulative density function

satn - phase saturation in differential form. For the cumulative saturation, just use cfd which is already normalized to 1.

bin_centers - the center point of each bin

bin_edges - locations of bin divisions, including 1 more value than the number of bins

bin_widths - useful for passing to the width argument of matplotlib.pyplot.bar

Return type:

named_tuple

Notes

(1) To ensure the returned values represent actual sizes be sure to scale the distance transform by the voxel size first (dt *= voxel_size)

plt.bar(psd.R, psd.satn, width=psd.bin_widths, edgecolor=’k’)

porespy.metrics.porosity(im)[source]

Calculates the porosity of an image assuming 1’s are void space and 0’s are solid phase.

All other values are ignored, so this can also return the relative fraction of a phase of interest.

Parameters:im (ND-array) – Image of the void space with 1’s indicating void space (or True) and 0’s indicating the solid phase (or False).
Returns:porosity – Calculated as the sum of all 1’s divided by the sum of all 1’s and 0’s.
Return type:float

See also

phase_fraction()

Notes

This function assumes void is represented by 1 and solid by 0, and all other values are ignored. This is useful, for example, for images of cylindrical cores, where all voxels outside the core are labelled with 2.

Alternatively, images can be processed with find_disconnected_voxels to get an image of only blind pores. This can then be added to the orignal image such that blind pores have a value of 2, thus allowing the calculation of accessible porosity, rather than overall porosity.

porespy.metrics.porosity_profile(im, axis)[source]

Returns a porosity profile along the specified axis

Parameters:
  • im (ND-array) – The volumetric image for which to calculate the porosity profile
  • axis (int) – The axis (0, 1, or 2) along which to calculate the profile. For instance, if axis is 0, then the porosity in each YZ plane is calculated and returned as 1D array with 1 value for each X position.
Returns:

result – A 1D-array of porosity along the specified axis

Return type:

1D-array

porespy.metrics.props_to_image(regionprops, shape, prop)[source]

Creates an image with each region colored according the specified prop, as obtained by regionprops_3d.

Parameters:
  • regionprops (list) – This is a list of properties for each region that is computed by PoreSpy’s regionprops_3D or Skimage’s regionsprops.
  • shape (array_like) – The shape of the original image for which regionprops was obtained.
  • prop (string) – The region property of interest. Can be a scalar item such as ‘volume’ in which case the the regions will be colored by their respective volumes, or can be an image-type property such as ‘border’ or ‘convex_image’, which will return an image composed of the sub-images.
Returns:

image – An ND-image the same size as the original image, with each region represented by the values specified in prop.

Return type:

ND-array

See also

props_to_DataFrame(), regionprops_3d()

porespy.metrics.props_to_DataFrame(regionprops)[source]

Returns a Pandas DataFrame containing all the scalar metrics for each region, such as volume, sphericity, and so on, calculated by regionprops_3D.

Parameters:regionprops (list) – This is a list of properties for each region that is computed by regionprops_3D. Because regionprops_3D returns data in the same list format as the regionprops function in Skimage you can pass in either.
Returns:DataFrame – A Pandas DataFrame with each region corresponding to a row and each column corresponding to a key metric. All the values for a given property (e.g. ‘sphericity’) can be obtained as val = df['sphericity']. Conversely, all the key metrics for a given region can be found with df.iloc[1].
Return type:Pandas DataFrame

See also

props_to_image(), regionprops_3d()

porespy.metrics.radial_density(im, bins=10, voxel_size=1)[source]

Computes radial density function by analyzing the histogram of voxel values in the distance transform. This function is defined by Torquato [1] as:

\[\int_0^\infty P(r)dr = 1.0\]

where P(r)dr is the probability of finding a voxel at a lying at a radial distance between r and dr from the solid interface. This is equivalent to a probability density function (pdf)

The cumulative distribution is defined as:

\[F(r) = \int_r^\infty P(r)dr\]

which gives the fraction of pore-space with a radius larger than r. This is equivalent to the cumulative distribution function (cdf).

Parameters:
  • im (ND-array) – Either a binary image of the pore space with True indicating the pore phase (or phase of interest), or a pre-calculated distance transform which can save time.
  • bins (int or array_like) – This number of bins (if int) or the location of the bins (if array). This argument is passed directly to Scipy’s histogram function so see that docstring for more information. The default is 10 bins, which reduces produces a relatively smooth distribution.
  • voxel_size (scalar) – The size of a voxel side in preferred units. The default is 1, so the user can apply the scaling to the returned results after the fact.
Returns:

result – A named-tuple containing several 1D arrays:

R - radius, equivalent to bin_centers

pdf - probability density function

cdf - cumulative density function

bin_centers - the center point of each bin

bin_edges - locations of bin divisions, including 1 more value than the number of bins

bin_widths - useful for passing to the width argument of matplotlib.pyplot.bar

Return type:

named_tuple

Notes

This function should not be taken as a pore size distribution in the explict sense, but rather an indicator of the sizes in the image. The distance transform contains a very skewed number of voxels with small values near the solid walls. Nonetheless, it does provide a useful indicator and it’s mathematical formalism is handy.

Torquato refers to this as the pore-size density function, and mentions that it is also known as the pore-size distribution function. These terms are avoided here since they have specific connotations in porous media analysis.

References

[1] Torquato, S. Random Heterogeneous Materials: Mircostructure and Macroscopic Properties. Springer, New York (2002) - See page 48 & 292

porespy.metrics.region_interface_areas(regions, areas, voxel_size=1, strel=None)[source]

Calculates the interfacial area between all pairs of adjecent regions

Parameters:
  • regions (ND-array) – An image of the pore space partitioned into individual pore regions. Note that zeros in the image will not be considered for area calculation.
  • areas (array_like) – A list containing the areas of each regions, as determined by region_surface_area. Note that the region number and list index are offset by 1, such that the area for region 1 is stored in areas[0].
  • voxel_size (scalar) – The resolution of the image, expressed as the length of one side of a voxel, so the volume of a voxel would be voxel_size-cubed. The default is 1.
  • strel (array_like) – The structuring element used to blur the region. If not provided, then a spherical element (or disk) with radius 1 is used. See the docstring for mesh_region for more details, as this argument is passed to there.
Returns:

result – A named-tuple containing 2 arrays. conns holds the connectivity information and area holds the result for each pair. conns is a N-regions by 2 array with each row containing the region number of an adjacent pair of regions. For instance, if conns[0, 0] is 0 and conns[0, 1] is 5, then row 0 of area contains the interfacial area shared by regions 0 and 5.

Return type:

named_tuple

porespy.metrics.region_surface_areas(regions, voxel_size=1, strel=None)[source]

Extracts the surface area of each region in a labeled image.

Optionally, it can also find the the interfacial area between all adjoining regions.

Parameters:
  • regions (ND-array) – An image of the pore space partitioned into individual pore regions. Note that zeros in the image will not be considered for area calculation.
  • voxel_size (scalar) – The resolution of the image, expressed as the length of one side of a voxel, so the volume of a voxel would be voxel_size-cubed. The default is 1.
  • strel (array_like) – The structuring element used to blur the region. If not provided, then a spherical element (or disk) with radius 1 is used. See the docstring for mesh_region for more details, as this argument is passed to there.
Returns:

result – A list containing the surface area of each region, offset by 1, such that the surface area of region 1 is stored in element 0 of the list.

Return type:

list

porespy.metrics.regionprops_3D(im)[source]

Calculates various metrics for each labeled region in a 3D image.

The regionsprops method in skimage is very thorough for 2D images, but is a bit limited when it comes to 3D images, so this function aims to fill this gap.

Parameters:im (array_like) – An imaging containing at least one labeled region. If a boolean image is received than the True voxels are treated as a single region labeled 1. Regions labeled 0 are ignored in all cases.
Returns:props – An augmented version of the list returned by skimage’s regionprops. Information, such as volume, can be found for region A using the following syntax: result[A-1].volume.

The returned list contains all the metrics normally returned by skimage.measure.regionprops plus the following:

’slice’: Slice indices into the image that can be used to extract the region

’volume’: Volume of the region in number of voxels.

’bbox_volume’: Volume of the bounding box that contains the region.

’border’: The edges of the region, found as the locations where the distance transform is 1.

’inscribed_sphere’: An image containing the largest sphere can can fit entirely inside the region.

’surface_mesh_vertices’: Obtained by applying the marching cubes algorithm on the region, AFTER first blurring the voxel image. This allows marching cubes more freedom to fit the surface contours. See also surface_mesh_simplices

’surface_mesh_simplices’: This accompanies surface_mesh_vertices and together they can be used to define the region as a mesh.

’surface_area’: Calculated using the mesh obtained as described above, using the porespy.metrics.mesh_surface_area method.

’sphericity’: Defined as the ratio of the area of a sphere with the same volume as the region to the actual surface area of the region.

’skeleton’: The medial axis of the region obtained using the skeletonize_3D method from skimage.

’convex_volume’: Same as convex_area, but translated to a more meaningful name.

Return type:list

See also

snow_partitioning()

Notes

This function may seem slow compared to the skimage version, but that is because they defer calculation of certain properties until they are accessed, while this one evalulates everything (inlcuding the deferred properties from skimage’s regionprops)

Regions can be identified using a watershed algorithm, which can be a bit tricky to obtain desired results. PoreSpy includes the SNOW algorithm, which may be helpful.

porespy.metrics.representative_elementary_volume(im, npoints=1000)[source]

Calculates the porosity of the image as a function subdomain size. This function extracts a specified number of subdomains of random size, then finds their porosity.

Parameters:
  • im (ND-array) – The image of the porous material
  • npoints (int) – The number of randomly located and sized boxes to sample. The default is 1000.
Returns:

result – A tuple containing the volume and porosity of each subdomain tested in arrays npoints long. They can be accessed as attributes of the tuple. They can be conveniently plotted by passing the tuple to matplotlib’s plot function using the * notation: plt.plot(*result, 'b.'). The resulting plot is similar to the sketch given by Bachmat and Bear [1]

Return type:

named_tuple

Notes

This function is frustratingly slow. Profiling indicates that all the time is spent on scipy’s sum function which is needed to sum the number of void voxels (1’s) in each subdomain.

Also, this function is a prime target for parallelization since the npoints are calculated independenlty.

References

[1] Bachmat and Bear. On the Concept and Size of a Representative Elementary Volume (Rev), Advances in Transport Phenomena in Porous Media (1987)

porespy.metrics.two_point_correlation_bf(im, spacing=10)[source]

Calculates the two-point correlation function using brute-force (see Notes)

Parameters:
  • im (ND-array) – The image of the void space on which the 2-point correlation is desired
  • spacing (int) – The space between points on the regular grid that is used to generate the correlation (see Notes)
Returns:

result – A tuple containing the x and y data for plotting the two-point correlation function, using the *args feature of matplotlib’s plot function. The x array is the distances between points and the y array is corresponding probabilities that points of a given distance both lie in the void space. The distance values are binned as follows: bins = range(start=0, stop=sp.amin(im.shape)/2, stride=spacing)

Return type:

named_tuple

Notes

The brute-force approach means overlaying a grid of equally spaced points onto the image, calculating the distance between each and every pair of points, then counting the instances where both pairs lie in the void space.

This approach uses a distance matrix so can consume memory very quickly for large 3D images and/or close spacing.

porespy.metrics.two_point_correlation_fft(im)[source]

Calculates the two-point correlation function using fourier transforms

Parameters:im (ND-array) – The image of the void space on which the 2-point correlation is desired
Returns:result – A tuple containing the x and y data for plotting the two-point correlation function, using the *args feature of matplotlib’s plot function. The x array is the distances between points and the y array is corresponding probabilities that points of a given distance both lie in the void space.
Return type:named_tuple

Notes

The fourier transform approach utilizes the fact that the autocorrelation function is the inverse FT of the power spectrum density. For background read the Scipy fftpack docs and for a good explanation see: http://www.ucl.ac.uk/~ucapikr/projects/KamilaSuankulova_BSc_Project.pdf