spam.label package#

Submodules#

spam.label.ITKwatershed module#

Library of wrapper functions for Simple ITK Copyright (C) 2020 SPAM Contributors

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

spam.label.ITKwatershed.watershed(binary, markers=None, watershedLevel=1)[source]#

This function runs an ITK watershed on a binary image and returns a labelled image. This function uses an interpixel watershed.

Parameters:
  • binary (3D numpy array) – This image which is non-zero in the areas which should be split by the watershed algorithm

  • markers (3D numpy array (optional, default = None)) – Not implemented yet, but try!

  • watershedLevel (int, optional) – Watershed level for merging maxima

Returns:

labelled – 3D array where each object is numbered

Return type:

3D numpy array of ints

spam.label.contacts module#

Library of SPAM functions for dealing with contacts between particles Copyright (C) 2020 SPAM Contributors

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

spam.label.contacts.contactingLabels(lab, labelsList=None, areas=False, boundingBoxes=None, centresOfMass=None)[source]#

This function returns contacting labels for a given label or list of labels, and optionally the number of voxels involved in the contact. This is designed for an interpixel watershed where there are no watershed lines, and so contact detection is done by dilating the particle in question once and seeing what other labels in comes in contact with.

Parameters:
  • lab (3D numpy array of ints ( of type spam.label.toolkit.labelType)) – An array of labels with 0 as the background

  • labels (int or a list of labels) – Labels for which we should compute neighbours

  • areas (bool, optional (default = False)) – Compute contact “areas”? i.e., number of voxels Careful, this is not a physically correct quantity, and is measured only as the overlap from label to its neightbours. See c    ontactPoint for something more meaningful

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

Returns:

  • For a single “labels”, a list of contacting labels.

  • if areas == True, then a list of “areas” is also returned.

  • For multiple labels, as above but a lsit of lists in the order of the labels

Note

Because of dilation, this function is not very safe at the edges, and will throw an exception

spam.label.contacts.contactPoints(lab, contactPairs, returnContactZones=False, boundingBoxes=None, centresOfMass=None, verbose=False)[source]#

Get the point, or area of contact between contacting particles. This is done by looking at the overlap of a 1-dilation of particle1 onto particle 2 and vice versa.

Parameters:
  • lab (3D numpy array of ints ( of type spam.label.toolkit.labelType)) – An array of labels with 0 as the background

  • contactPairs (list (or list of list) of labels) – Pairs of labels for which we should compute neighbours

  • returnContactZones (bool, optional (default = False)) – Output a labelled volume where contact zones are labelled according to the order of contact pairs? If False, the centres of mass of the contacts are returned

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

Returns:

  • if returnContactZones – lab : 3D numpy array of ints ( of type spam.label.toolkit.labelType)

  • else – Z Y X Centres of mass of contacts

spam.label.contacts.labelledContacts(lab, maximumCoordinationNumber=20)[source]#

Uniquely names contacts based on grains.

Parameters:
  • lab (3D numpy array of ints ( of type spam.label.toolkit.labelType)) – An array of labels with 0 as the background

  • maximumCoordinationNumber (int (optional, default = 20)) – Maximum number of neighbours to consider per grain

Returns:

contactVolumearray of ints

3D array where contact zones are uniquely labelled

Zarray of ints

Vector of length lab.max()+1 contaning the coordination number (number of touching labels for each label)

contactTablearray of ints

2D array of lab.max()+1 by 2*maximumCoordinationNumber containing, per grain: touching grain label, contact number

contactingLabelsarray of ints

2D array of numberOfContacts by 2 containing, per contact: touching grain label 1, touching grain label 2

Return type:

An array, containing

spam.label.contacts.fetchTwoGrains(volLab, labels, volGrey=None, boundingBoxes=None, padding=0, size_exclude=5)[source]#

Fetches the sub-volume of two grains from a labelled image

Parameters:
  • volLab (3D array of integers) – Labelled volume, with lab.max() labels

  • labels (1x2 array of integers) – the two labels that should be contained in the subvolume

  • volGrey (3D array) – Grey-scale volume

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • padding (integer) – padding of the subvolume for some purpose it might be benefitial to have a subvolume with a boundary of zeros

