Core objects

The Epithelium object

class tyssue.Epithelium(identifier, datasets, specs=None, coords=None, maxbackup=5)

Base class defining a connective tissue in 2D or 3D.

property Nc

The number of cells in the epithelium.

property Ne

The number of edges in the epithelium.

property Nf

The number of faces in the epithelium.

property Nv

The number of vertices in the epithelium.

backup()

Creates a copy of self and keeps a reference to it in the self._backups deque.

property cell_df

The cell pd.DataFrame containing cell associated data e.g. the position of their center or their volume

copy(deep_copy=True)

Returns a copy of the epithelium

deep_copybool, default True

if True, use a copy of the original object’s datasets to create the new object. If False, datasets are not copied

cut_out(bbox, coords=None)

Returns the index of edges with at least one vertex outside of the bounding box

bboxsequence of shape (dim, 2)

the bounding box as (min, max) pairs for each coordinates.

coordslist of str of len dim, default None

the coords corresponding to the bbox.

property edge_df

The edge pd.DataFrame containing edge associated data e.g. their length.This dataframe also contains the whole connexion of the epithelium through the “srce”, “trgt”, “face”, “cell” indices. In 2D, a half-edge is associated with a single (face, srce, trgt) positively oriented triangle. In 3D, the (cell, face, srce, trgt) positively oriented terahedron is also unique.

property face_df

The face pd.DataFrame containing face associated data e.g. the position of their center or their area

face_polygons(coords=None)

Returns a pd.Series of arrays with the coordinates the face polygons

Each element of the Series is a (num_sides, num_dims) array of points ordered counterclockwise.

Vertices are assumed to be ordered in a face. If you are not sure it is the case, you can run sheet.reset_index(order=True) before calling this function.

get_invalid()

Returns a mask over self.edge_df for invalid faces.

get_neighborhood(elem_id, order, elem='cell')

Returns elem_id neighborhood up to a degree of order

For example, if order is 2, it wil return the adjacent cells (or faces) and theses cells neighbors.

neighborspd.DataFrame with two colums, the index

of the neighboring cell (face), and it’s neighboring order

get_neighbors(elem_id, elem='cell')

Returns the indexes of the adjacent elements (cells or faces) of the element of index elem_id.

elem_idint

the index of the central element (a face or a cell)

element : {‘cell’ | ‘face’}, default ‘cell’

neghborsset

the cells (or faces) sharing an edge with the central cell (face)

get_opposite_faces()

Populates the ‘opposite’ column of self.face_df with the index of the opposite face or -1 if the face has no opposite.

get_orbits(center, periph)

Returns a dataframe with a (center, edge) MultiIndex with periph elements.

centerstr,

the name of the center element for example ‘face’, ‘srce’

periphstr,

the name of the periphery elements, for example ‘trgt’, ‘cell’

>>> cell_verts = sheet.get_orbits('face', 'srce')
>>> cell_verts.loc[45]
edge
218    75
219    78
220    76
221    81
222    90
223    87
Name: srce, dtype: int64
get_valid()

Set the ‘is_valid’ column to true if the faces are all closed polygons, and the cells closed polyhedra.

idx_lookup(elem_id, element)

returns the current index of the element with the “id” column equal to elem_id

elem_idint

id of the element to retrieve

element{“vert”|”edge”|”face”|”cell”}

the corresponding dataset.

remove(edge_out, trim_borders=False, order_edges=False)

Remove the edges indexed by edge_out associated with all the cells and faces containing those edges

If trim_borders is True (defaults to False), there will be a single border edge per border face.

reset_index(order=False)

Resets the datasets to have continuous indices

If order is True (the default), sorts the edges such that for each face, vertices are ordered clockwize

reset_topo()

Recomputes the number of sides for the faces and the number of faces for the cells.

restore()

Resets the eptithelium data to its last backed up state

A copy of the current state prior to restoring is kept in the _bad attribute for inspection. Calling this method multiple times (without calling backup) will go back in the epithelium backups.

sanitize(trim_borders=False, order_edges=False)

Removes invalid faces and associated vertices

If trim_borders is True (defaults to False), there will be a single border edge per border face.

set_bbox(margin=0.0)

Sets the attribute bbox with pairs of values bellow and above the min and max of the vert coords, with a margin.

property settings

Accesses the specs[‘settings’] dictionnary.

sum_cell(df)

Sums the values of the edge-indexed dataframe df grouped by the values of self.edge_df[“cell”]

summed : pd.DataFrame the summed data, indexed by the source vertices.

sum_face(df)

Sums the values of the edge-indexed dataframe df grouped by the values of self.edge_df[“face”]

summed : pd.DataFrame the summed data, indexed by the source vertices.

sum_srce(df)

Sums the values of the edge-indexed dataframe df grouped by the values of self.edge_df[“srce”]

summed : pd.DataFrame the summed data, indexed by the source vertices.

sum_trgt(df)

Sums the values of the edge-indexed dataframe df grouped by the values of self.edge_df[“trgt”]

summed : pd.DataFrame the summed data, indexed by the source vertices.

triangular_mesh(coords=None, return_mask=False)

Return a triangulation of an epithelial sheet (2D in a 3D space), with added edges between face barycenters and junction vertices.

coordslist of str:

pair of coordinates corresponding to column names for self.face_df and self.vert_df

return_maskbool, optional, default True

if True, returns face_mask

vertices(self.Nf+self.Nv, 3) ndarray

all the vertices’ coordinates

triangles(self.Ne, 3) ndarray of ints

triple of the vertices’ indexes forming the triangular elements. For each junction edge, this is simply the index (srce, trgt, face). This is correctly oriented.

face_mask: (self.Nf + self.Nv,) mask with 1 iff the vertex corresponds

to a face center

upcast_cell(df)

Reindexes input data to self.edge_df.index by repeating the values for each cell entry

dfpd.DataFrame, pd.Series np.ndarray or string

The data to be upcasted. If array like, should have self.Nc elements. If a string is passed it should be a column of self.cell_df

upcast_dfpd.DataFrame or pd.Series

The value repeated like the values of self.edge_df[“cell”]

upcast_cols(element, columns)

Syntactic sugar to upcast from the epithelium datasets.

element: {‘srce’|’trgt’|’face’|’cell’}

corresponding self.edge_df column over which to index if element is ‘srce’ or ‘trgt’, the upcast data will be taken form self.vert_df

columns: index

the column(s) to be taken from the input dataset.

upcast_face(df)

Reindexes input data to self.edge_df.index by repeating the values for each face entry

dfpd.DataFrame, pd.Series np.ndarray or string

The data to be upcasted. If array like, should have self.Nf elements. If a string is passed it should be a column of self.face_df

upcast_dfpd.DataFrame or pd.Series

The value repeated like the values of self.edge_df[“face”]

upcast_srce(df)

Reindexes input data to self.edge_df.index by repeating the values for each source entry.

dfpd.DataFrame, pd.Series np.ndarray or string

The data to be upcasted. If array like, should have self.Nv elements. If a string is passed it should be a column of self.vert_df

upcast_dfpd.DataFrame or pd.Series

The value repeated like the values of self.edge_df[“srce”]

upcast_trgt(df)

Reindexes input data to self.edge_df.index by repeating the values for each target entry

dfpd.DataFrame, pd.Series np.ndarray or string

The data to be upcasted. If array like, should have self.Nv elements. If a string is passed it should be a column of self.vert_df

upcast_dfpd.DataFrame or pd.Series

The value repeated like the values of self.edge_df[“trgt”]

update_num_faces()

Updates the number of faces around the cells. The data is registered in the “num_faces” column of self.cell_df.

update_num_sides()

Updates the number of half-edges around the faces. The data is registered in the “num_sides” column of self.face_df.

