tramway.tessellation package

tramway.tessellation.base module

class tramway.tessellation.base.Partition(points=None, tessellation=None, cell_index=None, location_count=None, bounding_box=None, param=None, locations=None, translocations=None)

Bases: tramway.core.lazy.Lazy

Container datatype for molecule location datasets partitioned using a tessellation.

A Partition instance conveniently stores the tessellation (tessellation) and the proper partition of the data (cell_index) together with the data itself (points) and a few more intermediate results frequently derivated from a data partition.

locations and translocations are aliases of points. No control is performed on whether translocations are actual translocations for example.

The partition cell_index may be in any of the following formats:

array
Cell index of size the number of data points. The element at index i is the cell index of the i th point or -1 if the i th point is not assigned to any cell.
pair of arrays
Point-cell association in the shape of a sparse representation (point_index, cell_index) such that for all i the point_index[i] point is in the cell_index[i] cell.
sparse matrix (scipy.sparse)
number_of_points * number_of_cells matrix with nonzero element wherever the corresponding point is in the corresponding cell.

Note

If the point coordinates are defined as a DataFrame, point indices are row indices and NOT row labels (see also iloc).

See also Tessellation.cell_index().

new in 0.6: the __len__() method returns the number of assigned points, counted as many times as the number of cells they are assigned to.

points

the original (trans-)location coordinates, unchanged.

Type:array-like
tessellation

The tessellation that defined the partition.

Type:Tessellation
cell_index

Point-cell association (or data partition).

Type:numpy.ndarray or pair of arrays or sparse matrix
location_count

point count per cell; location_count[i] is the number of points in cell i.

Type:numpy.ndarray, lazy
bounding_box

2 * D array with lower values in first row and upper values in second row, where D is the dimension of the point data.

Type:array-like, lazy
param

Arguments involved in the tessellation and the partition steps, as key-value pairs. Such information is maintained in Partition so that it can be stored in .rwa files and retrieve for traceability.

Type:dict

Functional dependencies:

  • setting tessellation unsets cell_index
  • setting points unsets cell_index and bounding_box
  • setting cell_index unsets location_count
bounding_box
cell_index
descriptors(*vargs, **kwargs)

Proxy method for Tessellation.descriptors().

freeze()

Proxy method for Tessellation.freeze().

location_count
locations
number_of_cells

Proxy property for Tessellation.number_of_cells.

param
points
tessellation
translocations
tramway.tessellation.base.CellStats

alias of tramway.tessellation.base.Partition

tramway.tessellation.base.point_adjacency_matrix(cells, symetric=True, cell_labels=None, adjacency_labels=None)

Adjacency matrix of data points such that a given pair of points is defined as adjacent iif they belong to adjacent and distinct cells.

Parameters:
  • cells (Partition) – Partition with both partition and tessellation defined.
  • symetric (bool) – If False, the returned matrix will not be symetric, i.e. wherever i->j is defined, j->i is not.
  • cell_labels (callable) – Takes an array of cell labels as input (see Tessellation.cell_label) and returns a bool array of equal shape.
  • adjacency_labels (callable) – Takes an array of edge labels as input (see Tessellation.adjacency_label) and returns a bool array of equal shape.
Returns:

Sparse square matrix with as many rows as data points.

Return type:

scipy.sparse.csr_matrix

class tramway.tessellation.base.Tessellation(scaler=None)

Bases: tramway.core.lazy.Lazy

Abstract class for tessellations.

The methods to be implemented are tessellate() and cell_index().

scaler

scaler.

Type:tramway.core.scaler.Scaler
cell_adjacency

square adjacency matrix for cells. If _adjacency_label is defined, _cell_adjacency should be sparse and the explicit elements should be indices in _adjacency_label.

Type:sparse matrix
cell_label

cell labels with as many elements as cells.

Type:numpy.ndarray
adjacency_label

inter-cell edge labels with as many elements as there are edges.

Type:numpy.ndarray
adjacency_label

Inter-cell edge labels, numpy.ndarray with as many elements as edges.

cell_adjacency

Square cell adjacency matrix. If adjacency_label is defined, cell_adjacency is sparse and the explicit elements are indices in adjacency_label.

cell_index(points, format=None, select=None, **kwargs)

Partition.