Returns:

  • subVolLab (3D array of integers) – labelled sub-volume containing the two input labels

  • subVolBin (3D array of integers) – binary subvolume

  • subVolGrey (3D array) – grey-scale subvolume

spam.label.contacts.localDetection(subVolGrey, localThreshold, radiusThresh=None)[source]#

local contact refinement checks whether two particles are in contact with a local threshold, that is higher than the global one used for binarisation

Parameters:
  • subVolGrey (3D array) – Grey-scale volume

  • localThreshold (integer or float, same type as the 3D array) – threshold for binarisation of the subvolume

  • radiusThresh (integer, optional) – radius for excluding patches that might not be relevant, e.g. noise can lead to patches that are not connected with the grains in contact the patches are calculated with equivalentRadii() Default is None and such patches are not excluded from the subvolume

Returns:

CONTACT – if True, that particles appear in contact for the local threshold

Return type:

boolean

Note

see https://doi.org/10.1088/1361-6501/aa8dbf for further information

spam.label.contacts.contactOrientations(volBin, volLab, watershed='ITK', peakDistance=1, markerShrink=True, verbose=False)[source]#

determines contact normal orientation between two particles uses the random walker implementation from skimage

Parameters:
  • volBin (3D array) – binary volume containing two particles in contact

  • volLab (3D array of integers) – labelled volume containing two particles in contact

  • watershed (string) – sets the basis for the determination of the orientation options are “ITK” for the labelled image from the input, “RW” for a further segmentation by the random walker Default = “ITK”

  • peakDistance (int, optional) – Distance in pixels used in skimage.feature.peak_local_max Default value is 1

  • markerShrink (boolean) – shrinks the markers to contain only one voxel. For large markers, the segmentation might yield strange results depending on the shape. Default = True

  • verbose (boolean (optional, default = False)) – True for printing the evolution of the process False for not printing the evolution of process

Returns:

  • contactNormal (1x3 array) – contact normal orientation in z,y,x coordinates the vector is flipped such that the z coordinate is positive

  • len(coordIntervox) (integer) – the number of points for the principal component analysis indicates the quality of the fit

  • notTreatedContact (boolean) – if False that contact is not determined if the contact consisted of too few voxels a fit is not possible or feasible

spam.label.contacts.localDetectionAssembly(volLab, volGrey, contactList, localThreshold, boundingBoxes=None, nProcesses=32, radiusThresh=4)[source]#

Local contact refinement of a set of contacts checks whether two particles are in contact with a local threshold using localDetection()

Parameters:
  • volLab (3D array) – Labelled volume

  • volGrey (3D array) – Grey-scale volume

  • contactList ((ContactNumber)x2 array of integers) – contact list with grain labels of each contact in a seperate row, as given by spam.label.contacts.labelledContacts in the list contactingLabels

  • localThreshold (integer or float, same type as the 3D array) – threshold for binarisation of the subvolume

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • nProcesses (integer (optional, default = nProcessesDefault)) – Number of processes for multiprocessing Default = number of CPUs in the system

  • radiusThresh (integer, optional) – radius in px for excluding patches that might not be relevant, e.g. noise can lead to patches that are not connected with the grains in contact the patches are calculated with equivalentRadii() Default is 4

Returns:

contactListRefined – refined contact list, based on the chosen local threshold

Return type:

(ContactNumber)x2 array of integers

Note

see https://doi.org/10.1088/1361-6501/aa8dbf for further information

spam.label.contacts.contactOrientationsAssembly(volLab, volGrey, contactList, watershed='ITK', peakDistance=5, boundingBoxes=None, nProcesses=32, verbose=False)[source]#

Determines contact normal orientation in an assembly of touching particles uses either directly the labelled image or the random walker implementation from skimage

Parameters:
  • volLab (3D array) – Labelled volume

  • volGrey (3D array) – Grey-scale volume

  • contactList ((ContactNumber)x2 array of integers) – contact list with grain labels of each contact in a seperate row, as given by spam.label.contacts.labelledContacts in the list contactingLabels or by spam.label.contacts.localDetectionAssembly()

  • watershed (string, optional) – sets the basis for the determination of the orientation options are “ITK” for the labelled image from the input, “RW” for a further segmentation by the random walker Default = “ITK”

  • peakDistance (int, optional) – Distance in pixels used in skimage.feature.peak_local_max Default value is 5

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • nProcesses (integer (optional, default = nProcessesDefault)) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean (optional, default = False)) – True for printing the evolution of the process False for not printing the evolution of process

