Segmentation of lnPi (segment)#

Routines to segment lnPi

  1. find max/peaks in lnPi

  2. segment lnPi about these peaks

  3. determine free energy difference between segments a. Merge based on low free energy difference

  4. combination of 1-3.

Functions:

peak_local_max_adaptive(data, *[, mask, ...])

Find local max with fall backs min_distance and filter.

Classes:

Segmenter([peak_kws, watershed_kws])

Data segmenter:

PhaseCreator(nmax[, nmax_peak, ref, ...])

Helper class to create phases

BuildPhasesBase(x, phase_creator)

Base class to build Phases objects from scalar values of lnz.

BuildPhases_mu(lnz, phase_creator)

create phases from scalar value of mu for fixed value of mu for other species

BuildPhases_dmu(dlnz, phase_creator)

Create phases from scalar value of mu at fixed value of dmu for other species

lnpy.segment.peak_local_max_adaptive(data, *, mask=None, min_distance=None, style='indices', threshold_rel=0.0, threshold_abs=0.2, num_peaks_max=None, connectivity=None, errors='warn', **kwargs)[source]#

Find local max with fall backs min_distance and filter.

This is an adaptation of peak_local_max(), which is called iteratively until the number of peaks is less than num_peaks_max.

Parameters:
  • data (array-like) – Image data to analyze

  • mask (None or ndarray of bool) – Mask using “image” convention. Where mask is True, values are included.

  • min_distance (int or sequence of int, optional) – min_distance parameter. If sequence, then call peak_local_max() until number of peaks <=num_peaks_max. Default value is (5, 10, 15, 20, 25).

  • style ({'indices', 'mask', 'marker'}) – Controls output style

    • indices : indices of peaks

    • mask : array of bool marking peak locations

    • marker : array of int

  • threshold_rel, threshold_abs (float) – thresholds parameters

  • num_peaks_max (int, optional) – Max number of maxima/peaks to find. If not specified, any number of peaks allowed.

  • connectivity (int, optional) – Maximum number of orthogonal hops to consider a pixel/voxel as a neighbor. Accepted values are ranging from 1 to input.ndim. If None, a full connectivity of input.ndim is used.

  • errors ({'ignore','raise','warn'}, default 'warn') –

    • If raise, raise exception if npeaks > num_peaks_max

    • If ignore, return all found maxima

    • If warn, raise warning if npeaks > num_peaks_max

  • **kwargs – Extra arguments to peak_local_max()

Returns:

out (ndarray of int or list of ndarray of bool) – Depending on the value of indices.

Notes

The option mask is passed as the value labels in peak_local_max()

See also

peak_local_max, label

class lnpy.segment.Segmenter(peak_kws=None, watershed_kws=None)[source]#

Bases: object

Data segmenter:

Parameters:

Methods:

peaks(data[, mask, num_peaks_max, style])

Interface to peak_local_max_adaptive() with default values from self.

watershed(data, markers, mask[, connectivity])

Interface to skimage.segmentation.watershed() function

segment_lnpi(lnpi[, markers, find_peaks, ...])

Perform segmentations of lnPi object using watershed on negative of lnPi data.

peaks(data, mask=None, *, num_peaks_max=None, style='marker', **kwargs)[source]#

Interface to peak_local_max_adaptive() with default values from self.

Parameters:
  • data (array-like) – Image data to analyze

  • num_peaks_max (int, optional) – Max number of maxima/peaks to find. If not specified, any number of peaks allowed.

  • style ({'indices', 'mask', 'marker'}) – Controls output style

    • indices : indices of peaks

    • mask : array of bool marking peak locations

    • marker : array of int

  • mask (None or ndarray of bool) – Mask using “image” convention. Where mask is True, values are included.

  • **kwargs – Extra arguments to peak_local_max_adaptive().

Returns:

ndarray of int or sequence of ndarray – If style=='marker', then return label array. Otherwise, return indices of peaks.

Notes

Any value not set will be inherited from self.peak_kws

watershed(data, markers, mask, connectivity=None, **kwargs)[source]#

Interface to skimage.segmentation.watershed() function

Parameters:
  • data (array-like) – Image data to analyze

  • markers (int, or ndarray of int, optional) – Same shape as image. The desired number of markers, or an array marking the basins with the values to be assigned in the label matrix. Zero means not a marker. If None (no markers given), the local minima of the image are used as markers.

  • mask (None or ndarray of bool) – Mask using “image” convention. Where mask is True, values are included.

  • connectivity (ndarray, optional) – An array with the same number of dimensions as image whose non-zero elements indicate neighbors for connection. Following the scipy convention, default is a one-connected array of the dimension of the image.

  • **kwargs – Extra arguments to watershed()

Returns:

labels (ndarray of int) – Each unique value i in labels indicates a mask. That is labels == i is a mask for feature i.