The returned value depends on the format input argument:

  • array: returns a vector v such that v[i] is cell index for
    point index i or -1.
  • pair: returns a pair of I-sized arrays (p, c) where, for each
    point-cell association i in range(I), p[i] is a point index and c[i] is a corresponding cell index.
  • matrix or coo or csr or csc:
    returns a sparse matrix with points as rows and cells as columns; non-zeros are all True or float weights.

By default with format undefined, any implementation may favor any format.

Note that array may not be an acceptable format and cell_index() may not comply with format='index' unless select is defined. When a location or a translocation is associated to several cells, select chooses a single cell among them.

The default implementation calls format_cell_index() on the result of an abstract _cell_index method that any Tessellation implementation can implement instead of cell_index().

See also format_cell_index().

Parameters:
  • points (pandas.DataFrame) – point (location) coordinates.
  • format (str) – preferred representation of the point-cell association (or partition).
  • select (callable) – takes the point index, an array of cell indices and the tessellation as arguments, and returns a cell index or -1 for no cell.
cell_label

Cell labels, numpy.ndarray with as many elements as there are cells.

contour(cell, distance=1, fallback=False, adjacency=None, **kwargs)

Select a close path around a cell.

This method may be moved out of Tessellation in the near future.

descriptors(points, asarray=False)

Keep the data columns that were involved in growing the tessellation.

Parameters:
  • points (pandas.DataFrame) – point coordinates.
  • asarray (bool) – returns a numpy.ndarray.
Returns:

selected coordinates; the data may not be copied.

Return type:

array-like

freeze()

Delete the data required to and only to incrementally update the tessellation.

This may save large amounts of memory but differs from unload() in that the subsequent data loss may not be undone.

neighbours(i)

Indices of neighbour cells.

Reminder: neighbours with False or negative adjacency labels
are also returned. Consider simplified_adjacency instead to ignore these neighbours.
Parameters:i (int) – cell index.
Returns:indices of the neighbour cells of cell i.
Return type:numpy.ndarray
number_of_cells

Number of cells.

Read-only property.

scaler
set_adjacency(i, j, label=True)

Set (or unset) the adjacency between two cells.

If adjacency_label is not defined, label should be True or False, otherwise adjacency_label is created with default value 1 and cell_adjacency is modified to encode indices in adjacency_label instead of labels.

Parameters:
  • i (int) – cell index.
  • j (int) – cell index.
  • label (scalar) – adjacency label.
simplified_adjacency(adjacency=None, label=None, format='coo')

Simplified copy of cell_adjacency as a scipy.sparse.spmatrix sparse matrix with no explicit zeros.

Non-zero values indicate adjacency and all these values are strictly positive.

In addition, cells with negative (or null) labels are also disconnected from their neighbours.

Labels are cell_label by default. Alternative labels can be provided as label.

To prevent label-based disconnection, set label to False.

Multiple arrays of labels can also be supplied as a tuple. Note that explicit labels always supersede cell_label and the later should be explicitely listed in the tuple so that it is applied in combination with other label arrays.

Parameters:
  • adjacency (scipy.sparse.spmatrix) – adjacency matrix (cell_adjacency is used if adjacency is None).
  • label (bool or array-like) – cell labels.
  • format (str) – any of ‘coo’, ‘csr’ and ‘csc’.
Returns:

simplified adjacency matrix.

Return type:

scipy.sparse.spmatrix