Returns:

  • contactOrientations ((numContacts)x6 array) – contact normal orientation for every contact [labelGrainA, labelGrainB, orientationZ, orientationY, orientationX, intervoxel positions for PCA]

  • notTreatedContact (boolean) – if False that contact is not determined if the contact consisted of too few voxels a fit is not possible or feasible

spam.label.label module#

Library of SPAM functions for dealing with labelled images Copyright (C) 2020 SPAM Contributors

This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see <http://www.gnu.org/licenses/>.

spam.label.label.boundingBoxes(lab)[source]#

Returns bounding boxes for labelled objects using fast C-code which runs a single time through lab

Parameters:

lab (3D array of integers) – Labelled volume, with lab.max() labels

Returns:

boundingBoxes – This array contains, for each label, 6 integers:

  • Zmin, Zmax

  • Ymin, Ymax

  • Xmin, Xmax

Return type:

lab.max()x6 array of ints

Note

Bounding boxes are not slices and so to extract the correct bounding box from a numpy array you should use:

lab[ Zmin:Zmax+1, Ymin:Ymax+1, Xmin:Xmax+1 ]

Otherwise said, the bounding box of a single-voxel object at 1,1,1 will be:

1,1,1,1,1,1

Also note: for labelled images where some labels are missing, the bounding box returned for this case will be obviously wrong: e.g., Zmin = (z dimension-1) and Zmax = 0

spam.label.label.centresOfMass(lab, boundingBoxes=None, minVol=None)[source]#

Calculates (binary) centres of mass of each label in labelled image

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • minVol (int, optional) – The minimum volume in vx to be treated, any object below this threshold is returned as 0

Returns:

centresOfMass – This array contains, for each label, 3 floats, describing the centre of mass of each label in Z, Y, X order

Return type:

lab.max()x3 array of floats

spam.label.label.volumes(lab, boundingBoxes=None)[source]#

Calculates (binary) volumes each label in labelled image, using potentially slow numpy.where

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

Returns:

volumes – This array contains the volume in voxels of each label

Return type:

lab.max()x1 array of ints

spam.label.label.equivalentRadii(lab, boundingBoxes=None, volumes=None)[source]#

Calculates (binary) equivalent sphere radii of each label in labelled image

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • volumes (lab.max()x1 array of ints) – Vector contining volumes, if this is passed, the others are ignored

Returns:

equivRadii – This array contains the equivalent sphere radius in pixels of each label

Return type:

lab.max()x1 array of floats

spam.label.label.momentOfInertia(lab, boundingBoxes=None, minVol=None, centresOfMass=None)[source]#

Calculates (binary) moments of inertia of each label in labelled image

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

  • minVol (int, optional) – The minimum volume in vx to be treated, any object below this threshold is returned as 0 Default = default for spam.label.centresOfMass

Returns:

  • eigenValues (lab.max()x3 array of floats) – The values of the three eigenValues of the moment of inertia of each labelled shape

  • eigenVectors (lab.max()x9 array of floats) – 3 x Z,Y,X components of the three eigenValues in the order of the eigenValues

spam.label.label.ellipseAxes(lab, volumes=None, MOIeigenValues=None, enforceVolume=True, twoD=False)[source]#

Calculates length of half-axes a,b,c of the ellipitic fit of the particle. These are half-axes and so are comparable to the radius – and not the diameter – of the particle.

See appendix of for inital work:

“Three-dimensional study on the interconnection and shape of crystals in a graphic granite by X-ray CT and image analysis.”, Ikeda, S., Nakano, T., & Nakashima, Y. (2000).

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels Note: This is not strictly necessary if volumes and MOI is given

  • volumes (1D array of particle volumes (optional, default = None)) – Volumes of particles (length of array = lab.max())

  • MOIeigenValues (lab.max()x3 array of floats, (optional, default = None)) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • True) (enforceVolume = bool (default =) – Should a, b and c be scaled to enforce the fitted ellipse volume to be the same as the particle? This causes eigenValues are no longer completely consistent with fitted ellipse

  • twoD (bool (default = False)) – Are these in fact 2D ellipses? Not implemented!!

Returns:

ABCaxes – a, b, c lengths of particle in pixels

Return type:

lab.max()x3 array of floats

Note

Our elliptic fit is not necessarily of the same volume as the original particle, although by default we scale all axes linearly with enforceVolumes to enforce this condition.

Reminder: volume of an ellipse is (4/3)*pi*a*b*c

Useful check from TM: Ia = (4/15)*pi*a*b*c*(b**2+c**2)

Function contributed by Takashi Matsushima (University of Tsukuba)

spam.label.label.convertLabelToFloat(lab, vector)[source]#

Replaces all values of a labelled array with a given value. Useful for visualising properties attached to labels, e.g., sand grain displacements.

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • vector (a lab.max()x1 vector with values to replace each label with)

Returns:

relabelled

Return type:

3D array of converted floats

spam.label.label.makeLabelsSequential(lab)[source]#

This function fills gaps in labelled images, by relabelling them to be sequential integers. Don’t forget to recompute all your grain properties since the label numbers will change

Parameters:

lab (3D numpy array of ints ( of type spam.label.toolkit.labelType)) – An array of labels with 0 as the background

Returns:

lab – An array of labels with 0 as the background

Return type:

3D numpy array of ints ( of type spam.label.toolkit.labelType)

spam.label.label.getLabel(labelledVolume, label, boundingBoxes=None, centresOfMass=None, margin=0, extractCube=False, extractCubeSize=None, maskOtherLabels=True, labelDilate=0, labelDilateMaskOtherLabels=False)[source]#

Helper function to extract labels from a labelled image/volume. A dictionary is returned with the a subvolume around the particle. Passing boundingBoxes and centresOfMass is highly recommended.

Parameters:
  • labelVolume (3D array of ints) – 3D Labelled volume

  • label (int) – Label that we want information about

  • boundingBoxes (nLabels*2 array of ints, optional) – Bounding boxes as returned by boundingBoxes. Optional but highly recommended. If unset, bounding boxes are recalculated for every call.

  • centresOfMass (nLabels*3 array of floats, optional) – Centres of mass as returned by centresOfMass. Optional but highly recommended. If unset, centres of mass are recalculated for every call.

  • extractCube (bool, optional) – Return label subvolume in the middle of a cube? Default = False

  • extractCubeSize (int, optional) – half-size of cube to extract. Default = calculate minimum cube

  • margin (int, optional) – Extract a margin pixel margin around bounding box or cube. Default = 0

  • maskOtherLabels (bool, optional) – In the returned subvolume, should other labels be masked? If true, the mask is directly returned. Default = True

  • labelDilate (int, optional) – Number of times label should be dilated before returning it? This can be useful for catching the outside/edge of an image. margin should at least be equal to this value. Requires maskOtherLabels. Default = 0

  • labelDilateMaskOtherLabels (bool, optional) – Strictly cut the other labels out of the dilated image of the requested label? Only pertinent for positive labelDilate values. Default = False

Returns:

Keys:
subvol3D array of bools or ints

subvolume from labelled image

slicetuple of 3*slices

Slice used to extract subvol for the bounding box mode

sliceCubetuple of 3*slices

Slice used to extract subvol for the cube mode, warning, if the label is near the edge, this is the slice up to the edge, and so it will be smaller than the returned cube

boundingBox1D numpy array of 6 ints

Bounding box, including margin, in bounding box mode. Contains: [Zmin, Zmax, Ymin, Ymax, Xmin, Xmax] Note: uses the same convention as spam.label.boundingBoxes, so if you want to use this to extract your subvolume, add +1 to max

boundingBoxCube1D numpy array of 6 ints

Bounding box, including margin, in cube mode. Contains: [Zmin, Zmax, Ymin, Ymax, Xmin, Xmax] Note: uses the same convention as spam.label.boundingBoxes, so if you want to use this to extract your subvolume, add +1 to max

centreOfMassABS3*float

Centre of mass with respect to labelVolume

centreOfMassREL3*float

Centre of mass with respect to subvol

volumeInitialint

Volume of label (before dilating)

volumeDilatedint

Volume of label (after dilating, if requested)

Return type:

Dictionary containing

spam.label.label.getImagettesLabelled(lab1, label, Phi, im1, im2, searchRange, boundingBoxes, centresOfMass, margin=0, labelDilate=0, maskOtherLabels=True, applyF='all', volumeThreshold=100)[source]#

This function is responsible for extracting correlation windows (“imagettes”) from two larger images (im1 and im2) with the help of a labelled im1. This is generally to do image correlation, this function will be used for spam-ddic and pixelSearch modes.

Parameters:
  • lab1 (3D numpy array of ints) – Labelled image containing nLabels

  • label (int) – Label of interest

  • Phi (4x4 numpy array of floats) – Phi matrix representing the movement of imagette1, if not equal to I, imagette1 is deformed by the non-translation parts of Phi (F) and the displacement is added to the search range (see below)

  • im1 (3D numpy array) – This is the large input reference image of greyvalues

  • im2 (3D numpy array) – This is the large input deformed image of greyvalues

  • searchRange (6-component numpy array of ints) – This defines where imagette2 should be extracted with respect to imagette1’s position in im1. The 6 components correspond to [ Zbot Ztop Ybot Ytop Xbot Xtop ]. If Z, Y and X values are the same, then imagette2 will be displaced and the same size as imagette1. If ‘bot’ is lower than ‘top’, imagette2 will be larger in that dimension

  • boundingBoxes (nLabels*2 array of ints) – Bounding boxes as returned by boundingBoxes

  • centresOfMass (nLabels*3 array of floats) – Centres of mass as returned by centresOfMass

  • margin (int, optional) – Margin around the grain to extract in pixels Default = 0

  • labelDilate (int, optional) – How much to dilate the label before computing the mask? Default = 0

  • maskOtherLabels (bool, optional) – In the returned subvolume, should other labels be masked? If true, the mask is directly returned. Default = True

  • applyF (string, optional) – If a non-identity Phi is passed, should the F be applied to the returned imagette1? Options are: ‘all’, ‘rigid’, ‘no’ Default = ‘all’ Note: as of January 2021, it seems to make more sense to have this as ‘all’ for pixelSearch, and ‘no’ for local DIC

  • volumeThreshold (int, optional) – Pixel volume of labels that are discarded Default = 100

Returns:

‘imagette1’ : 3D numpy array,

’imagette1mask’: 3D numpy array of same size as imagette1 or None,

’imagette2’: 3D numpy array, bigger or equal size to imagette1

’returnStatus’: int,

Describes success in extracting imagette1 and imagette2. If == 1 success, otherwise negative means failure.

’pixelSearchOffset’: 3-component list of ints

Coordinates of the top of the pixelSearch range in im1, i.e., the displacement that needs to be added to the raw pixelSearch output to make it a im1 -> im2 displacement

Return type:

Dictionary

spam.label.label.labelsOnEdges(lab)[source]#

Return labels on edges of volume

Parameters:

lab (3D numpy array of ints) – Labelled volume

Returns:

uniqueLabels – List of labels on edges

Return type:

list of ints

spam.label.label.removeLabels(lab, listOfLabelsToRemove)[source]#

Resets a list of labels to zero in a labelled volume.

Parameters:
  • lab (3D numpy array of ints) – Labelled volume

  • listOfLabelsToRemove (list-like of ints) – Labels to remove

Returns:

lab – Labelled volume with desired labels blanked

Return type:

3D numpy array of ints

Note

You might want to use makeLabelsSequential after using this function, but don’t forget to recompute all your grain properties since the label numbers will change

spam.label.label.setVoronoi(lab, poreEDT=None, maxPoreRadius=10)[source]#

This function computes an approximate set Voronoi for a given labelled image. This is a voronoi which does not have straight edges, and which necessarily passes through each contact point, so it is respectful of non-spherical grains.

See: Schaller, F. M., Kapfer, S. C., Evans, M. E., Hoffmann, M. J., Aste, T., Saadatfar, M., … & Schroder-Turk, G. E. (2013). Set Voronoi diagrams of 3D assemblies of aspherical particles. Philosophical Magazine, 93(31-33), 3993-4017. https://doi.org/10.1080/14786435.2013.834389

and

Weis, S., Schonhofer, P. W., Schaller, F. M., Schroter, M., & Schroder-Turk, G. E. (2017). Pomelo, a tool for computing Generic Set Voronoi Diagrams of Aspherical Particles of Arbitrary Shape. In EPJ Web of Conferences (Vol. 140, p. 06007). EDP Sciences.

Parameters:
  • lab (3D numpy array of labelTypes) – Labelled image

  • poreEDT (3D numpy array of floats (optional, default = None)) – Euclidean distance map of the pores. If not given, it is computed by scipy.ndimage.distance_transform_edt

  • maxPoreRadius (int (optional, default = 10)) – Maximum pore radius to be considered (this threshold is for speed optimisation)

Returns:

lab – Image labelled with set voronoi labels

Return type:

3D numpy array of labelTypes

spam.label.label.labelTetrahedra(dims, points, connectivity, nThreads=1)[source]#

Labels voxels corresponding to tetrahedra according to a connectivity matrix and node points

Parameters:
  • dims (tuple representing z,y,x dimensions of the desired labelled output)

  • points (number of points x 3 array of floats) – List of points that define the vertices of the tetrahedra in Z,Y,X format. These points are referred to by line number in the connectivity array

  • connectivity (number of tetrahedra x 4 array of integers) – Connectivity matrix between points that define tetrahedra. Each line defines a tetrahedron whose number is the line number + 1. Each line contains 4 integers that indicate the 4 points in the nodePos array.

  • nThreads (int (optional, default=1)) – The number of threads used for the cpp parallelisation.

Returns:

Labelled 3D volume where voxels are numbered according to the tetrahedron number they fall inside of # WARNING: Voxels outside of the mesh get a value of #tetrahedra + 1

Return type:

3D array of ints, shape = dims

spam.label.label.labelTetrahedraForScipyDelaunay(dims, delaunay)[source]#

Labels voxels corresponding to tetrahedra coming from scipy.spatial.Delaunay Apparently the cells are not well-numbered, which causes a number of zeros when using labelledTetrahedra

Parameters:
  • dims (tuple) – represents z,y,x dimensions of the desired labelled output

  • delaunay ("delaunay" object) – Object returned by scipy.spatial.Delaunay( centres ) Hint: If using label.toolkit.centresOfMass( ), do centres[1:] to remove the position of zero.

Returns:

lab – Labelled 3D volume where voxels are numbered according to the tetrahedron number they fall inside of

Return type:

3D array of ints, shape = dims

spam.label.label.filterIsolatedCells(array, struct, size)[source]#

Return array with completely isolated single cells removed

Parameters:
  • array (3-D (labelled or binary) array) – Array with completely isolated single cells

  • struct (3-D binary array) – Structure array for generating unique regions

  • size (integer) – Size of the isolated cells to exclude (Number of Voxels)

Returns:

filteredArray – Array with minimum region size > size

Return type:

3-D (labelled or binary) array

Notes

function from: http://stackoverflow.com/questions/28274091/removing-completely-isolated-cells-from-python-array

spam.label.label.trueSphericity(lab, boundingBoxes=None, centresOfMass=None, gaussianFilterSigma=0.75, minVol=256)[source]#

Calculates the degree of True Sphericity (psi) for all labels, as per: “Sphericity measures of sand grains” Rorato et al., Engineering Geology, 2019 and originlly proposed in: “Volume, shape, and roundness of rock particles”, Waddell, The Journal of Geology, 1932.

True Sphericity (psi) = Surface area of equivalent sphere / Actual surface area

The actual surface area is computed by extracting each particle with getLabel, a Gaussian smooth of 0.75 is applied and the marching cubes algorithm from skimage is used to mesh the surface and compute the surface area.

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

  • gaussianFilterSigma (float, optional) – Sigma of the Gaussian filter used to smooth the binarised shape Default = 0.75

  • minVol (int, optional) – The minimum volume in vx to be treated, any object below this threshold is returned as 0 Default = 256 voxels

Returns:

trueSphericity – The values of the degree of true sphericity for each particle

Return type:

lab.max() array of floats

Notes

Function contributed by Riccardo Rorato (UPC Barcelona)

Due to numerical errors, this value can be >1, it should be clipped at 1.0

spam.label.label.convexVolume(lab, boundingBoxes=None, centresOfMass=None, volumes=None, nProcesses=32, verbose=True)[source]#

This function compute the convex hull of each label of the labelled image and return a list with the convex volume of each particle.

Parameters:
  • lab (3D array of integers) – Labelled volume, with lab.max() labels

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

  • volumes (lab.max()x1 array of ints) – Volumes in format returned by volumes If not defined (Default = None), it is recomputed by running volumes

  • nProcesses (integer (optional, default = nProcessesDefault)) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean (optional, default = False)) – True for printing the evolution of the process False for not printing the evolution of process