See also

watershed

segment_lnpi(lnpi, markers=None, find_peaks=True, num_peaks_max=None, connectivity=None, peaks_kws=None, watershed_kws=None)[source]#

Perform segmentations of lnPi object using watershed on negative of lnPi data.

Parameters:
  • lnpi (lnPiMasked) – Object to be segmented

  • markers (int, or ndarray of int, optional) – Same shape as image. The desired number of markers, or an array marking the basins with the values to be assigned in the label matrix. Zero means not a marker. If None (no markers given), the local minima of the image are used as markers.

  • find_peaks (bool, default True) – If True, use peak_local_max_adaptive() to construct markers.

  • num_peaks_max (int, optional) – Max number of maxima/peaks to find. If not specified, any number of peaks allowed.

  • connectivity (ndarray, optional) – An array with the same number of dimensions as image whose non-zero elements indicate neighbors for connection. Following the scipy convention, default is a one-connected array of the dimension of the image.

  • peak_kws (mapping, optional) – Optional parameters to peak_local_max_adaptive()

  • watershed_kws (mapping, optional) – Optional parameters to watershed()

Returns:

labels (ndarray of int) – Each unique value i in labels indicates a mask. That is labels == i is a mask for feature i.

class lnpy.segment.PhaseCreator(nmax, nmax_peak=None, ref=None, segmenter=None, segment_kws=None, tag_phases=None, phases_factory=None, free_energy_kws=None, merge_kws=None)[source]#

Bases: object

Helper class to create phases

Parameters:
  • nmax (int) – Maximum number of phases to allow

  • nmax_peak (int, optional) – if specified, the allowable number of peaks to locate. This can be useful for some cases. These phases will be merged out at the end.

  • ref (lnPiMasked, optional) – Reference object.

  • segmenter (Segmenter, optional) – segmenter object to create labels/masks. Defaults to using base segmenter.

  • segment_kws (mapping, optional) – Optional arguments to be passed to Segmenter.segment_lnpi().

  • tag_phases (callable(), optional) – Optional function which takes a list of lnPiMasked objects and returns on integer label for each object.

  • phases_factory (callable(), optional) – Factory function for returning Collection from a list of lnPiMasked object. Defaults to lnPiCollection.from_list().

  • free_energy_kws (mapping, optional) – Optional arguments to …

  • merge_kws (mapping, optional) – Optional arguments to merge_regions()

Methods:

build_phases([lnz, ref, efac, nmax, ...])

Construct 'phases' for a lnPi object.

build_phases_mu(lnz)

Factory constructor at fixed values of mu

build_phases_dmu(dlnz)

Factory constructor at fixed values of dmu.

build_phases(lnz=None, ref=None, *, efac=None, nmax=None, nmax_peak=None, connectivity=None, reweight_kws=None, merge_phase_ids=True, merge_phases=True, phases_factory=True, phase_kws=None, segment_kws=None, free_energy_kws=None, merge_kws=None, tag_phases=None)[source]#

Construct ‘phases’ for a lnPi object.

This is quite an involved process. The steps are

  • Optionally find the location of the maxima in lnPi.

  • Segment lnpi using watershed

  • Merge phases which are energetically similar

  • Optionally merge phases which have the same phase_id

Parameters:
  • lnz (int or sequence of int, optional) – lnz value to evaluate ref at. If not specified, use ref.lnz

  • ref (lnPiMasked) – Object to be segmented

  • efac (float, optional) – Optional value to use in energetic merging of phases.

  • nmax (int, optional) – Maximum number of phases. Defaults to self.nmax

  • nmax_peak (int, optional) – Maximum number of peaks to allow in peak_local_max_adaptive(). Note that this value can be larger than nmax. Defaults to self.nmax_peak.

  • connectivity (int, optional) – Maximum number of orthogonal hops to consider a pixel/voxel as a neighbor. Accepted values are ranging from 1 to input.ndim. If None, a full connectivity of input.ndim is used.

  • reweight_kws (mapping, optional) – Extra arguments to ref.reweight

  • merge_phase_ids (bool, default True) – If True and calling tag_phases routine, merge phases with same phase_id.

  • phases_factory (callable() or bool, default True) – Function to convert list of phases into Phases object. If phases_factory True, revert to self.phases_factory. If phases_factory is False, do not apply a factory, and return list of lnpy.lnpidata.lnPiMasked and array of phase indices.

  • phase_kws (mapping, optional) – Extra arguments to phases_factory

  • segment_kws (mapping, optional) – Extra arguments to self.segmenter.segment_lnpi

  • free_energy_kws (mapping, optional) – Extra arguments to free energy calculation

  • merge_kws (mapping, optional) – Extra arguments to merge

  • tag_phases (callable(), optional) – Function to tag phases. Defaults to self.tag_phases

Returns:

output (list of lnPiMasked and ndarray, or lnPiCollection) – If no phase creator, return list of lnPiMasked objects and array of phase indices. Otherwise, lnPiCollection object.

build_phases_mu(lnz)[source]#

Factory constructor at fixed values of mu

Parameters:

{lnz_buildphases_mu}

See also

BuildPhases_mu

Examples

>>> import lnpy.examples
>>> e = lnpy.examples.hsmix_example()

The default build phases from this multicomponent system requires specifies the activity for both species. For example:

>>> e.phase_creator.build_phases([0.1, 0.2])
<class lnPiCollection>
lnz_0  lnz_1  phase
0.1    0.2    0        [0.1, 0.2]
              1        [0.1, 0.2]
dtype: object

But if we want to creat phases at a fixed value of either lnz_0 or lnz_1, we can do the following:

>>> b = e.phase_creator.build_phases_mu([None, 0.5])

Note the syntax [None, 0.5]. This means that calling b(lnz_0) will create a new object at [lnz_0, 0.5].

>>> b(0.1)
<class lnPiCollection>
lnz_0  lnz_1  phase
0.1    0.5    0        [0.1, 0.5]
              1        [0.1, 0.5]
dtype: object

Likewise, we can fix lnz_0 with

>>> b = e.phase_creator.build_phases_mu([0.5, None])
>>> b(0.1)
<class lnPiCollection>
lnz_0  lnz_1  phase
0.5    0.1    0        [0.5, 0.1]
              1        [0.5, 0.1]
dtype: object

To create an object at fixed value of dmu_i = lnz_i - lnz_fixed, we use the following:

>>> b = e.phase_creator.build_phases_dmu([None, 0.5])

Now any phase created will have lnz = [lnz_0, 0.5 + lnz_0]

>>> b(0.5)
<class lnPiCollection>
lnz_0  lnz_1  phase
0.5    1.0    0        [0.5, 1.0]
              1        [0.5, 1.0]
dtype: object
build_phases_dmu(dlnz)[source]#

Factory constructor at fixed values of dmu.

Parameters:

{dlnz_buildphases_dmu}

class lnpy.segment.BuildPhasesBase(x, phase_creator)[source]#

Bases: object

Base class to build Phases objects from scalar values of lnz.

Attributes:

index

Index number which varies

Methods:

__call__(lnz_index, *[, phases_factory])

Build phases from scalar value of lnz.

property index#

Index number which varies

__call__(lnz_index, *, phases_factory=True, **kwargs)[source]#

Build phases from scalar value of lnz.

Parameters:
Returns:

output (list of lnPiMasked and ndarray, or lnPiCollection) – If no phase creator, return list of lnPiMasked objects and array of phase indices. Otherwise, lnPiCollection object.

class lnpy.segment.BuildPhases_mu(lnz, phase_creator)[source]#

Bases: BuildPhasesBase

create phases from scalar value of mu for fixed value of mu for other species

Parameters:
  • lnz (list of float or None) – list with one element equal to None. This is the component which will be varied For example, lnz=[lnz0, None, lnz2] implies use values of lnz0,lnz2 for components 0 and 2, and vary component 1.

  • phase_creator (PhaseCreator) – Factory method to create phases collection object. For example, lnPiCollection.from_list().

Methods:

__call__(lnz_index, *[, phases_factory])

Build phases from scalar value of lnz.

Attributes:

index

Index number which varies

__call__(lnz_index, *, phases_factory=True, **kwargs)[source]#

Build phases from scalar value of lnz.

Parameters:
Returns:

output (list of lnPiMasked and ndarray, or lnPiCollection) – If no phase creator, return list of lnPiMasked objects and array of phase indices. Otherwise, lnPiCollection object.

property index#

Index number which varies

class lnpy.segment.BuildPhases_dmu(dlnz, phase_creator)[source]#

Bases: BuildPhasesBase

Create phases from scalar value of mu at fixed value of dmu for other species

Parameters:
  • dlnz (list of float or None) – list with one element equal to None. This is the component which will be varied For example, dlnz=[dlnz0,None,dlnz2] implies use values of dlnz0,dlnz2 for components 0 and 2, and vary component 1. dlnz_i = lnz_i - lnz_index, where lnz_index is the value varied.

  • phase_creator (PhaseCreator) – Factory method to create phases collection object. For example, lnPiCollection.from_list().

Methods:

__call__(lnz_index, *[, phases_factory])

Build phases from scalar value of lnz.

Attributes:

index

Index number which varies

__call__(lnz_index, *, phases_factory=True, **kwargs)[source]#

Build phases from scalar value of lnz.

Parameters:
Returns:

output (list of lnPiMasked and ndarray, or lnPiCollection) – If no phase creator, return list of lnPiMasked objects and array of phase indices. Otherwise, lnPiCollection object.

property index#

Index number which varies