tessellate(points, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

class tramway.tessellation.base.Delaunay(scaler=None)

Bases: tramway.tessellation.base.Tessellation

Delaunay graph.

A cell is represented by a centroid and an edge of the graph represents a neighour relationship between two cells.

Delaunay implements the nearest neighour feature and support for cell overlap.

cell_centers

coordinates of the cell centers.

Type:numpy.ndarray
cell_centers

Unscaled coordinates of the cell centers (numpy.ndarray).

cell_index(points, format=None, select=None, knn=None, radius=None, min_location_count=None, metric='euclidean', filter=None, filter_descriptors_only=False, fast_min_nn=True, **kwargs)

See Tessellation.cell_index().

A single array representation of the point-cell association may not be possible with knn defined, because a point can be associated to multiple cells. If such a case happens the default output format will be pair.

In addition to the values allowed by Tessellation.cell_index(), format admits value force array that acts like format='array', select=nearest_cell(...). The implementation however is more straight-forward and simply ignores the minimum number of nearest neighbours if provided.

new in 0.6: fast_min_nn fixes a bug and should be set to False; default is True for backward-compatibility.

Parameters:
  • points – see Tessellation.cell_index().
  • format – see Tessellation.cell_index(); additionally admits force array.
  • select – see Tessellation.cell_index().
  • knn (int or tuple or callable) – If int: minimum number of points per cell (or of nearest neighours to the cell center). Cells may overlap and the returned cell index may be a sparse point-cell association. If tuple (pair of ints): minimum and maximum number of points per cell respectively. If callable: takes a cell index and returns the minimum and maximum number of points.
  • radius (float or tuple or callable) – If float: distance from the cell center; smaller cells may include locations from neighbour cells and larger cells may include only part of their associated locations. If tuple (pair of floats): minimum and maximum radius of a cell respectively; any of these values can be None. If callable: takes a cell index and returns the minimum and maximum radius.
  • min_location_count (int) – minimum number of points for a cell to be included in the labeling. This argument applies before knn. The points in these cells, if not associated with another cell, are labeled -1. The other cell labels do not change.
  • metric (str) – any metric name understandable by cdist().
  • filter (callable) – takes the calling instance, a cell index and the corresponding subset of points; returns True if the corresponding cell should be included in the labeling.
  • filter_descriptors_only (bool) – whether filter should get points as descriptors only.
  • fast_min_nn (bool) – new in 0.6 if True, the initial assignment for small cells is overriden by the extension ball; if False, the initial assignment is included even if the extension ball is smaller.
Returns:

see Tessellation.cell_index().

tessellate(points)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

class tramway.tessellation.base.Voronoi(scaler=None)

Bases: tramway.tessellation.base.Delaunay

Voronoi graph.

Voronoi explicitly represents the cell boundaries, as a Voronoi graph, on top of the Delaunay graph that connects the cell centers. It implements the construction of this additional graph using scipy.spatial.Voronoi. This default implementation is lazy. If vertices and ridges are available, they are stored in private attributes _vertices, _vertex_adjacency and _cell_vertices. Otherwise, when vertices, vertex_adjacency or cell_vertices properties are called, the attributes are transparently made available calling the _postprocess() private method. Memory space can thus be freed again, setting vertices, vertex_adjacency and cell_vertices to None. Note however that subclasses may override these on-time calculation mechanics.

vertices

coordinates of the Voronoi vertices.

Type:numpy.ndarray
vertex_adjacency

adjacency matrix for Voronoi vertices.

Type:scipy.sparse
cell_vertices

mapping of cell indices to their associated vertices as indices in vertices.

Type:dict of array-like
cell_volume

cell volume (or surface area in 2D); only 2D is supported for now. Important note: if time is treated just as an extra dimension like in TimeLattice, then cell volume may actually be spatial volume multiplied by segment duration.

Type:numpy.ndarray
cell_adjacency

Square cell adjacency matrix. If adjacency_label is defined, cell_adjacency is sparse and the explicit elements are indices in adjacency_label.

cell_vertices
cell_volume
delete_cells(cell_indices, adjacency_label=True, pack_indices=True, _delaunay_adjacency=False, exclude_neighbours=False)

Delete cells.

Both the Delaunay and Voronoi graphs are modified.

As of version 0.4.*, all not-Delaunay-compatible adjacency links are discarded anyway.

Parameters:
  • cell_indices (numpy.ndarray) – indices of the cells to delete.
  • adjacency_label (scalar) – label for the newly-adjacent cells if adjacency labels are defined; for integer labels, True is translated as label_max+1, and False as label_min-1; passing None prevents any extra adjacency link.
  • pack_indices (bool) – cell indices are shifted down.
Returns:

index mapping (useful if pack_indices is True).

Return type:

numpy.ndarray

See also: pack_indices().

pack_indices(_delete_cell=True, _delete_vertex=True)
vertex_adjacency
vertices

Unscaled coordinates of the Voronoi vertices (numpy.ndarray).

tramway.tessellation.base.format_cell_index(K, format=None, select=None, shape=None, copy=False, **kwargs)

Convert from any valid index format to any other.

Converting an array index to any other format assumes that the point indices are in a contiguous range from 0 to the number of elements in the index.

Parameters:
  • K (any) – original point-cell association representation.
  • format (str) – either array, pair, matrix, coo, csr or csc. See also Tessellation.cell_index().
  • select (callable) – called only if format == 'array' and points are associated to multiple cells; select takes the point index as first argument, the corresponding cell indices (numpy.ndarray) as second argument and the extra keyword arguments given to format_cell_index().
  • shape (int, int) – number of points, number of cells.
  • copy (bool) – if True, ensures that a copy of K is returned if K is already in the requested format.
Returns:

point-cell association in the requested format.

Return type:

any

See also Tessellation.cell_index() and nearest_cell().

tramway.tessellation.base.nearest_cell(locations, cell_centers)

Generate a function suitable for use as format_cell_index()’s argument select.

The returned function takes a point index and cell indices as arguments and returns the index of the nearest cell.

Parameters:
  • locations (numpy.ndarray) – location coordinates.
  • cell_centers (numpy.ndarray) – cell center coordinates.
Returns:

select function.

Return type:

callable

tramway.tessellation.base.dict_to_sparse(cell_vertex, shape=None)

Convert cell-vertex association dict to sparse matrices.

tramway.tessellation.base.sparse_to_dict(cell_vertex)

Convert cell-vertex associations sparse matrices to dict.

tramway.tessellation.base.cell_index_by_radius(tessellation, points, radius, format=None, select=None, metric='euclidean', **kwargs)

See Delaunay.cell_index().

Specialized routine to assign locations to cells which center is no further than radius.

tramway.tessellation.grid plugin

class tramway.tessellation.grid.RegularMesh(scaler=None, lower_bound=None, upper_bound=None, count_per_dim=None, min_probability=None, max_probability=None, avg_probability=None, min_distance=None, avg_distance=None, **kwargs)

Bases: tramway.tessellation.base.Voronoi

Regular k-D grid.

lower_bound

(scaled)

Type:pandas.Series or numpy.ndarray
upper_bound

(scaled)

Type:pandas.Series or numpy.ndarray
count_per_dim
Type:list of ints
min_probability

(not used)

Type:float
avg_probability
Type:float
max_probability

(not used)

Type:float
min_distance

minimum distance between adjacent cell centers; ignored if avg_distance is defined.

Type:float
avg_distance

average distance between adjacent cell centers; ignored if avg_probability is defined.

Type:float
cell_adjacency

Square cell adjacency matrix. If adjacency_label is defined, cell_adjacency is sparse and the explicit elements are indices in adjacency_label.

cell_centers

Unscaled coordinates of the cell centers (numpy.ndarray).

contour(cell, distance=1, **kwargs)

See also tramway.feature.adjacency.contour().

tessellate(points, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

vertices

Unscaled coordinates of the Voronoi vertices (numpy.ndarray).

tramway.tessellation.gwr plugin

class tramway.tessellation.gwr.GasMesh(scaler=None, min_distance=None, avg_distance=None, max_distance=None, min_probability=None, avg_probability=None, **kwargs)

Bases: tramway.tessellation.base.Voronoi

GWR based tessellation.

gas

internal graph representation of the gas.

Type:Gas
min_probability

minimum probability of a point to be in any given cell.

Type:float
_min_distance

scaled minimum distance between adjacent cell centers.

Type:float
_avg_distance

upper bound on the average scaled distance between adjacent cell centers.

Type:float
_max_distance

scaled maximum distance between adjacent cell centers.

Type:float
freeze()

Delete the data required to and only to incrementally update the tessellation.

This may save large amounts of memory but differs from unload() in that the subsequent data loss may not be undone.

tessellate(points, pass_count=(), residual_factor=0.7, error_count_tol=0.005, min_growth=0.0001, collapse_tol=0.01, stopping_criterion=0, verbose=False, plot=False, alpha_risk=1e-15, grab=None, max_frames=None, max_batches=None, axes=None, complete_delaunay=False, topology='approximate density', **kwargs)

Grow the tessellation.

Parameters:
  • points – see tessellate().
  • pass_count (float or pair of floats) – number of points to sample (with replacement) from data points, as a multiple of the size of the data. If pass_count is a pair of numbers, they are the lower and the upper bounds on the number of samples. If pass_count is a single number, it is interpreted as the lower bound, and the upper bound is set equal to 2 * pass_count.
  • residual_factor (float) – multiplies with _max_distance to determine residual_max in train().
  • error_count_tol (float) – (see train())
  • min_growth (float) – (see train())
  • collapse_tol (float) – (see train())
  • stopping_criterion (int) – (see train())
  • verbose (bool) – verbose output.
  • batch_size (int) – (see Gas)
  • tau (float) – (see Gas)
  • trust (float) – (see Gas)
  • lifetime (int) – (see Gas)
  • alpha_risk (float) – location distributions of potential neighbor cells are compared with a t-test
  • complete_delaunay (bool) – complete the Delaunay graph
  • topology (str) – any of ‘approximate density’ (default), ‘displacement length’
Returns:

See tessellate().

See also

tramway.tessellation.gwr.gas.Gas and tramway.tessellation.gwr.gas.Gas.train().

tramway.tessellation.hexagon plugin

class tramway.tessellation.hexagon.HexagonalMesh(scaler=None, tilt=0.0, min_probability=None, max_probability=None, avg_probability=None, min_distance=None, avg_distance=None, **kwargs)

Bases: tramway.tessellation.base.Voronoi

2D hexagonal mesh.

tilt

angle of the “main” axis, in pi/6 radians; 0= the main axis is the x-axis; 1= the main axis is the y-axis.

Type:float
hexagon_radius

hexagon radius.

Type:float
hexagon_count

number of hexagons along the “main” axis, number of hexagons along the perpendical axis.

Type:pair of ints
min_probability

(not used)

Type:float
avg_probability

inverse of the average number of points per cell.

Type:float
max_probability

(not used)

Type:float
min_distance

minimum distance between adjacent cell centers; ignored if avg_distance is defined.

Type:float
avg_distance

average distance between adjacent cell centers; ignored if avg_probability is defined.

Type:float
avg_distance
avg_probability
cell_adjacency

Square cell adjacency matrix. If adjacency_label is defined, cell_adjacency is sparse and the explicit elements are indices in adjacency_label.

cell_centers

Unscaled coordinates of the cell centers (numpy.ndarray).

cell_volume
hexagon_count
hexagon_radius
max_probability
min_distance
min_probability
tessellate(points, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

tilt
vertices

Unscaled coordinates of the Voronoi vertices (numpy.ndarray).

tramway.tessellation.kdtree plugin

class tramway.tessellation.kdtree.KDTreeMesh(scaler=None, min_distance=None, avg_distance=None, min_probability=None, max_probability=None, max_level=None, **kwargs)

Bases: tramway.tessellation.base.Voronoi

k-dimensional tree (quad tree in 2D) based tessellation.

scaler

see Tessellation.

min_probability

minimum probability of a point to be in a given cell.

Type:float
max_probability

maximum probability of a point to be in a given cell.

Type:float
max_level

maximum level, considering that the smallest cells are at level 0 and the level increments each time the cell size doubles.

Type:float
_min_distance

scaled minimum distance between neighbor cell centers.

Type:float, private
_avg_distance

scaled average distance between neighbor cell centers.

Type:float, private
cell_index(points, format=None, select=None, knn=None, min_location_count=None, metric='chebyshev', filter=None, filter_descriptors_only=False, **kwargs)

See Tessellation.cell_index().

A single array representation of the point-cell association may not be possible with knn defined, because a point can be associated to multiple cells. If such a case happens the default output format will be pair.

In addition to the values allowed by Tessellation.cell_index(), format admits value force array that acts like format='array', select=nearest_cell(...). The implementation however is more straight-forward and simply ignores the minimum number of nearest neighbours if provided.

new in 0.6: fast_min_nn fixes a bug and should be set to False; default is True for backward-compatibility.

Parameters:
  • points – see Tessellation.cell_index().
  • format – see Tessellation.cell_index(); additionally admits force array.
  • select – see Tessellation.cell_index().
  • knn (int or tuple or callable) – If int: minimum number of points per cell (or of nearest neighours to the cell center). Cells may overlap and the returned cell index may be a sparse point-cell association. If tuple (pair of ints): minimum and maximum number of points per cell respectively. If callable: takes a cell index and returns the minimum and maximum number of points.
  • radius (float or tuple or callable) – If float: distance from the cell center; smaller cells may include locations from neighbour cells and larger cells may include only part of their associated locations. If tuple (pair of floats): minimum and maximum radius of a cell respectively; any of these values can be None. If callable: takes a cell index and returns the minimum and maximum radius.
  • min_location_count (int) – minimum number of points for a cell to be included in the labeling. This argument applies before knn. The points in these cells, if not associated with another cell, are labeled -1. The other cell labels do not change.
  • metric (str) – any metric name understandable by cdist().
  • filter (callable) – takes the calling instance, a cell index and the corresponding subset of points; returns True if the corresponding cell should be included in the labeling.
  • filter_descriptors_only (bool) – whether filter should get points as descriptors only.
  • fast_min_nn (bool) – new in 0.6 if True, the initial assignment for small cells is overriden by the extension ball; if False, the initial assignment is included even if the extension ball is smaller.
Returns:

see Tessellation.cell_index().

delete_cell(i, adjacency_label=True, metric='chebyshev', pack_indices=True)
freeze()

Delete the data required to and only to incrementally update the tessellation.

This may save large amounts of memory but differs from unload() in that the subsequent data loss may not be undone.

tessellate(points, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

tramway.tessellation.kmeans plugin

class tramway.tessellation.kmeans.KMeansMesh(scaler=None, min_probability=None, avg_probability=None, min_distance=None, initial='grid', **kwargs)

Bases: tramway.tessellation.base.Voronoi

K-Means and Voronoi based tessellation.

avg_probability

probability of a point to be in a given cell (controls the number of cells and indirectly their size).

Type:float

Other Attributes:

_min_distance (float): scaled minimum distance between adjacent cell centers;
not used.
tessellate(points, tol=1e-06, prune=2.5, plot=False, **kwargs)

Grow the tessellation.

points

see tessellate().

tol

error tolerance. Passed as thresh to scipy.cluster.vq.kmeans().

Type:float
prune

prunes the Voronoi and removes the edges which length is greater than prune times the median edge length; True is translated to the default value.

Type:bool or float

tramway.tessellation.nesting module

class tramway.tessellation.nesting.NestedTessellations(scaler=None, parent=None, factory=None, parent_index_arguments={}, **kwargs)

Bases: tramway.tessellation.base.Tessellation

Tessellation of tessellations.

When nesting, the parent tessellation should have been grown already.

tessellation grows all the nested tessellations.

In __init__, tessellation and cell_index, the scaler (if any), args (if any) and kwargs arguments only apply to the nested tessellations.

adjacency_label

Inter-cell edge labels, numpy.ndarray with as many elements as edges.

cell_adjacency

Square cell adjacency matrix. If adjacency_label is defined, cell_adjacency is sparse and the explicit elements are indices in adjacency_label.

cell_centers
cell_index(points, *args, **kwargs)

Partition.

The returned value depends on the format input argument:

  • array: returns a vector v such that v[i] is cell index for
    point index i or -1.
  • pair: returns a pair of I-sized arrays (p, c) where, for each
    point-cell association i in range(I), p[i] is a point index and c[i] is a corresponding cell index.
  • matrix or coo or csr or csc:
    returns a sparse matrix with points as rows and cells as columns; non-zeros are all True or float weights.

By default with format undefined, any implementation may favor any format.

Note that array may not be an acceptable format and cell_index() may not comply with format='index' unless select is defined. When a location or a translocation is associated to several cells, select chooses a single cell among them.

The default implementation calls format_cell_index() on the result of an abstract _cell_index method that any Tessellation implementation can implement instead of cell_index().

See also format_cell_index().

Parameters:
  • points (pandas.DataFrame) – point (location) coordinates.
  • format (str) – preferred representation of the point-cell association (or partition).
  • select (callable) – takes the point index, an array of cell indices and the tessellation as arguments, and returns a cell index or -1 for no cell.
cell_label

Cell labels, numpy.ndarray with as many elements as there are cells.

child_cell_indices(u)
Parameters:u (int) – child index.
Returns:child’s cell indices.
Return type:slice
child_factory
child_factory_arguments
children
freeze()

Delete the data required to and only to incrementally update the tessellation.

This may save large amounts of memory but differs from unload() in that the subsequent data loss may not be undone.

parent
parent_index_arguments
tessellate(points, *args, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

tramway.tessellation.random plugin

class tramway.tessellation.random.RandomMesh(scaler=None, avg_probability=None, lower_bound=None, upper_bound=None, cell_count=None, **kwargs)

Bases: tramway.tessellation.base.Voronoi

avg_probability
cell_centers

Unscaled coordinates of the cell centers (numpy.ndarray).

lower_bound
tessellate(points, allow_empty_cells=False, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

upper_bound

tramway.tessellation.time module

class tramway.tessellation.time.TimeLattice(scaler=None, segments=None, time_label=None, mesh=None, time_dimension=None)

Bases: tramway.tessellation.base.Tessellation

Proxy Tessellation for time lattice expansion.

If time_lattice contains integers, these elements are regarded as frame indices, whereas if it contains floats, the elements represent time.

The time_edge attribute (time_label init argument) drives the encoding of temporal adjacency. It is a two-element sequence respectively representing past and future relationships. If time_label is supplied as a single scalar value, past and future relationships are encoded with this same label. If a label is False, then the corresponding relationship is not represented. If a label is True, then it is translated into an integer value that is not used yet.

The time_dimension attribute may be useful combined with a defined spatial_mesh to include time in the representation of cell centers and the calculation of cell volumes.

If time_label is not None and time_dimension is None, time_dimension will default to True. If time_dimension is True and time_label is None, then time_label will default to True. None to both arguments will be treated as False. These rules are resolved in tessellate().

Functional dependencies:

  • setting time_lattice unsets cell_label, cell_adjacency and adjacency_label
  • setting spatial_mesh unsets cell_centers, cell_label, cell_adjacency and adjacency_label, cell_volume
  • cell_centers, cell_volume and split_frames are available only when spatial_mesh
    is defined
adjacency_label

Inter-cell edge labels, numpy.ndarray with as many elements as edges.

cell_adjacency

Square cell adjacency matrix. If adjacency_label is defined, cell_adjacency is sparse and the explicit elements are indices in adjacency_label.

cell_centers
cell_index(points, *args, **kwargs)

In addition to the arguments to cell_index() for the underlying spatial tessellation:

Parameters:
  • time_col (int or str) – column name; default is ‘t’; keyword-argument only.
  • time_knn (int or tuple) – min_nn minimum number of nearest “neighbours” of the time segment center, or (min_nn, max_nn) minimum and maximum numbers of nearest neighbours in time; keyword-argument only.

knn and time_knn are not compatible.

cell_label

Cell labels, numpy.ndarray with as many elements as there are cells.

cell_volume

If time_dimension is True, then the product of spatial volume and segment duration is returned. Otherwise, cell_volume is only the spatial volume.

descriptors(points, *args, **kwargs)

Keep the data columns that were involved in growing the tessellation.

Parameters:
  • points (pandas.DataFrame) – point coordinates.
  • asarray (bool) – returns a numpy.ndarray.
Returns:

selected coordinates; the data may not be copied.

Return type:

array-like

freeze()

Delete the data required to and only to incrementally update the tessellation.

This may save large amounts of memory but differs from unload() in that the subsequent data loss may not be undone.

future_edge
number_of_cells

Number of cells.

Read-only property.

past_edge
simplified_adjacency(adjacency=None, label=None, format='coo', distinguish_time=None)

distinguish_time allows to keep temporal relationships distinct from the spatial relationships. As a consequence of encoding time, the simplified adjacency matrix is not boolean.

If distinguish_time is None, then distinguish_time will default to True if time_edge is defined, False otherwise.

spatial_mesh
split_frames(df, return_times=False)
split_segments(spt_data, return_times=False)
tessellate(points, time_only=False, **kwargs)

Grow the tessellation.

Parameters:points (pandas.DataFrame) – point coordinates.

Admits keyword arguments.

time_dimension
time_edge
time_lattice
tramway.tessellation.time.with_time_lattice(cells, frames, exclude_cells_by_location_count=None, **kwargs)

tramway.tessellation.window plugin

class tramway.tessellation.window.SlidingWindow(scaler=None, duration=None, shift=None, frames=False, time_label=None, time_dimension=None, start_time=None)

Bases: tramway.tessellation.time.TimeLattice

cell_index(points, *args, **kwargs)

In addition to the arguments to cell_index() for the underlying spatial tessellation:

Parameters:
  • time_col (int or str) – column name; default is ‘t’; keyword-argument only.
  • time_knn (int or tuple) – min_nn minimum number of nearest “neighbours” of the time segment center, or (min_nn, max_nn) minimum and maximum numbers of nearest neighbours in time; keyword-argument only.

knn and time_knn are not compatible.

duration
shift
start_time