Returns:

convexVolume

Return type:

lab.max()x1 array of floats with the convex volume.

Note

convexVolume can only be computed for particles with volume greater than 3 voxels. If it is not the case, it will return 0.

spam.label.label.moveLabels(lab, PhiField, returnStatus=None, boundingBoxes=None, centresOfMass=None, margin=3, PhiCOM=True, threshold=0.5, labelDilate=0, nProcesses=32)[source]#

This function applies a discrete Phi field (from DDIC?) over a labelled image.

Parameters:
  • lab (3D numpy array) – Labelled image

  • PhiField ((multidimensional x 4 x 4 numpy array of floats)) – Spatial field of Phis

  • returnStatus (lab.max()x1 array of ints, optional) – Array with the return status for each label (usually returned by spam-ddic) If not defined (Default = None), all the labels will be moved If returnStatus[i] == 2, the label will be moved, otherwise is omitted and erased from the final image

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

  • margin (int, optional) – Margin, in pixels, to take in each label. Default = 3

  • PhiCOM (bool, optional) – Apply Phi to centre of mass of particle?, otherwise it will be applied in the middle of the particle’s bounding box. Default = True

  • threshold (float, optional) – Threshold to keep interpolated voxels in the binary image. Default = 0.5

  • labelDilate (int, optional) – Number of times label should be dilated/eroded before returning it. If labelDilate > 0 a dilated label is returned, while labelDilate < 0 returns an eroded label. Default = 0

  • nProcesses (integer (optional, default = nProcessesDefault)) – Number of processes for multiprocessing Default = number of CPUs in the system