update_rank()
update_specs(new, reset=False)

Recursively updates the self.specs nested dictionnary, and set the new values to the corresponding columns in the datasets. If reset is True, existing values will be overwritten.

validate()

returns True if the mesh is validated

e.g. has only closed polygons and polyhedra

validate_closed_cells()

Returns True if all cells of the epithelium are closed.

property vert_df

The face pd.DataFrame containing vertex associated data e.g. the position of each vertex.

vertex_mesh(coords, vertex_normals=True)

Returns the vertex coordinates and a list of vertex indices for each face of the tissue. If vertex_normals is True, also returns the normals of each vertex (set as the average of the vertex’ edges), suitable for .OBJ export

Vertices are assumed to be ordered in a face. If you are not sure it is the case, you can run sheet.reset_index() before calling this function.

Sheet object

class tyssue.Sheet(identifier, datasets, specs=None, coords=None)

An epithelial sheet, i.e a 2D mesh in a 2D or 3D space, akin to a HalfEdge data structure in CGAL.

The geometric properties are defined in tyssue.geometry.sheet_geometry A dynamical model derived from Fahradifar et al. 2007 is provided in tyssue.dynamics.sheet_vertex_model

extract(face_mask, coords=['x', 'y', 'z'])

Extract a new sheet from the sheet that correspond to a key word that define a face.

face_mask : column name in face composed by boolean value coords :

sheet_fold_patch_extract :

subsheet corresponding to the fold patch area.

extract_bounding_box(x_boundary=None, y_boundary=None, z_boundary=None, coords=['x', 'y', 'z'])

Extracts a new sheet from the embryo sheet

that correspond to boundary coordinate define by the user.

x_boundary : pair of floats y_boundary : pair of floats z_boundary : pair of floats coords : list of strings, default [‘x’, ‘y’, ‘z’]

coordinates over which to crop the sheet

subsheet : a new Sheet object

get_extra_indices()

Computes extra indices:

  • self.free_edges: half-edges at the epithelium boundary

  • self.dble_edges: half-edges inside the epithelium, with an opposite

  • self.east_edges: half of the dble_edges, pointing east (figuratively)

  • self.west_edges: half of the dble_edges, pointing west

    (the order of the east and west edges is conserved, so that the ith west half-edge is the opposite of the ith east half-edge)

  • self.sgle_edges: joint index over free and east edges, spanning

    the entire graph without double edges

  • self.wrpd_edges: joint index over free edges followed by the
    east edges twice, such that a vector over the whole half-edge

    dataframe is wrapped over the single edges

  • self.srtd_edges: index over the whole half-edge sorted such that

    the free edges come first, then the east, then the west

Also computes: - self.Ni: the number of inside full edges

(i.e. len(self.east_edges))

  • self.No: the number of outside full edges (i.e. len(self.free_edges))

  • self.Nd: the number of double half edges (i.e. len(self.dble_edges))

  • self.anti_sym: pd.Series with shape (self.Ne,) with 1 at the free and east half-edges and -1 at the opposite half-edges.

  • East and west is resepctive to some orientation at the moment the indices are computed the partition stays valid as long as there are no changes in the topology, so due to vertex displacement, ‘east’ and ‘west’ might not stay valid. This is just a practical naming convention.

  • As the name suggest, this method is not working for edges in 3D pointing exactly north or south, ie iff edge[‘dx’] == edge[‘dy’] == 0. Until we need or find a better solution, we’ll just assert it worked.

get_force_inference(coords=None, column=None, free_border_edges=False)

Measure force based on Brodland method.

g_gamma_matrix*tension_vector = 0 Equation is homogenous and to avoid tension_vectors = 0, Construction and solve the constrained least-squares equation system [[g_gamma_matrix.T g_gamma_matrix, C^T_1],

[C_1, 0]] où C1 = {1….1}

shape of g_gamma_matrix = (Ne/2, Nv*len(coords))

Note

Results might not be consistens for highly curved epithelium

coords: coordinates column: None, specify a column name in edge_df to put tension value free_border_edges: bool, default False, take into account edges in the border of the tissue if True

edges_tensions: tension values array if column not define

get_neighborhood(face, order, elem='face')

Returns face neighborhood up to a degree of order.

For example, if order is 2, it wil return the adjacent, faces and theses faces neighbors.

neighborspd.DataFrame with two colums, the index

of the neighboring face, and it’s neighboring order

get_neighbors(face, elem='face')

Returns the faces adjacent to face.

get_opposite()
classmethod planar_sheet_2d(identifier, nx, ny, distx, disty, noise=None)

Creates a planar sheet from an hexagonal grid of cells.

identifier : string nx, ny : int

number of cells in the x and y axes

distx, distyfloat,

the distances in x and y between the cells

noisefloat, default None

position noise on the hexagonal grid

planar_sheet: a 2D Sheet instance

in the (x, y) plane

classmethod planar_sheet_3d(identifier, nx, ny, distx, disty, noise=None)

Creates a planar sheet from an hexagonal grid of cells.

identifier : string nx, ny : int

number of cells in the x and y axes

distx, distyfloat,

the distances in x and y between the cells

noisefloat, default None

position noise on the hexagonal grid

flat_sheet: a 2.5D Sheet instance

reset_topo()

Recomputes the number of sides for the faces and the number of faces for the cells.

sort_edges_eastwest()

reorder edges such the free edges are first, then the first half of the double edges, then the other half of the double edges, this way, each subset of the edges dataframe are contiguous.

Monolayer object

class tyssue.Monolayer(name, datasets, specs=None, coords=None)

3D monolayer epithelium

property apical_edges
property apical_faces
property apical_verts
property basal_edges
property basal_faces
property basal_verts
classmethod from_flat_sheet(name, apical_sheet, specs, thickness=1)
get_sub_sheet(segment)

Returns a Sheet object of the corresponding segment

segment: str, the corresponding segment, wether ‘apical’ or ‘basal’

guess_face_segment(face)

Infers the face segment.

guess_vert_segment(vert)

Infers the vertex segment from its surrounding edges.

property lateral_edges
property lateral_faces
property lateral_verts
segment_index(segment, element)

Drawing functions

Matplotlib based

Matplotlib based plotting

tyssue.draw.plt_draw.create_gif(history, output, num_frames=None, interval=None, draw_func=None, margin=5, **draw_kwds)

Creates an animated gif of the recorded history.

You need imagemagick on your system for this function to work.

history : a tyssue.History object output : path to the output gif file num_frames : int, the number of frames in the gif interval : tuples, define begin and end frame of the gif draw_func : a drawing function

this function must take a sheet object as first argument and return a fig, ax pair. Defaults to quick_edge_draw (aka sheet_view with quick mode)

marginint, the graph margins in percents, default 5

if margin is -1, let the draw function decide

**draw_kwds are passed to the drawing function

tyssue.draw.plt_draw.curved_view(sheet, radius_cutoff=1000.0)
tyssue.draw.plt_draw.draw_edge(sheet, coords, ax, **draw_spec_kw)
tyssue.draw.plt_draw.draw_face(sheet, coords, ax, **draw_spec_kw)

Draws epithelial sheet polygonal faces in matplotlib Keyword values can be specified at the element level as columns of the sheet.face_df

tyssue.draw.plt_draw.draw_vert(sheet, coords, ax, **draw_spec_kw)

Draw junction vertices in matplotlib.

tyssue.draw.plt_draw.get_arc_data(sheet)
tyssue.draw.plt_draw.parse_face_specs(face_draw_specs, sheet)
tyssue.draw.plt_draw.plot_forces(sheet, geom, model, coords, scaling, ax=None, approx_grad=None, **draw_specs_kw)

Plot the net forces at each vertex, with their amplitudes multiplied by scaling. To be clear, this is the oposite of the gradient - grad E.