Returns:

labOut – New labelled image with the labels moved by the deformations established by the PhiField.

Return type:

3D numpy array

spam.label.label.erodeLabels(lab, erosion=1, boundingBoxes=None, centresOfMass=None, nProcesses=32)[source]#

This function erodes a labelled image.

Parameters:
  • lab (3D numpy array) – Labelled image

  • erosion (int, optional) – Erosion level

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

  • nProcesses (integer (optional, default = nProcessesDefault)) – Number of processes for multiprocessing Default = number of CPUs in the system

Returns:

erodeImage – New labelled image with the eroded labels.

Return type:

3D numpy array

Note

The function makes use of spam.label.moveLabels() to generate the eroded image.

spam.label.label.convexFillHoles(lab, boundingBoxes=None, centresOfMass=None)[source]#

This function fills the holes computing the convex volume around each label.

Parameters:
  • lab (3D numpy array) – Labelled image

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

Returns:

labOut – New labelled image.

Return type:

3D numpy array

Note

The function works nicely for convex particles. For non-convex particles, it will alter the shape.

spam.label.label.getNeighbours(lab, listOfLabels, method='getLabel', neighboursRange=None, centresOfMass=None, boundingBoxes=None)[source]#

This function computes the neighbours for a list of labels.

Parameters:
  • lab (3D numpy array) – Labelled image

  • listOfLabels (list of ints) – List of labels to which the neighbours will be computed

  • method (string) – Method to compute the neighbours. ‘getLabel’ : The neighbours are the labels inside the subset obtained through spam.getLabel() ‘mesh’ : The neighbours are computed using a tetrahedral connectivity matrix Default = ‘getLabel’

  • neighboursRange (int) – Parameter controlling the search range to detect neighbours for each method. For ‘getLabel’, it correspond to the size of the subset. Default = meanRadii For ‘mesh’, it correspond to the size of the alpha shape used for carving the mesh. Default = 5*meanDiameter.

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running centresOfMass

Returns:

neighbours – List with the neighbours for each label in listOfLabels.

Return type:

list

spam.label.label.detectUnderSegmentation(lab, nProcesses=32, verbose=True)[source]#

This function computes the coefficient of undersegmentation for each particle, defined as the ratio of the convex volume and the actual volume.

Parameters:
  • lab (3D numpy array) – Labelled image

  • nProcesses (integer (optional, default = nProcessesDefault)) – Number of processes for multiprocessing Default = number of CPUs in the system

  • verbose (boolean (optional, default = False)) – True for printing the evolution of the process

Returns:

underSegCoeff – An array of float values that suggests the respective labels are undersegmentated.

Return type:

lab.max() array of floats

Note

For perfect convex particles, any coefficient higher than 1 should be interpreted as a particle with undersegmentation problems. However, for natural materials the threshold to define undersegmentation varies. It is suggested to plot the histogram of the undersegmentation coefficient and select the threshold accordingly.

spam.label.label.detectOverSegmentation(lab)[source]#

This function computes the coefficient of oversegmentation for each particle, defined as the ratio between a characteristic lenght of the maximum contact area and a characteristic length of the particle.

Parameters:

lab (3D numpy array) – Labelled image

Returns:

  • overSegCoeff (lab.max() array of floats) – An array of float values with the oversegmentation coefficient

  • sharedLabel (lab.max() array of floats) – Array of floats with the the oversegmentation coefficient neighbours - label that share the maximum contact area