tyssue.draw.plt_draw.plot_junction(eptm, edge_index, coords=['x', 'y'])

Plots local graph around a junction, for debugging purposes.

tyssue.draw.plt_draw.plot_scaled_energies(sheet, geom, model, scales, ax=None)

Plot scaled energies

sheet: a:class: Sheet object geom: a Geometry class model: a :class:’Model’ scales: np.linspace of float

fig: a :class:matplotlib.figure.Figure instance ax: :class:matplotlib.Axes instance, default None

tyssue.draw.plt_draw.quick_edge_draw(sheet, coords=['x', 'y'], ax=None, **draw_spec_kw)
tyssue.draw.plt_draw.sheet_view(sheet, coords=['x', 'y'], ax=None, **draw_specs_kw)

Base view function, parametrizable through draw_secs

The default sheet_spec specification is:

{‘edge’: {

‘visible’: True, ‘width’: 0.5, ‘head_width’: 0.2, # arrow head width for the edges ‘length_includes_head’: True, # see matplotlib Arrow artist doc ‘shape’: ‘right’, ‘color’: ‘#2b5d0a’, # can be an array ‘alpha’: 0.8, ‘zorder’: 1, ‘colormap’: ‘viridis’},

‘vert’: {

‘visible’: True, ‘s’: 100, ‘color’: ‘#000a4b’, ‘alpha’: 0.3, ‘zorder’: 2},

‘face’: {

‘visible’: False, ‘color’: ‘#8aa678’, ‘alpha’: 1.0, ‘zorder’: -1} }

Ipyvolume based

3D visualisation inside the notebook.

tyssue.draw.ipv_draw.browse_history(history, coords=['x', 'y', 'z'], **draw_specs_kw)
tyssue.draw.ipv_draw.edge_mesh(sheet, coords, **edge_specs)

Creates a ipyvolume Mesh of the edge lines to be displayed in Jupyter Notebooks

mesh: a ipyvolume.widgets.Mesh mesh widget

tyssue.draw.ipv_draw.face_mesh(sheet, coords, **face_draw_specs)

Creates a ipyvolume Mesh of the face polygons

tyssue.draw.ipv_draw.sheet_view(sheet, coords=['x', 'y', 'z'], **draw_specs_kw)

Creates a javascript renderer of the edge lines to be displayed in Jupyter Notebooks

fig: a ipyvolume.widgets.Figure widget mesh: a ipyvolume.widgets.Mesh mesh widget

tyssue.draw.ipv_draw.update_view(fig, meshes)
tyssue.draw.ipv_draw.view_ipv(sheet, coords=['x', 'y', 'z'], **edge_specs)

Creates a javascript renderer of the edge lines to be displayed in Jupyter Notebooks

fig: a ipyvolume.widgets.Figure widget mesh: a ipyvolume.widgets.Mesh mesh widget

Geometry classes

Planar

class tyssue.PlanarGeometry

Geomtetry methods for 2D planar cell arangements

static face_projected_pos(sheet, face, psi)

returns the sheet vertices position translated to center the face face at (0, 0) and rotated in the (x, y) plane by and angle psi radians

classmethod get_phis(sheet)
classmethod update_all(sheet)

Updates the sheet geometry by updating: * the edge vector coordinates * the edge lengths * the face centroids * the normals to each edge associated face * the face areas

static update_areas(sheet)

Updates the normal coordinate of each (srce, trgt, face) face.

static update_normals(sheet)

Sheet (2D 1/2)

class tyssue.SheetGeometry

Geometry definitions for 2D sheets in 3D

static face_projected_pos(sheet, face, psi=0)

Returns the position of a face vertices projected on a plane perpendicular to the face normal, and translated so that the face center is at the center of the coordinate system

sheet: a :class:Sheet object face: int,

the index of the face on which to rotate the sheet

psi: float,

Optional angle giving the rotation along the z axis

rot_pos: pd.DataFrame

The rotated, relative positions of the face’s vertices

static face_rotation(sheet, face, psi=0)

Returns a 3D rotation matrix such that the face normal points in the z axis

sheet: a :class:Sheet object face: int,

the index of the face on which to rotate the sheet

psi: float,

Optional angle giving the rotation along the new z axis

rotation: (3, 3) np.ndarray

The rotation matrix

classmethod face_rotations(sheet, method='normal', output_as='edge')

Returns the (sheet.Ne, 3, 3) array of rotation matrices such that each rotation returns a coordinate system (u, v, w) where the face vertices are mostly in the u, v plane.

If method is ‘normal’, face is oriented with it’s normal along w if method is ‘svd’, the u, v, w is determined through singular value decompostion of the face vertices relative positions.

svd is slower but more effective at reducing face dimensionality.

output_as: string, default ‘edge’ Return the (sheet.Ne, 3, 3),

else ‘face’ Return the (sheet.Nf, 3, 3)

classmethod get_phis(sheet, method='normal')

Returns the angle of the vertices in the plane perpendicular to each face normal. For not-too-deformed faces, sorting vertices by this gives clockwize orientation.

I think not-too-deformed means starconvex here.

The ‘method’ argument is passed to face_rotations

static normal_rotations(sheet, output_as='edge')

Returns the (sheet.Ne, 3, 3) array of rotation matrices such that each rotation aligns the coordinate system along each face normals

output_as: string, default ‘edge’ Return the (sheet.Ne, 3, 3),

else ‘face’ Return the (sheet.Nf, 3, 3)

classmethod reset_scafold(sheet)

Re-centers and (in the case of a rod sheet) resets the a-b parameters and tip masks

static svd_rotations(sheet, output_as='edge')

Returns the (sheet.Ne, 3, 3) array of rotation matrices such that each rotation aligns the coordinate system according to each face vertex SVD

output_as: string, default ‘edge’ Return the (sheet.Ne, 3, 3),

else ‘face’ Return the (sheet.Nf, 3, 3)

classmethod update_all(sheet)

Updates the sheet geometry by updating: * the edge vector coordinates * the edge lengths * the face centroids * the normals to each edge associated face * the face areas * the vertices heights (depends on geometry) * the face volumes (depends on geometry)

static update_areas(sheet)

Updates the normal coordniate of each (srce, trgt, face) face.

classmethod update_height(sheet)

Update the height of the sheet vertices, based on the geometry specified in the sheet settings:

sheet.settings[‘geometry’] can be set to

  • cylindrical: the vertex height is

    measured with respect to the distance to the the axis specified in sheet.settings[‘height_axis’] (e.g z)

  • flat: the vertex height is

    measured with respect to the position on the axis specified in sheet.settings[‘height_axis’]

  • ‘spherical’: the vertex height is measured with respect to its

    distance to the coordinate system centers

  • ‘rod’: the vertex height is measured with respect to its

    distance to the coordinate height axis if between the focii, and from the closest focus otherwise. The focii positions are updated before the height update.

In all the cases, this distance is shifted by an amount of sheet.vert_df[‘basal_shift’]

static update_normals(sheet)

Updates the face_df coords columns as the face’s vertices center of mass.

static update_vol(sheet)

Note that this is an approximation of the sheet geometry module.

class tyssue.ClosedSheetGeometry

Geometry for a closed 2.5D sheet.

Apart from the geometry update from a normal sheet, the enclosed volume is also computed. The value is stored in sheet.settings[“lumen_vol”]

classmethod update_all(sheet)

Updates the sheet geometry by updating: * the edge vector coordinates * the edge lengths * the face centroids * the normals to each edge associated face * the face areas * the vertices heights (depends on geometry) * the face volumes (depends on geometry)

static update_lumen_vol(sheet)
class tyssue.geometry.sheet_geometry.EllipsoidGeometry
static scale(eptm, scale, coords)

Scales the coordinates coords by a factor delta

static update_height(eptm)

Update the height of the sheet vertices, based on the geometry specified in the sheet settings:

sheet.settings[‘geometry’] can be set to

  • cylindrical: the vertex height is

    measured with respect to the distance to the the axis specified in sheet.settings[‘height_axis’] (e.g z)

  • flat: the vertex height is

    measured with respect to the position on the axis specified in sheet.settings[‘height_axis’]

  • ‘spherical’: the vertex height is measured with respect to its

    distance to the coordinate system centers

  • ‘rod’: the vertex height is measured with respect to its

    distance to the coordinate height axis if between the focii, and from the closest focus otherwise. The focii positions are updated before the height update.

In all the cases, this distance is shifted by an amount of sheet.vert_df[‘basal_shift’]

Bulk (3D)

class tyssue.BulkGeometry

Geometry functions for 3D cell arangements

classmethod update_all(eptm)

Updates the eptm geometry by updating: * the edge vector coordinates * the edge lengths * the face centroids * the normals to each edge associated face * the face areas * the cell areas * the vertices heights (depends on geometry) * the face volumes (depends on geometry)

static update_areas(eptm)

Updates the normal coordniate of each (srce, trgt, face) face.

static update_centroid(eptm)

Updates the face_df coords columns as the face’s vertices center of mass. Also updates the edge_df fx, fy, fz columns with their upcasted values

static update_dcoords(eptm)

Update the edge vector coordinates on the coords basis (default_coords by default). Modifies the corresponding columns (i.e [‘dx’, ‘dy’, ‘dz’]) in sheet.edge_df.

Also updates the upcasted coordinates of the source and target vertices

static update_ucoords(eptm)
static update_vol(eptm)
static validate_face_norms(eptm)
class tyssue.MonolayerGeometry
static basal_apical_axis(eptm, cell)

Returns a unit vector allong the apical-basal axis of the cell

classmethod cell_projected_pos(eptm, cell, psi=0)

Returns the positions of the cell vertices transformed such that the cell center sits at the coordinate system’s origin and the basal-apical axis is the new z axis.

class tyssue.ClosedMonolayerGeometry
classmethod update_all(eptm)

Updates the eptm geometry by updating: * the edge vector coordinates * the edge lengths * the face centroids * the normals to each edge associated face * the face areas * the cell areas * the vertices heights (depends on geometry) * the face volumes (depends on geometry)

static update_lumen_vol(eptm)

Topology functions

Base

tyssue.topology.base_topology.add_vert(eptm, edge)

Adds a vertex in the middle of the edge,

which is split as is its opposite(s)

eptm : a Epithelium instance edge : int the index of one of the half-edges to split

new_vert : int the index to the new vertex new_edges : int or list of ints index to the new edge(s). For a sheet, returns a single index, for a 3D epithelium, returns the list of all the new parallel edges new_opp_edges : int or list of ints index to the new opposite edge(s). For a sheet, returns a single index, for a 3D epithelium, returns the list of all the new parallel edges

In the simple case whith two half-edge, returns indices to the new edges, with the following convention:

s e t

——>

  • <—— *

oe

s e ne t

—— —–>

  • <—– * —— *

    oe nv noe

where “e” is the passed edge as argument, “s” its source “t” its target and “oe” its opposite. The returned edges are the ones between the new vertex and the input edge’s original target.

tyssue.topology.base_topology.close_face(eptm, face)

Closes the face if a single edge is missing.

This function does not close the adjacent and opposite faces.

tyssue.topology.base_topology.collapse_edge(sheet, edge, reindex=True, allow_two_sided=False)

Collapses edge and merges it’s vertices, creating (or increasing the rank of) a rosette structure.

If reindex is True (the default), resets indexes and topology data. The edge is collapsed on the smaller of the srce, trgt indexes (to minimize reindexing impact)

tyssue.topology.base_topology.condition_4i(eptm)

Return an index over the faces violating condition 4 i in Okuda et al 2013, that is edges (from the same face) sharing two vertices simultaneously.

tyssue.topology.base_topology.condition_4ii(eptm)

Return an array of face pairs sharing more than two half-edges, as defined in Okuda et al. 2013 condition 4 ii

An indication way to solve this:

faces = condition_4ii(eptm)

pairs = set(frozenset(p) for p in faces)


cols = ['srce', 'trgt', 'face', 'cell', 'length', 'sub_area']
edges = eptm.edge_df[eptm.edge_df["face"].isin(faces[0])][cols].sort_values("face")
all_edges = eptm.edge_df.loc[eptm.edge_df["face"].isin(set(faces.ravel())), cols].sort_values("face")

all_edges['single'] = all_edges[["srce", "trgt"]].apply(frozenset, axis=1)

ufaces = set(faces.ravel())

com_vs = set(all_edges.srce)
for face in ufaces:
    com_vs = com_vs.intersection(all_edges.loc[all_edges["face"] == face, "srce"])
tyssue.topology.base_topology.drop_two_sided_faces(eptm)

Removes all the two (or one?) sided faces from the epithelium

Note that they are not collapsed, but simply eliminated Does not reindex

tyssue.topology.base_topology.get_neighbour_face_pairs(eptm)

Returns a pandas Series of neighboring face pairs (as forzen sets of 2 indexes)

tyssue.topology.base_topology.get_num_common_edges(eptm)

Returns the number of common edges between two neighboring faces this number is set to -1 if those faces are opposite and share the same edges.

tyssue.topology.base_topology.merge_border_edges(sheet, drop_two_sided=True)

Merge edges at the border of a sheet such that no vertex has only one incoming and one outgoing edge.

tyssue.topology.base_topology.merge_vertices(sheet, vert0, vert1, reindex=True)

Merge the two vertices vert0 and vert1 iff they are linked by an edge

If reindex is True (the default), resets indexes and topology data

It is more efficient to call directly collapse_edge

tyssue.topology.base_topology.remove_face(sheet, face)

Removes a face from the mesh.

tyssue.topology.base_topology.split_vert(sheet, vert, face, to_rewire, epsilon, recenter=False)

Creates a new vertex and moves it towards the center of face.

The edges in to_rewire will be connected to the new vertex.

sheet : a tyssue.Sheet instance vert : int, the index of the vertex to split face : int, the index of the face where to move the vertex to_rewire : pd.DataFrame a subset of sheet.edge_df

where all the edges pointing to (or from) the old vertex will point to (or from) the new.

This will leave opened faces and cells

Sheet

tyssue.topology.sheet_topology.cell_division(sheet, mother, geom, angle=None)

Causes a cell to divide

sheet : a ‘Sheet’ instance mother : face index of target dividing cell geom : a 2D geometry angle : division angle for newly formed edge

daughter: face index of new cell

  • Function checks for perodic boundaries if there are, it checks if dividing cell rests on an edge of the periodic boundaries if so, it displaces the boundaries by a half a period and moves the target cell in the bulk of the tissue. It then performs cell division normally and reverts the periodic boundaries to the original configuration

tyssue.topology.sheet_topology.face_division(sheet, mother, vert_a, vert_b)

Divides the face associated with edges indexed by edge_a and edge_b, splitting it in the middle of those edes.

tyssue.topology.sheet_topology.get_division_edges(sheet, mother, geom, angle=None, axis='x')
tyssue.topology.sheet_topology.resolve_t1s(sheet, geom, model, solver, max_iter=60)
tyssue.topology.sheet_topology.split_vert(sheet, vert, face=None, multiplier=1.5, reindex=True, recenter=False, epsilon=None)

Splits a vertex towards the center of the face.

This operation removes the face face from the neighborhood of the vertex.

tyssue.topology.sheet_topology.type1_transition(sheet, edge01, *, remove_tri_faces=True, multiplier=1.5)

Performs a type 1 transition around the edge edge01

See ../../doc/illus/t1_transition.png for a sketch of the definition of the vertices and cells letterings See Finegan et al. for a description of the algotithm https://doi.org/10.1101/704932

sheet : a Sheet instance edge_01 : int

index of the edge around which the transition takes place

epsilonfloat, optional, deprecated

default 0.1, the initial length of the new edge, in case “threshold_length” is not in the sheet.settings

remove_tri_facesbool, optional

if True (the default), will remove triangular cells after the T1 transition is performed

multiplierfloat, optional

default 1.5, the multiplier to the threshold length, so that the length of the new edge is set to multiplier * threshold_length

Bulk and Monolayer

tyssue.topology.bulk_topology.HI_transition(eptm, face)

H → I transition as defined in Okuda et al. 2013 (DOI 10.1007/s10237-012-0430-7). See tyssue/doc/illus/IH_transition.png for the algorithm

tyssue.topology.bulk_topology.IH_transition(eptm, edge)

I → H transition as defined in Okuda et al. 2013 (DOI 10.1007/s10237-012-0430-7). See tyssue/doc/illus/IH_transition.png for the algorithm

tyssue.topology.bulk_topology.cell_division(eptm, mother, geom, vertices=None, mother_verts=None, daughter_verts=None)
tyssue.topology.bulk_topology.close_cell(eptm, cell)

Closes the cell by adding a face. Assumes a single face is missing

tyssue.topology.bulk_topology.find_HIs(eptm, shorts=None)
tyssue.topology.bulk_topology.find_IHs(eptm, shorts=None)
tyssue.topology.bulk_topology.find_rearangements(eptm)

Finds the candidates for IH and HI transitions Returns ——- edges_HI: set of indexes of short edges faces_IH: set of indexes of small triangular faces

tyssue.topology.bulk_topology.fix_pinch(eptm)

Due to rearangements, some faces in an epithelium will have more than one opposite face.

This method fixes the issue so we can have a valid epithelium back.

tyssue.topology.bulk_topology.get_division_edges(eptm, mother, plane_normal, plane_center=None, return_verts=False)

Returns an index of the mother cell edges crossed by the division plane, ordered clockwize around the division plane normal.

tyssue.topology.bulk_topology.get_division_vertices(eptm, division_edges=None, mother=None, plane_normal=None, plane_center=None, return_all=False)
tyssue.topology.bulk_topology.remove_cell(eptm, cell)

Removes a tetrahedral cell from the epithelium.

tyssue.topology.bulk_topology.split_vert(eptm, vert, face=None, multiplier=1.5)

Splits a vertex towards a face.

eptm : a tyssue.Epithelium instance vert : int the vertex to split face : int, optional, the face to split

if face is None, one face will be chosen at random

multiplier: float, default 1.5

length of the new edge(s) in units of eptm.settings[“threshold_length”]

For a given face, we look for the adjacent cell with the lowest number of faces converging on the vertex. If this number is higher than 4 we raise a ValueError

If it’s 3, we do a OI transition, resulting in a new edge but no new faces If it’s 4, we do a IH transition, resulting in a new face and 2 ne edges.

see ../doc/illus/IH_transition.png

tyssue.topology.monolayer_topology.cell_division(monolayer, mother, orientation='vertical', psi=None)

Divides the cell mother in the monolayer.

  • monolayer: a Monolayer instance

  • mother: int, the index of the cell to devide

  • orientation: str, {“vertical” | “horizontal” | “apical”} if “horizontal”, performs a division in the equatorial plane of the cell. If “vertical” (the default), performs a division along the basal-apical axis of the cell. If “apical”, performs a division cutting the apical face perpendicularly to its principal axis

  • psi: float, default 0 extra rotation angle of the division plane around the basal-apical plane

  • daughter: int, the index of the daughter cell

tyssue.topology.monolayer_topology.find_basal_edge(monolayer, apical_edge)

Returns the basal edge parallel to the apical edge passed in argument.

monolayer: a Monolayer instance

Dynamic models definitions

Gradients

tyssue.dynamics.base_gradients.length_grad(sheet)

returns -(dx/l, dy/l, dz/l), ie grad_i(l_ij))

tyssue.dynamics.planar_gradients.area_grad(sheet)
tyssue.dynamics.planar_gradients.lumen_area_grad(eptm)

Base gradients for sheet like geometries

tyssue.dynamics.sheet_gradients.area_grad(sheet)
tyssue.dynamics.sheet_gradients.height_grad(sheet)

Base gradients for bulk geometry

tyssue.dynamics.bulk_gradients.lumen_volume_grad(eptm)

Calculates the gradient for the volume enclosed by the epithelium.

For a monolayer, it will by default compute the volume enclosed by the basal side (edges whose ‘segment’ column is “basal”). If the polarity is reversed and the apical side faces the lumen, this can be changed by setting eptm.settings[“lumen_side”] to ‘apical’

tyssue.dynamics.bulk_gradients.volume_grad(eptm)

Effectors and Model factory

tyssue.dynamics.factory.model_factory(effectors, ref_effector=None)

Produces a Model class with the provided effectors.

effectors : list of effectors.AbstractEffectors classes. ref_effector : optional, default None

if passed, will be used for normalization, by default, the last effector in the list is used

NewModela Model derived class with compute_enregy and compute_gradient

methods

Generic forces and energies

class tyssue.dynamics.effectors.AbstractEffector

The effector class is used by model factories to construct a model.

dimensions = None
element = None
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Abstract effector'
magnitude = None
spatial_ref = (None, None)
specs = {'cell': {}, 'edge': {}, 'face': {}, 'vert': {}}
temporal_ref = (None, None)
class tyssue.dynamics.effectors.BarrierElasticity

Barrier use to maintain the tissue integrity.

dimensions = array(1.) * fJ/um**2
element = 'vert'
static energy(eptm)
static gradient(eptm)
label = 'Barrier elasticity'
magnitude = 'barrier_elasticity'
specs = {'vert': {'barrier_elasticity': 1.0, 'delta_rho': 0.0, 'is_active': 1}}
class tyssue.dynamics.effectors.BorderElasticity
dimensions = array(1.) * fJ/um**2
element = 'edge'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Border edges elasticity'
magnitude = 'border_elasticity'
spatial_ref = ('prefered_length', UnitLength('micrometer', 0.001 * mm, 'um'))
specs = {'edge': {'border_elasticity': 1.0, 'is_active': 1, 'is_border': 1.0, 'length': 1.0, 'prefered_length': 1.0}}
class tyssue.dynamics.effectors.CellAreaElasticity
dimensions = array(1.) * fJ/um**4
element = 'cell'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Area elasticity'
magnitude = 'area_elasticity'
spatial_ref = ('prefered_area', array(1.) * um**2)
specs = {'cell': {'area': 1.0, 'area_elasticity': 1.0, 'is_alive': 1, 'prefered_area': 1.0}}
class tyssue.dynamics.effectors.CellVolumeElasticity
dimensions = array(1.) * fJ/um**6
element = 'cell'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Volume elasticity'
magnitude = 'vol_elasticity'
spatial_ref = ('prefered_vol', array(1.) * um**3)
specs = {'cell': {'is_alive': 1, 'prefered_vol': 1.0, 'vol': 1.0, 'vol_elasticity': 1.0}}
class tyssue.dynamics.effectors.FaceAreaElasticity
dimensionless = False
dimensions = array(1.) * fJ/um**4
element = 'face'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Area elasticity'
magnitude = 'area_elasticity'
spatial_ref = ('prefered_area', array(1.) * um**2)
specs = {'edge': {'sub_area': 0.16666666666666666}, 'face': {'area': 1.0, 'area_elasticity': 1.0, 'is_alive': 1, 'prefered_area': 1.0}}
class tyssue.dynamics.effectors.FaceContractility
dimensions = array(1.) * fJ/um**2
element = 'face'
static energy(eptm)
static gradient(eptm)
label = 'Contractility'
magnitude = 'contractility'
spatial_ref = ('mean_perimeter', UnitLength('micrometer', 0.001 * mm, 'um'))
specs = {'face': {'contractility': 1.0, 'is_alive': 1, 'perimeter': 1.0}}
class tyssue.dynamics.effectors.FaceVolumeElasticity
dimensions = array(1.) * fJ/um**6
element = 'face'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Volume elasticity'
magnitude = 'vol_elasticity'
spatial_ref = ('prefered_vol', array(1.) * um**3)
specs = {'edge': {'sub_area': 0.16666666666666666}, 'face': {'is_alive': 1, 'prefered_vol': 1.0, 'vol': 1.0, 'vol_elasticity': 1.0}, 'vert': {'height': 1.0}}
class tyssue.dynamics.effectors.LengthElasticity

Elastic half edge

dimensions = array(1.) * fJ/um**2
element = 'edge'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Length elasticity'
magnitude = 'length_elasticity'
spatial_ref = ('prefered_length', UnitLength('micrometer', 0.001 * mm, 'um'))
specs = {'edge': {'is_active': 1, 'length': 1.0, 'length_elasticity': 1.0, 'prefered_length': 1.0, 'ux': 0.5773502691896257, 'uy': 0.5773502691896257, 'uz': 0.5773502691896257}}
class tyssue.dynamics.effectors.LineTension
dimensions = array(1.) * fJ/um
element = 'edge'
static energy(eptm)
static gradient(eptm)
label = 'Line tension'
magnitude = 'line_tension'
spatial_ref = ('mean_length', UnitLength('micrometer', 0.001 * mm, 'um'))
specs = {'edge': {'is_active': 1, 'line_tension': 1.0}}
class tyssue.dynamics.effectors.LineViscosity
dimensions = array(1.) * s*nN/um
element = 'edge'
static gradient(eptm)
label = 'Linear viscosity'
magnitude = 'edge_viscosity'
spatial_ref = ('mean_length', UnitLength('micrometer', 0.001 * mm, 'um'))
specs = {'edge': {'edge_viscosity': 1.0, 'is_active': 1}}
temporal_ref = ('dt', UnitTime('second', 's'))
class tyssue.dynamics.effectors.LumenAreaElasticity

..math:

rac{K_Y}{2}(A_{mathrm{lumen}} - A_{0,mathrm{lumen}})^2

dimensions = array(1.) * fJ/um**4
element = 'settings'
static energy(eptm)
static gradient(eptm)
label = 'Lumen volume constraint'
magnitude = 'lumen_elasticity'
spatial_ref = ('lumen_prefered_vol', array(1.) * um**2)
specs = {'settings': {'lumen_elasticity': 1.0, 'lumen_prefered_vol': 1.0, 'lumen_vol': 1.0}}
class tyssue.dynamics.effectors.LumenVolumeElasticity

Global volume elasticity of the object. For example the volume of the yolk in the Drosophila embryo

dimensions = array(1.) * fJ/um**6
element = 'settings'
static energy(eptm)
static get_nrj_norm(specs)
static gradient(eptm)
label = 'Lumen volume elasticity'
magnitude = 'lumen_vol_elasticity'
spatial_ref = ('lumen_prefered_vol', array(1.) * um**3)
specs = {'settings': {'lumen_prefered_vol': 1.0, 'lumen_vol': 1.0, 'lumen_vol_elasticity': 1.0}}
class tyssue.dynamics.effectors.PerimeterElasticity

From Mapeng Bi et al. https://doi.org/10.1038/nphys3471

dimensions = array(1.) * fJ/um**2
element = 'face'
static energy(eptm)
static gradient(eptm)
label = 'Perimeter Elasticity'
magnitude = 'perimeter_elasticity'
spatial_ref = ('prefered_perimeter', UnitLength('micrometer', 0.001 * mm, 'um'))
specs = {'face': {'is_alive': 1, 'perimeter': 1.0, 'perimeter_elasticity': 0.1, 'prefered_perimeter': 3.81}}
class tyssue.dynamics.effectors.RadialTension

Apply a tension perpendicular to a face.

dimensions = array(1.) * fJ/um
element = 'face'
static energy(eptm)
static gradient(eptm)
label = 'Apical basal tension'
magnitude = 'radial_tension'
specs = {'face': {'height': 1.0, 'radial_tension': 1.0}}
class tyssue.dynamics.effectors.SurfaceTension
dimensions = array(1.) * fJ/um**2
element = 'face'
static energy(eptm)
static gradient(eptm)
label = 'Surface tension'
magnitude = 'surface_tension'
spatial_ref = ('prefered_area', array(1.) * um**2)
specs = {'face': {'area': 1.0, 'is_active': 1, 'surface_tension': 1.0}}
tyssue.dynamics.effectors.dimensionalize(nondim_specs, dim_specs, effector, ref_effector)
tyssue.dynamics.effectors.elastic_energy(element_df, var, elasticity, prefered)
tyssue.dynamics.effectors.elastic_force(element_df, var, elasticity, prefered)
tyssue.dynamics.effectors.normalize(dim_specs, nondim_specs, effector, ref_effector)
tyssue.dynamics.effectors.scaler(nondim_specs, dim_specs, effector, ref_effector)

Predefined models

Vertex model for an Epithelial sheet (see definitions).

Depends on the sheet vertex geometry functions.

Specific functions for apoptosis vertex model

class tyssue.dynamics.apoptosis_model.ApicoBasalTension

Effector for the apical-basal tension.

The energy is proportional to the heigth of the cell

dimensions = array(1.) * fJ/um**2
element = 'vert'
static energy(sheet)
static gradient(sheet)
label = 'Apical-basal tension'
magnitude = 'radial_tension'
specs = {'vert': {'height': 1.0, 'is_active': 1, 'radial_tension': 0.0}}

Dynamical models for monlayer and bulk epithelium.

class tyssue.dynamics.bulk_model.LaminaModel

Not implemented yet

Quasistatic solver

Epithelium generation utilities

The generation module provides utilities to easily create Epithelium objects.

class tyssue.generation.shapes.AnnularSheet(identifier, datasets, specs=None, coords=None)

2D annular model of a cylinder-like monolayer.

Provides syntactic sugar to access the apical, basal and lateral segments of the epithlium

property apical_edges
property apical_verts
property basal_edges
property basal_verts
property lateral_edges
reset_topo()

Recomputes the number of sides for the faces and the number of faces for the cells.

segment_index(segment, element)
tyssue.generation.shapes.Lloyd_relaxation(sheet, geom, steps=10, coords=None, update_method=None)

Performs Lloyd relaxation on the sheet.

tyssue.generation.shapes.ellipse_rho(theta, a, b)
tyssue.generation.shapes.ellipsoid_sheet(a, b, c, n_zs, **kwargs)

Creates an ellipsoidal apical mesh.

a, b, cfloats

Size of the ellipsoid half axes in the x, y, and z directions, respectively

n_zsint

The (approximate) number of faces along the z axis.

kwargs are passed to get_ellipsoid_centers

eptm : a Epithelium object

The mesh returned is an Epithelium and not a simpler Sheet so that a unique cell data can hold information on the whole volume of the ellipsoid.

tyssue.generation.shapes.generate_ring(Nf, R_in, R_out, R_vit=None, apical='in')

Generates a 2D tyssue object aranged in a ring of Nf tetragonal cells with inner diameter R_in and outer diameter R_out

Nfint

The number of cells in the tissue

R_infloat

The inner ring diameter

R_outfloat

The outer ring diameter

R_vitfloat

The vitelline membrane diameter (a non strechable membrane around the annulus)

apicalstr {‘in’ | ‘out’}

The side of the apical surface if “in”, the apical surface is inside the annulus, facing the lumen as in an organoid; if ‘out’: the apical side is facing the exterior of the tissue, as in an embryo

eptmAnnularSheet

2D annular tissue. The R_in and R_out parameters are stored in the class settings attribute.

tyssue.generation.shapes.get_ellipsoid_centers(a, b, c, n_zs, pos_err=0.0, phase_err=0.0)

Creates hexagonaly organized points on the surface of an ellipsoid

a, b, c: float

ellipsoid radii along the x, y and z axes, respectively i.e the ellipsoid boounding box will be [[-a, a], [-b, b], [-c, c]]

n_zsfloat

number of cells on the z axis, typical

pos_errfloat, default 0.

normaly distributed noise of std. dev. pos_err is added to the centers positions

phase_errfloat, default 0.

normaly distributed noise of std. dev. phase_err is added to the centers angle ϕ

tyssue.generation.shapes.sheet_from_cell_centers(points, noise=0, interp_s=0.0001)

Returns a Sheet object from the Voronoï tessalation of the cell centers.

The strategy is to project the points on a sphere, get the Voronoï tessalation on this sphere and reproject the vertices on the original (implicit) surface through linear interpolation of the cell centers.

Works for relatively smooth surfaces (at the very minimum star convex).

pointsnp.ndarray of shape (Nf, 3)

the x, y, z coordinates of the cell centers

noisefloat, default 0.0

addiditve normal noise stdev

interp_sfloat, default 1e-4

interpolation smoothing factor (might need to set higher)

sheet : a Sheet object with Nf faces

tyssue.generation.shapes.spherical_monolayer(R_in, R_out, Nc, apical='out', Lloyd_relax=False)

Returns a spherical monolayer with the given inner and outer radii, and approximately the gieven number of cells.

The apical argument can be ‘in’ out ‘out’ to specify wether the apical face of the cells faces inward or outward, reespectively.

tyssue.generation.shapes.spherical_sheet(radius, Nf, Lloyd_relax=False, **kwargs)

Returns a spherical sheet with the given radius and (approximately) the given number of cells

tyssue.generation.shapes.update_on_sphere(sheet)

This module provides utlities to modify an input tissue through extrusion or subdivision

tyssue.generation.modifiers.create_anchors(sheet)

Adds an edge linked to every vertices at the boundary and create anchor vertices

tyssue.generation.modifiers.extrude(apical_datasets, method='homotecy', scale=0.3, vector=[0, 0, - 1])

Extrude a sheet to form a monlayer epithelium

  • apical_datasets: dictionnary of three DataFrames,

‘vert’, ‘edge’, ‘face’ * method: str, optional {‘homotecy’|’translation’|’normals’} default ‘homotecy’ * scale: float, optional the scale factor for homotetic scaling, default 0.3. * vector: sequence of three floats, optional, used for the translation

default [0, 0, -1]

if method == ‘homotecy’, the basal layer is scaled down from the apical one homoteticaly w/r to the center of the coordinate system, by a factor given by scale

if method == ‘translation’, the basal vertices are translated from the apical ones by the vector vect

if method == ‘normals’, basal vertices are translated from the apical ones along the normal of the surface at each vertex, by a vector whose size is given by scale

tyssue.generation.modifiers.subdivide_faces(eptm, faces)

Adds a vertex at the center of each face, and returns a new dataset

eptm: a Epithelium instance faces: list,

indices of the faces to be subdivided

new_dset: dict

a dataset with the new faces devided

Hexagonal grids

tyssue.generation.hexagonal_grids.circle(num_t, radius=1.0, phase=0.0)

Returns x and y positions of num_t points regularly placed around a circle of radius radius, shifted by phase radians.

num_tint

the number of points around the circle

radiusfloat, default 1.

the radius of the circle

phasefloat, default 0.0

angle shift w/r to the x axis in radians

points : np.Ndarray of shape (num_t, 2), the x, y positions of the points

tyssue.generation.hexagonal_grids.hexa_cylinder(num_t, num_z, radius=1.0, capped=False, noise=0, orientation='transverse')

Returns an arrays of x, y positions of points evenly spread on a cylinder with num_t points on the periphery and num_z points on its length.

num_tint,

The number of points on the periphery

num_zint,

The number of points along the z axis (the length of the cylinder)

radiusfloat, default 1

The radius of the cylinder

cappedbool, default False

If True, the tips of the cylinder are capped by a disk of point as generated by the hexa_disk function.

noisefloat, default 0

normaly distributed position noise around the cell points

orientation{‘transverse’ | ‘longitudinal’}, default ‘transverse’

the orientation of the cells (with the longueur axis perpendicular or along the length of the cylinder)

tyssue.generation.hexagonal_grids.hexa_disk(num_t, radius=1)

Returns an arrays of x, y positions of points evenly spread on a disk with num_t points on the periphery.

num_tint

the number of poitns on the disk periphery, the rest of the disk is filled automaticaly

radiusfloat, default 1.

the radius of the disk

tyssue.generation.hexagonal_grids.hexa_grid2d(nx, ny, distx, disty, noise=None)

Creates an hexagonal grid of points

tyssue.generation.hexagonal_grids.hexa_grid3d(nx, ny, nz, distx=1.0, disty=1.0, distz=1.0, noise=None)

Creates an hexagonal grid of points

tyssue.generation.hexagonal_grids.three_faces_sheet(zaxis=False)

Creates the apical junctions mesh of three packed hexagonal faces. If zaxis is True (defaults to False), adds a z coordinates, with z = 0.

Faces have a side length of 1.0 +/- 1e-3.

face_df: the faces DataFrame indexed from 0 to 2 vert_df: the junction vertices DataFrame edge_df: the junction edges DataFrame

tyssue.generation.hexagonal_grids.three_faces_sheet_array()

Creates the apical junctions mesh of three packed hexagonal faces. If zaxis is True (defaults to False), adds a z coordinates, with z = 0.

Faces have a side length of 1.0 +/- 1e-3.

points: (13, ndim) np.array of floats

the positions, where ndim is 2 or 3 depending on zaxis

edges: (15, 2) np.array of ints

indices of the edges

(Nc, Nv, Ne): triple of ints

number of faces, vertices and edges (3, 13, 15)

tyssue.generation.from_voronoi.from_2d_voronoi(voro, specs=None)

Creates 2D (sheet geometry) datasets from a Voronoï tessalation

voro: a scipy.spatial.Voronoi object

datasets: dict

datasets suitable for Epithelium implementation

tyssue.generation.from_voronoi.from_3d_voronoi(voro)

Creates 3D (bulk geometry) datasets from a Voronoï tessalation

voro: a scipy.spatial.Voronoi object

datasets: dict

datasets suitable for Epithelium implementation

Input/output

tyssue.io.hdf5.load_datasets(h5store, data_names=['face', 'vert', 'edge', 'cell', 'settings'])
tyssue.io.hdf5.save_datasets(h5store, eptm)
tyssue.io.csv.write_storm_csv(filename, points, coords=['x', 'y', 'z'], split_by=None, **csv_args)

Saves a point cloud array in the storm format

tyssue.io.obj.save_junction_mesh(filename, eptm)
tyssue.io.obj.save_splitted_cells(fname, sheet, epsilon=0.1)
tyssue.io.obj.save_triangulated(filename, eptm)
tyssue.io.obj.write_splitted_cells(*args, **kwargs)

Predefined datasets

Available predefined datasets:small_hexagonal_snaped.hf5 small_ellipsoid.hf5 rod_sheet.hf5 sheet6x5.hf5 with_4sided_cell.hf5 planar_periodic8x8.hf5 __pycache__ __init__.py small_hexagonal.hf5 before_apoptosis.hf5 15_cells_patch.hf5 with_pinch.hf5

tyssue.stores.load_datasets(store, **kwargs)

This will soon be deprecated, use the leaner `tyssue.stores.stores_dir and tyssue.stores.stores_list

Collision detection and correction

Detection

tyssue.collisions.intersection.self_intersections(sheet)

Checks for self collisions for the sheet

sheeta Sheet object

This object must have a triangular_mesh method returning a valid triangular mesh.

edge_pairs: np.ndarray of indices

Array of shape (n_intersections, 2) with the indices of the pairs of intersecting edges

Resolution

class tyssue.collisions.solvers.CollidingBoxes(sheet, position_buffer, intersecting_edges)

Utility class to manage collisions

get_limits(shyness=1e-10)

Iterator over the position boundaries avoiding the collisions.

shynessfloat

the extra distance between two colliding vertices, on each side of the collision plane.

lower, uppertwo pd.Series

those Series are indexed by the vertices of the colliding faces giving the lower and upper bounds for the vertices

solve_collisions(shyness=1e-10)

Solves the collisions by finding the collision plane.

Modifies the sheet vertex positions inplace such that they rest at a distance shyness apart on each side of the collision plane.

shynessfloat, default 1e-10

the extra distance between two colliding vertices, on each side of the collision plane.

Based on Liu, J.-D., Ko, M.-T., & Chang, R.-C. (1998), A simple self-collision avoidance for cloth animation. Computers & Graphics, 22(1), 117–128. DOI

tyssue.collisions.solvers.auto_collisions(fun)

Decorator to solve collisions detections after the execution of the decorated function.

It is assumed that the two first arguments of the decorated function are a Sheet object and a geometry class

The function is re-executed with the updated geometry

tyssue.collisions.solvers.revert_positions(sheet, position_buffer, intersecting_edges)
tyssue.collisions.solvers.solve_bulk_collisions(eptm, position_buffer)

Corrects the auto-collisions for the outer surface(s) of a 3D epithelium.

eptm : a Epithelium object position_buffer : np.array of shape (eptm.Nv, eptm.dim):

positions of the vertices prior to the collisions

changedbool

True if the positions of some vertices were changed

tyssue.collisions.solvers.solve_sheet_collisions(sheet, position_buffer)

Corrects the auto-collisions for the outer surface(s) of a 2.5D sheet.

sheet : a Sheet object position_buffer : np.array of shape (sheet.Nv, sheet.dim):

positions of the vertices prior to the collisions

changedbool

True if the positions of some vertices were changed

Miscellanous utils

tyssue.utils.utils.ar_calculation(sheet, coords=['x', 'y'])

Calculates the aspect ratio of each face of the sheet

eptm : a Sheet object coords : list of str, optional, default [‘x’, ‘y’]

the coordinates on which to compute the aspect ratio

AR: pandas series of aspect ratio for all faces.

As is the case in ImageJ, the returned aspect ratio is always higher than 1

tyssue.utils.utils.cell_centered_patch(eptm, cell, neighbour_order)

Return subsheet centered on cell with a distance of neighbour order around the cell

eptm : a Epithelium instance face : int, id of the center face neighbour_order: int, number of neighbour around the center face

patch: an object of the same class as the input object

tyssue.utils.utils.combine_specs(*specs)
tyssue.utils.utils.data_at_opposite(sheet, edge_data, free_value=None)

Returns a pd.DataFrame with the values of the input edge_data at the opposite edges. For free edges, optionaly replaces Nan values with free_value

sheet: a Sheet instance edge_data: dataframe contain value of edge

opposite: pandas series contain value of opposite edge

tyssue.utils.utils.elem_centered_patch(eptm, elem_idx, neighbour_order, elem)

Return subeptm centered on the element (cell or face) with index elem_idx with neighbour_order neighbours around it.

eptm : a Epithelim instance index : int, id of the center element neighbour_order: int,

neighbourhood ‘degree’ around the center element

patch: an object with the same class as eptm

tyssue.utils.utils.face_centered_patch(sheet, face, neighbour_order)

Return subsheet centered on face with a distance of neighbour order around the face

sheet : a Sheet object face : int, id of the center face neighbour_order: int, number of neighbour around the center face

patch: an object of the same class as the input object

tyssue.utils.utils.get_next(eptm)

Returns the indices of the next edge for each edge

tyssue.utils.utils.get_sub_eptm(eptm, edges, copy=False)

Define sub-epithelium corresponding to the edges.

eptm: a Epithelium instance edges: list of edges includes in the sub-epithelium

sub_eptm: a Epithelium instance

tyssue.utils.utils.modify_segments(eptm, modifiers)

Modifies the datasets of a segmented epithelium according to the passed modifiers.

eptm : tyssue.Epithelium modifiers : nested dictionnary

This functions assumes that the epithelium has a segment_index method as implemented in the tyssue.Monolayer.

>>> modifiers = {
>>>     'apical' : {
>>>         'edge': {'line_tension': 1.},
>>>         'face': {'prefered_area': 0.2},
>>>     },
>>>     'basal' : {
>>>         'edge': {'line_tension': 3.},
>>>         'face': {'prefered_area': 0.1},
>>>     }
>>> modify_segments(monolayer, modifiers)
>>> monolayer.ver_df.loc[monolayer.apical_edges,
>>>                      'line_tension'].unique()[0] == 1.
True
tyssue.utils.utils.scaled_unscaled(func, scale, eptm, geom, args=(), kwargs={}, coords=None)

Scales the epithelium by an homotetic factor scale, applies the function func, and scales back to original size.

func: the function to apply to the scaled epithelium scale: float, the scale to apply eptm: a Epithelium instance geom: a Geometry class args: sequence, the arguments to pass to func kwargs: dictionary, the keywords arguments

to pass to func

coords: the coordinates on which the scaling applies

If the execution of function fails, the scaling is still reverted

res: the result of the function func

tyssue.utils.utils.set_data_columns(datasets, specs, reset=False)

Sets the columns of the dataframes in the datasets dictionnary to the uniform values in the specs sub-dictionnaries.

datasets : dict of dataframes specs : dict of dicts reset : bool, default False

For each key in specs, the value is a dictionnary whose keys are column names for the corresponding dataframe in datasets. If there is no such column in the dataframe, it is created. If the columns allready exists and reset is True, the new value is used.

tyssue.utils.utils.single_cell(eptm, cell, copy=False)

Define epithelium instance for all element to a define cell.

eptm : a Epithelium instance cell : identifier of a cell copy : bool, default False

sub_etpm: class:’Epithelium’ instance corresponding to the cell

tyssue.utils.utils.spec_updater(specs, new)

Add element to the new dictionary to the specs dictionary. Update value if the key already exist.

specs: specification that will be modified new: dictionary of new specification

tyssue.utils.utils.swap_apico_basal(organo)

Swap apical and basal segments of an organoid.

tyssue.utils.utils.to_nd(df, ndim)

Give a new shape to an input data by duplicating its column.

df : input data that will be reshape ndim : dimension of the new reshape data.

df_nd : return array reshaped in ndim.

Decorators

tyssue.utils.decorators.cell_lookup(func)
tyssue.utils.decorators.do_undo(func)

Decorator that creates a copy of the first argument (usually an epithelium object) and restores it if the function fails.

The first argument in *args should have backup() and restore() methods.

tyssue.utils.decorators.face_lookup(func)
tyssue.utils.decorators.time_exe(func)
tyssue.utils.decorators.validate(func)

Decorator that validate the epithelium after the decorated function was applied. the first argument of func should be an epithelium instance, and is at least assumed to have a validate method.