Note

The threshold to define oversegmentation is dependent on each material and conditions of the test. It is suggested to plot the histogram of the oversegmentation coefficient and select the threshold accordingly.

spam.label.label.fixUndersegmentation(lab, imGrey, targetLabels, underSegCoeff, boundingBoxes=None, centresOfMass=None, imShowProgress=False, verbose=True)[source]#

This function fixes undersegmentation problems, by performing a watershed with a higher local threshold for the problematic labels.

Parameters:
  • lab (3D numpy array) – Labelled image

  • imGrey (3D numpy array) – Normalised greyscale of the labelled image, with a greyscale range between 0 and 1 and with void/solid peaks at 0.25 and 0.75, respectively. You can use helpers.histogramTools.findHistogramPeaks and helpers.histogramTools.histogramNorm to obtain a normalized greyscale image.

  • targetLabels (int or a list of labels) – List of target labels to solve undersegmentation

  • underSegCoeff (lab.max() array of floats) – Undersegmentation coefficient as returned by detectUnderSegmentation

  • boundingBoxes (lab.max()x6 array of ints, optional) – Bounding boxes in format returned by boundingBoxes. If not defined (Default = None), it is recomputed by running boundingBoxes

  • centresOfMass (lab.max()x3 array of floats, optional) – Centres of mass in format returned by centresOfMass. If not defined (Default = None), it is recomputed by running ``centresOfMass``boundingBoxes : 3D numpy array Labelled image

  • imShowProgress (bool, optional) – Graphical interface to observe the process for each label. Default = False

  • verbose (boolean (optional, default = False)) – True for printing the evolution of the process

Returns:

lab – Labelled image after running makeLabelsSequential

Return type:

3D numpy array

spam.label.label.fixOversegmentation(lab, targetLabels, sharedLabel, verbose=True, imShowProgress=False)[source]#

This function fixes oversegmentation problems, by merging each target label with its oversegmentation coefficient neighbour.

Parameters:
  • lab (3D numpy array) – Labelled image

  • targetLabels (int or a list of labels) – List of target labels to solve oversegmentation

  • sharedLabel (lab.max() array of floats) – List ofoversegmentation coefficient neighbour as returned by detectOverSegmentation

  • imShowProgress (bool, optional) – Graphical interface to observe the process for each label. Default = False

  • verbose (boolean (optional, default = False)) – True for printing the evolution of the process

Returns:

lab – Labelled image after running makeLabelsSequential

Return type:

3D numpy array

spam.label.label.shuffleLabels(lab)[source]#

This function re-assigns randomly the labels of a label image, usually for visualisation purposes.

Parameters:

lab (3D array of integers) – Labelled volume, with lab.max() labels

Returns:

lab – Labelled volume, with lab.max() labels

Return type:

3D array of integers

Module contents#