14.2. steps.API_2.geom

The geom module contains classes related to meshes and geometries: compartments, patches, region of interest, diffusion boundaries, etc. It also contains utility classes for handling element lists.

14.2.1. Detailed documentation

  • Containers

  • References

  • Reference lists

  • Convenience classes


class Geometry(*args, **kwargs)[source]

Bases: NamedObject

Top-level geometry container

A number of compartment objects and patches objects may be grouped. Should be used as a context manager for the declaration of compartments etc:

geom = Geometry()
with geom:
    comp1 = Compartment.Create(vsys, volume)

    ... # Declare other objects in geom

After having declared children objects in the with mesh: block, they can be accessed as attributes of geom with their name (see steps.API_2.utils.NamedObject and steps.API_2.utils.NamedObject.Create()):

geom.comp1
class Compartment(*args, **kwargs)[source]

Bases: NamedObject

Base class for compartment objects

The same class is used to declare compartments in well-mixed models and in 3D tetrahedral meshes.

Compartments in well-mixed geometry can be declared in the following ways:

with geom:
    wmComp1 = Compartment.Create()
    wmComp2 = Compartment.Create(vsys)
    wmComp3 = Compartment.Create(vol=volume)
    wmComp4 = Compartment.Create(vsys, volume)
Parameters

Compartments in a tetrahedral mesh can be declared in the following ways:

with mesh:
    tetComp1 = Compartment.Create(tetLst)
    tetComp2 = Compartment.Create(tetLst, vsys)
Parameters
  • tetLst (TetList or any argument that can be used to build one) – List of tetrahedrons associated with the compartment.

  • vsys (steps.API_2.model.VolumeSystem or str) – Optional, the volume system associated with this compartment.

The volume of a compartment in a tetrahedral mesh is the total volume of the encapsulated tetrahedrons.

Compartments in a distributed tetrahedral mesh can be declared in the following ways:

with mesh:
    distComp1 = Compartment.Create()
    distComp2 = Compartment.Create(vsys, conductivity=g)
    distComp3 = Compartment.Create(physicalTag=10)
    distComp4 = Compartment.Create(tetLst, vsys)
Parameters
  • vsys (steps.API_2.model.VolumeSystem or str) – Optional, the volume system associated with this compartment.

  • conductivity (float) – Optional, the volume conductivity of the compartment.

  • physicalTag (int) – Optional, the index of the physical tag used to build the compartment (see DistMesh).

  • tetLst (TetList or any argument that can be used to build one) – Optional, List of tetrahedrons associated with the compartment. If it is ommited or an empty list is given, the compartment name will be used to find tetrahedrons tagged with it.

Some of the methods documented below are only available for tetrahedral compartments, they are grouped accordingly.

addSystem(vsys)[source]

Add a volume system to the compartment

Parameters

vsys (steps.API_2.model.VolumeSystem or str) – The volume system to be added, or its name.

Returns

None

property Vol

Volume of the compartment

Type

Union[float, steps.API_2.utils.Parameter], read-only for tetrahedral compartments

Only available for compartments defined in tetrahedral meshes:
property surface

All triangles in the surface of the compartment

Type

TriList, read-only

property tets

All tetrahedrons in the compartment

Type

TetList, read-only

property bbox

The bounding box of the compartment

Type

steps.API_2.geom.BoundingBox, read-only

Only available for compartments defined in distributed tetrahedral meshes:
property Conductivity

Conductivity of the compartment

Type

Union[float, steps.API_2.utils.Parameter]

class Patch(*args, **kwargs)[source]

Bases: NamedObject

Base class for patch objects

A patch is a piece of 2D surface surrounding (part of) a 3D compartment, which may be connected to another compartment. The same class is used to declare patches in well-mixed models and in 3D tetrahedral meshes.

Patches in well-mixed geometry can be declared in the following ways:

with geom:
    patch1 = Patch.Create(inner)
    patch2 = Patch.Create(inner, outer)
    patch3 = Patch.Create(inner, outer, ssys)
    patch4 = Patch.Create(inner, None , ssys)
    patch5 = Patch.Create(inner, None , None, area)
    patch6 = Patch.Create(inner, outer, None, area)
    patch7 = Patch.Create(inner, outer, ssys, area)
Parameters

Patches in a tetrahedral mesh can be declared in the following ways:

with mesh:
    patch1 = Patch.Create(triLst, inner)
    patch2 = Patch.Create(triLst, inner, outer)
    patch3 = Patch.Create(triLst, inner, None , ssys)
    patch4 = Patch.Create(triLst, inner, outer, ssys)
Parameters

Patches in a distributed tetrahedral mesh can be declared in the following ways:

with mesh:
    comp1 = Patch.Create(inner)
    comp2 = Patch.Create(inner, ssys=ssys)
    comp3 = Patch.Create(inner, outer, ssys)
    comp4 = Patch.Create(inner, physicalTag=10)
    comp5 = Patch.Create(triList, inner, ssys=ssys)
Parameters
  • inner (steps.API_2.geom.Compartment) – Inner compartment, see below

  • outer (steps.API_2.geom.Compartment) – Optional, Outer compartment, see below

  • ssys (steps.API_2.model.SurfaceSystem or str) – Optional, the surface system associated with this patch.

  • physicalTag (int) – Optional, the index of the physical tag used to build the patch (see DistMesh).

  • triLst (TriList or any argument that can be used to build one) – Optional, list of triangles associated with the patch. If an empy list is given the patch will be created from triangles that are tagged with the patch name.

Relationship between Compartments and Patches

It is necessary to explain the inner/outer relationship between compartments and patches. When a patch object is created. it is necessary to arbitrarily label the compartment(s) “inner” and “outer” (if a patch is connected to only one compartment then the compartment must be labelled “inner” by convention). This is necessary in order to fully describe the surface reaction rules. Accordingly, compartments also store a list of connections, “inner” patches and “outer” patches. So if a patch1 is created with comp1 as it’s “inner” compartment, comp1 knows patch1 as an “outer” patch. The labelling is purely defined when creating the Patch objects, bearing in mind the stoichiometry defined in the surface reaction objects. This may seem a little confusing at first, but will become clearer when experience is gained with these objects.

For patches in tetrahedral meshes, the area is the total area of the encapsulated triangles.

Some of the methods documented below are only available for tetrahedral patches, they are grouped accordingly.

addSystem(ssys)[source]

Add a surface system to the patch

Parameters

ssys (steps.API_2.model.SurfaceSystem or str) – The surface system to be added, or its name.

Returns

None

property Area

Surface area of the patch

Type

Union[float, steps.API_2.utils.Parameter], read-only for tetrahedral patches

property innerComp

Inner compartment

Type

Compartment

property outerComp

Outer compartment

Type

Compartment

Only available for patches defined in tetrahedral meshes:
property bbox

The bounding box of the patch

Type

steps.API_2.geom.BoundingBox, read-only

property tris

All triangles in the patch

Type

TriList, read-only

property edges

Perimeter of the triangles in the patch

I.e. all bars that are not shared by two triangles in the patch.

Type

BarList, read-only

class Membrane(*args, **kwargs)[source]

Bases: NamedObject

A set of patches on which membrane potential will be computed

This class provides annotation for a group of triangles that comprise a surface describing a membrane in a Tetmesh. This may be the same description of one or several TmPatches in order for voltage-dependent transitions, currents and so on to be inserted in the membrane. A Membrane object must be available if membrane potential calculation is to be performed.

By default STEPS only adds the inner compartments of the patches to the conduction volume, if other compartments should also be added to the conduction volume, they should be supplied through the supplementary_comps keyword parameter.

Parameters

patches (list) – List of patches whose triangles will be added to the membrane (only one patch for distributed meshes)

The following parameters are only available for non-distributed meshes:

Parameters
  • verify (bool) – Perform some checks on suitability of membrane (see below)

  • opt_method (int) – Optimization method (see below)

  • search_percent (float) – See below

  • opt_file_name (str) – Full path to a membrane optimization file

  • supplementary_comps (List[Compartment]) – Supplementary compartments to be considered part of the conduction volume

Perform some checks on suitability of membrane if verify is True: these checks will print warnings if membrane forms an open surface or if any triangle is found to have more than 3 neighbours. Specify optimization method with opt_method:

  1. principal axis ordering (quick to set up but usually results in slower simulation than method 2).

  2. breadth first search (can be time-consuming to set up, but usually faster simulation.

If breadth first search is chosen then argument search_percent can specify the number of starting points to search for the lowest bandwidth. If a filename (with full path) is given in optional argument opt_file_name the membrane optimization will be loaded from file, which was saved previously for this membrane with steps.API_2.sim.Simulation.saveMembOpt()

property open

Whether the membrane is open (contains holes) or not

Type

bool, read-only

property tris

List of triangles associated with the membrane

Type

steps.API_2.geom.TriList, read-only

property Area

Surface area of the membrane

Type

float, read-only

Only available for membranes defined in distributed tetrahedral meshes:
property Capacitance

Capacitance of the compartment

Type

Union[float, steps.API_2.utils.Parameter]

class EndocyticZone(tris, *args, **kwargs)[source]

Bases: NamedObject

A set of triangles that model a zone in which endocytosis reactions can happen

Endocytosis reactions declared in surface systems do not happen at any point of the corresponding patches, they only happen in endocytic zones. By default all endocytosis events are active in this zone, but they can be deactivated during simulation.

Parameters

tris (TriList) – The list of triangle

Endocytic zones are declared by using a patch as a context manager:

with patch1:

zone1 = EndocyticZone.Create(triLst)

Note that the triangles in triLst all need to be part of patch1.

property tris

All triangles in the endocytic zone

Type

TriList, read-only

class ROI(lst=None, *args, **kwargs)[source]

Bases: NamedObject

Region of Interest class

Parameters

lst (TetList, or TriList, or VertList) – Element list

Should be declared using the context manager syntax with TetMesh():

with mesh:
    tetROI1 = ROI.Create(mesh.tets[10:30])

The ROI is then treated as a physical location, like patches or compartment, and can be accessed from the mesh with its name:

mesh.tetROI1
property Vol

Total volume of the ROI, if applicable

Type

float, read-only

property Area

Total area of the ROI, if applicable

Type

float, read-only

property tets

All tetrahedrons in the ROI

Raises a TypeError if the ROI was not created from a TetList.

Type

TetList, read-only

property tris

All triangles in the ROI

Raises a TypeError if the ROI was not created from a TriList.

Type

TriList, read-only

property verts

All vertices in the ROI

Raises a TypeError if the ROI was not created from a VertList.

Type

VertList, read-only

class TetMesh(*args, **kwargs)[source]

Bases: Geometry

Container class for static tetrahedral meshes

This class stores the vertices points, 3D tetrahedral and 2D triangular elements that comprise the mesh. The indices of the elements will be stored as unsigned integers (a positive integer or zero) beginning at zero and incremented by 1. For example, if there are ntets number of tetrahedrons in the mesh, the indices of the tetrahedrons will be [0,1,2,…, (ntets-1)].

A TetMesh object should not be created from scratch but should instead be created with one of the dedicated class methods:

TetMesh objects can be used as a context manager in the same way as Geometry:

mesh = TetMesh.Load('/path/to/file')
with mesh:
    comp1 = Comp.Create(vsys)

    ... # Declare other objects in mesh
classmethod FromData(verts, tets, tris=[], name=None)[source]

Create a mesh from geometrical data

Construct a Tetmesh container with the folowing method: Supply a list of all vertices verts (by Cartesian coordinates), supply a list of all tetrahedrons tets (by indices of the 4 vertices) and supply a full or partial list of triangles tris (by indices of the 3 vertices). Indexing in STEPS begins at zero, so the first 3 coordinates in verts will describe the zeroth vertex, t he next 3 coordinates will describe the 1st vertex and so on. Labelling of the vertices in tets and tris should follow this indexing. Lists must be one-dimensional. Length of verts = nverts*3 where nverts is the total number of vertices; length of tets = ntets*4 where ntets is the total number of tetrahedrons; maximum length of tris ntris*3 where ntris is the total number of triangles. For example, if we have just three tetrahedrons; tet0=[0,1,2,3], tet1=[0,1,3,4] and tet2=[1,3,4,5] then the required one-dimensional list tets=[0,1,2,3,0,1,3,4,1,3,4,5].

Parameters
  • verts (List[int]) – Vertex data

  • tets (List[int]) – Tetrahedron data

  • tris (List[int]) – Optional triangle data

  • name (str) – Optional name for the mesh

Returns

The constructed TetMesh

Return type

TetMesh

classmethod Load(path, scale=1, strict=False, name=None)[source]

Load a mesh in STEPS

The mesh is loaded from an XML file that was previously generated by the TetMesh.Save() method.

Parameters
  • path (str) – The root of the path where the file(s) are stored. e.g. with ‘meshes/spine1’ this function will look for the file /meshes/spine1.xml

  • scale (float) – Optionally rescale the mesh on loading by given factor

  • strict (bool) – Apply strict(-er) checking to the input XML

  • name (str) – Optional name for the mesh

Returns

The loaded TetMesh

Return type

TetMesh

classmethod LoadAbaqus(filename, scale=1, ebs=None, shadow_mesh=None, name=None)[source]

Load a mesh from an ABAQUS-formated mesh file

If blocks or groups of elements are present in the file, they will also be loaded.

Parameters
  • filename (str or Tuple[str, str]) – The Abaqus filename (or path) including any suffix. Can also take a 2-tuple containing separate paths for tetrahedron and triangle data (in this order).

  • scale (float) – Optionally rescale the mesh on loading by given factor

  • ebs (List[str]) – Names of selected element blocks which are included in the mesh (does not apply if a 2-tuple is given as filename)

  • name (str) – Optional name for the mesh

Returns

The loaded TetMesh

Return type

TetMesh

classmethod LoadGmsh(filename, scale=1, name=None)[source]

Load a mesh from a Gmsh (2.2 ASCII)-formated mesh file

If blocks or groups of elements are present in the file, they will also be loaded.

Parameters
  • filename (str) – The Gmsh filename (or path) including any suffix.

  • scale (float) – Optionally rescale the mesh on loading by given factor

  • name (str) – Optional name for the mesh

Returns

The loaded TetMesh

Return type

TetMesh

classmethod LoadVTK(filename, scale=1, name=None)[source]

Load a mesh from a VTK (legacy ASCII)-formated mesh file

If blocks or groups of elements are present in the file, they will also be loaded.

Parameters
  • filename (str) – The VTK filename (or path) including any suffix.

  • scale (float) – Optionally rescale the mesh on loading by given factor

  • name (str) – Optional name for the mesh

Returns

The loaded TetMesh

Return type

TetMesh

classmethod LoadTetGen(pathroot, scale, name=None)[source]

Load a mesh from a TetGen-formated set of files

If blocks or groups of elements are present, they will also be loaded.

Parameters
  • pathroot (str) – The root of the path to the mesh files. E.g. mesh/torus should be given to read files mesh/torus.node, mesh/torus.ele, and optionally mesh/torus.face

  • scale (float) – Optionally rescale the mesh on loading by given factor

  • name (str) – Optional name for the mesh

Returns

The loaded TetMesh

Return type

TetMesh

The TetGen file format for describing a mesh actually consists of multiple files with varying suffixes. Currently, this class method only reads meshes consisting of three files:

  • <input>.node: describing the tetrahedral mesh node points.

  • <input>.ele: describing tetrahedral elements, each of which

    consists of 4 pointers into the <input>.node file. (TetGen also supports 10-node elements; these 6 extra nodes are obviously not read by STEPS.)

  • <input>.face: describing triangular faces, each of which

    consists of 3 pointers into the <input>.node file. This file is optional.

Other files are .vol (list of maximum volumes), .var (variant constraints) .neigh (list of neighbours), .smesh (simple PLC descriptions) and .edge (list of boundary edges) files. None of these seem relevant for our use cases, so we don’t load them even when they are there. In particular, the .neigh file is computed by STEPS itself.

Please refer to the TetGen manual (pages 31-40 in the last edition) for more information on these file formats

tetgen.berlios.de/files/tetgen-manual.pdf

Save(path)[source]

Save the mesh to an XML file

Parameters

path (str) – The root of the path to store the files. e.g. ‘meshes/spine1’ will save data in /meshes/spine1.xml

This file stores the basic information about the mesh which tends to be common information for any software that supports tetrahedral meshes.

  • NODES are stored by cartesian coordinates.

  • TRIANGLES are stored by the indices of their 3 nodes.

  • TETRAHEDRONS are stored by the indices of their 4 nodes.

The XML file also stores information about any compartments or patches created in STEPS (class steps.API_2.geom.Comp steps.API_2.geom.Patch respectively).

  • COMPARTMENT(S) are stored by:

    • their string identification.

    • a list of any volume systems added to the compartment at time of saving.

    • a list of tetrahedrons belonging to the compartment

  • PATCH(ES) are stored by:

    • their string identification.

    • a list of any surface systems added to the patch at time of saving.

    • the inner compartment id.

    • the outer compartment id (if it exists).

    • a list of triangles belonging to this patch.

ConvertToMetis(path)[source]

Convert the mesh data to metis connectivity data

Parameters

path (str) – The path to the file in which to save the Metis data.

property vertGroups

Groups of vertices defined in the loaded mesh

Type

Dict[str, VertList], read-only

Note

Both blocks and groups are accessed through the same property.

property triGroups

Groups of triangles defined in the loaded mesh

Type

Dict[str, TriList], read-only

Note

Both blocks and groups are accessed through the same property.

property tetGroups

Groups of tetrahedrons defined in the loaded mesh

Type

Dict[str, TetList], read-only

Note

Both blocks and groups are accessed through the same property.

property vertProxy

Element proxy object for vertices (see steps.API_1.utilities.meshio.ElementProxy)

Type

steps.API_1.utilities.meshio.ElementProxy

Note

Only available with some import methods.

property triProxy

Element proxy object for triangles (see steps.API_1.utilities.meshio.ElementProxy)

Type

steps.API_1.utilities.meshio.ElementProxy

Note

Only available with some import methods.

property tetProxy

Element proxy object for tetrahedrons (see steps.API_1.utilities.meshio.ElementProxy)

Type

steps.API_1.utilities.meshio.ElementProxy

Note

Only available with some import methods.

property Vol

Total volume of the mesh

Type

float, read-only

property bars

All bars in the mesh

Type

BarList, read-only

property bbox

The bounding box of the mesh

Type

steps.API_2.geom.BoundingBox, read-only

intersect(points, sampling=- 1)

Computes the intersection of the current mesh and line segment(s), given their vertices

Parameters
  • points (numpy.ndarray) – A 2-D NumPy array (/memview) of 3D points where each element contains the 3 point coordinates

  • sampling (int) – if > 0, use montecarlo method with sampling points, otherwise use deterministic method.

Returns

A list of lists of tuples representing the intersected tetrahedrons, one element per line segment. Each tuple is made of 2 elements, a tetrahedron global identifier, and its respective intersection ratio.

Return type

List[List[Tuple[TetReference, float]]]

intersectIndependentSegments(points, sampling=- 1)

Similar to the intersect method but here we deal with independent segments, i.e. every two points we have a segment not related to previous or following ones. E.g. seg0 = (points[0], points[1]), seg1 = (points[2], points[3]), etc.

Parameters
  • points (numpy.ndarray) – A 2-D NumPy array (/memview) of 3D points where each element contains the 3 point coordinates

  • sampling (int) – if > 0, use montecarlo method with sampling points, otherwise use deterministic method.

Returns

A list of lists of tuples representing the intersected tetrahedrons, one element per line segment. Each tuple is made of 2 elements, a tetrahedron global identifier, and its respective intersection ratio.

Return type

List[List[Tuple[TetReference, float]]]

property surface

All triangles in the surface of the mesh

Type

TriList, read-only

property tets

All tetrahedrons in the mesh

Type

TetList, read-only

property tris

All triangles in the mesh

Type

TriList, read-only

property verts

All vertices in the mesh

Type

VertList, read-only

class DistMesh(filename, scale=1, mpiComm=None, **kwargs)[source]

Bases: Geometry

Container class for distributed tetrahedral meshes

In contrast to TetMesh, only gmsh mesh files are accepted.

Parameters
  • filename (str) – The Gmsh filename (or path) including any suffix.

  • scale (float) – Optionally rescale the mesh on loading by given factor

  • mpiComm (mpi4py.MPI.Intracomm) – mpi4py MPI communicator, defaults to COMM_WORLD

filename should either be a full path to a single gmsh mesh or a path prefix for pre-partitioned meshes. If a single gmsh file is provided, it will be partitioned automatically. Pre-partitioned meshes can be loaded by providing a path prefix e.g. 'path/to/meshName' and individual partition files should be named 'path/to/meshName_1.msh', 'path/to/meshName_2.msh', etc. The numbering of files starts at 1 and goes up to steps.API_2.sim.MPI.nhosts.

Note that, for pre-partitioned meshes, elements are not necessarily numbered between 0 and n (with n the number of elements).

asGlobal()[source]

Context manager method to use the mesh globally

This should be used with the with statement in the following way:

with mesh.asGlobal():
    lst1 = mesh.tets    # Returns a global list of tetrahedron indices
    lst2 = mesh.surface # Returns a global list of boundary triangles

It is equivalent to using only with mesh: since this defaults to using global indices. with mesh.asGlobal(): is however preferable since it is more explicit.

Return type

DistMesh

asLocal(owned=True)[source]

Context manager method to use the mesh locally

This should be used with the with statement in the following way:

with mesh.asLocal():
    lst1 = mesh.tets    # Returns a local list of tetrahedron indices
    lst2 = mesh.surface # Returns a local list of boundary triangles

All properties and methods that would usually return lists of global indices will return local lists instead when wrapped with this context manager.

In the above example, lst1 would end up having the same value as mesh.tets.toLocal() but the latter is less efficient since it first queries the full list of all tetrahedrons and then restricts it to only the local ones.

Note

This does not affect existing lists or newly created lists. For example, a tetrahedron list created in a with mesh.asLocal() block will still default to containing global elements (unless the local=True keyword argument is given at creation). This context manager only affects the behavior of the mesh object, and its subparts like compartments and patches. For example comp.tets will return a list of local indices in this context manager.

Parameters

owned (bool) – Whether only owned elements should be considered (defaults to True)

Return type

DistMesh

intersect(points, sampling=- 1, raw=False, local=True)[source]

Computes the intersection of the current mesh and line segment(s), given their vertices

Parameters
  • points (numpy.ndarray) – A 2-D NumPy array (/memview) of 3D points where each element contains the 3 point coordinates

  • sampling (int) – if > 0, use montecarlo method with sampling points, otherwise use deterministic method.

  • raw (bool) – If True, return raw integer tetrahedron indices.

  • local (bool) – if False, return global tetrahedron identifiers instead of local ones.

This method can only be called when mesh.asLocal() is used, and it will return local tetrahedron references by default.

Returns

A list of lists of tuples representing the intersected tetrahedrons, one element per line segment. Each tuple is made of 2 elements, a tetrahedron identifier (local or global depending on local), and its respective intersection ratio.

Return type

List[List[Tuple[Union[TetReference, int], float]]]

intersectIndependentSegments(points, sampling=- 1, raw=False, local=True)[source]

Similar to the intersect method but here we deal with independent segments, i.e. every two points we have a segment not related to previous or following ones. E.g. seg0 = (points[0], points[1]), seg1 = (points[2], points[3]), etc.

Parameters
  • points (numpy.ndarray) – A 2-D NumPy array (/memview) of 3D points where each element contains the 3 point coordinates

  • sampling (int) – if > 0, use montecarlo method with sampling points, otherwise use deterministic method.

  • raw (bool) – If True, return raw integer tetrahedron indices.

  • local (bool) – if False, return global tetrahedron identifiers instead of local ones.

This method can only be called when mesh.asLocal() is used, and it will return local tetrahedron references by default.

Returns

A list of lists of tuples representing the intersected tetrahedrons, one element per line segment. Each tuple is made of 2 elements, a tetrahedron identifier (local or global depending on local), and its respective intersection ratio.

Return type

List[List[Tuple[Union[TetReference, int], float]]]

property tets

All tetrahedrons in the mesh

Type

TetList, read-only

property tris

All triangles in the mesh

Type

TriList, read-only

property verts

All vertices in the mesh

Type

VertList, read-only

property tetGroups
property triGroups
property vertGroups
property Vol

Total volume of the mesh

Type

float, read-only

property bars

All bars in the mesh

Type

BarList, read-only

property bbox

The bounding box of the mesh

Type

steps.API_2.geom.BoundingBox, read-only

property surface

All triangles in the surface of the mesh

Type

TriList, read-only

class DiffBoundary(lst, *args, **kwargs)[source]

Bases: NamedObject

Diffusion Boundary class

Annotation of a group of triangles in a Tetmesh. The triangles form a boundary between two compartments, that may allow diffusion of some specified species.

Parameters

lst (TriList) – List of triangles forming the boundary

__call__(direc=None)[source]

Get a version of the diffusion boundary with added information

Parentheses notation (function call) can be used on a diffusion boundary to represent a diffusion boundary further specified with additional information. This should never be needed during model declaration but can become handy for simulation control and data saving (see steps.API_2.sim.SimPath).

Parameters

direc – A direction for the diffusion boundary (i.e. a target compartment for diffusion boundaries and a target patch for surface diffusion boundary)

Returns

A version of the diffusion boundary that is specific to direction direc.

property tris

All triangles in the boundary, if applicable

Type

TriList, read-only

class SDiffBoundary(lst, patch1, patch2, *args, **kwargs)[source]

Bases: DiffBoundary

Surface Diffusion Boundary class

Annotation of a group of bars in a Tetmesh. The bars form a boundary between two patches, that may allow diffusion of some specified species. The patches to be connected by this surface diffusion boundary need to be specified to avoid potential ambiguity.

Parameters
  • lst (BarList) – List of bars forming the boundary

  • patch1 (Patch) – One of the patches connected to the boundary

  • patch2 (Patch) – The other patch connected to the boundary

property bars

All bars in the boundary, if applicable

Type

BarList, read-only

class MeshPartition(mesh, tet_hosts=[], tri_hosts={}, wm_hosts=[])[source]

Bases: Params

Mesh element partition for parallel simulations

Wrapper class that groups tetrahedrons, triangles and well-mixed partitions.

Parameters
  • mesh (TetMesh) – Mesh being partitioned

  • tet_hosts (List[int]) – List of host process indices for each tetrahedron

  • tri_hosts (Mapping[int, int]) – Mapping between triangle index and host process index

  • wm_hosts (List[int]) – List of host process indices for each well mixed compartment

property tetPart

The partition data for tetrahedrons

Type

List[int], read-only

property triPart

The partition data for triangles

Type

Dict[int, int], read-only

property wmPart

The partition data for well mixed locations

Type

List[int], read-only

property tetTable

A dictionary in the format of {host0: TetList([tet0, tet1, …]), …}

Type

Dict[int, TetList], read-only

property triTable

A dictionary in the format of {host0: TriList([tri0, tri1, …]), …}

Type

Dict[int, TriList], read-only

validate()[source]

Validate the partitioning of the mesh

printStats()[source]

Print out partitioning stastics

Returns

(tet_stats, tri_stats, wm_stats, num_hosts, min_degree, max_degree, mean_degree) [tet/tri/wm]_stats contains the number of tetrahedrons/triangles/well-mixed volumes in each hosting process, num_hosts provide the number of hosting processes, [min/max/mean]_degree provides the minimum/maximum/average connectivity degree of the partitioning.

LinearMeshPartition(mesh, xbin, ybin, zbin, default_tris=None)[source]

Partition the mesh by linearly binning along the 3 axes

First partition the tetrahedrons and then compute a triangle partition that matches the tetrahedron one.

Parameters
  • mesh (TetMesh) – The mesh to be partitioned

  • xbin (int) – Number of bins on the x axis

  • ybin (int) – Number of bins on the y axis

  • zbin (int) – Number of bins on the z axis

  • default_tris (TriList) – Optional list of triangles that should be partitioned even if they are not part of any patch

Returns

The partition object

Return type

MeshPartition

Note

The triangle partition only takes into account triangles that are part of a Patch.

MetisPartition(mesh, path, default_tris=None)[source]

Partition the mesh using a Metis .epart file

First partition the tetrahedrons according to the Metis file and then compute a triangle partition that matches the tetrahedron one.

Parameters
  • mesh (TetMesh) – The mesh to be partitioned

  • path (str) – Path to the .epart file generated by Metis

  • default_tris (TriList) – Optional list of triangles that should be partitioned even if they are not part of any patch

Returns

The partition object

Return type

MeshPartition

Note

The triangle partition only takes into account triangles that are part of a Patch.

GmshPartition(mesh, default_tris=None)[source]

Retrive partition information stored in a pre-partitioned gmsh.

The mesh needs to be pre-partitioned and stored in a single Gmsh file with format version 2.2 ascii.

Parameters
  • mesh (TetMesh) – The mesh to be partitioned

  • default_tris (TriList) – Optional list of triangles that should be partitioned even if they are not part of any patch

Returns

The partition object

Return type

MeshPartition

Note

The triangle partition only takes into account triangles that are part of a Patch.

MorphPartition(mesh, morph, scale=1e-06, default_tris=None)[source]

Partition the mesh using morphological sectioning data

First partition the tetrahedrons according to the morph data and then compute a triangle partition that matches the tetrahedron one.

Parameters
  • mesh (TetMesh) – The mesh to be partitioned

  • morph (Morph) – Morphological sectioning data

  • scale (float) – Scaling factor from morphological sectioning data to Tetmesh data, the default value is 1e-6, that is 1 unit (usually micron) in the morph file equals 1e-6 unit (meter) in Tetmesh

  • default_tris (TriList) – Optional list of triangles that should be partitioned even if they are not part of any patch

Returns

The partition object and a mapping between partition index and morphological section name

Return type

Tuple[MeshPartition, Dict[int, str]]

class Reference(idx, mesh=None, **kwargs)[source]

Bases: NamedObject

Base class for all element references.

Parameters
  • idx (int) – Index of the element

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the element

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

The actual index of the element can be accessed with the idx attribute.

toList()[source]

Get a list holding only this element

Returns

A list with only this element

Return type

RefList

property idx

Index of the element

Type

int, read-only

Only available for references defined in distributed tetrahedral meshes:
isLocal()

Return whether the reference uses a local index

Return type

bool

toLocal(owned=True)

Get the local version of the reference

If the reference is local, it returns itself. If the reference is global, and if the element exists in the current MPI process, a local reference to it is returned. Otherwise, None is returned. The returned reference is thus different for each MPI process.

Parameters

owned (bool) – Whether only owned elements should be considered (defaults to True)

Return type

Reference

toGlobal(**_kwargs)

Get the global version of the reference

If the reference is not local, it returns itself. Otherwise, the local index of the original reference is converted to global mesh index and the corresponding global reference is returned.

Return type

Reference

class TetReference(idx, mesh=None, **kwargs)[source]

Bases: Reference

Convenience class for accessing properties of a tetrahedron

Parameters
  • idx (int) – Index of the tetrahedron

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the tetrahedron

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

containsPoint(point)[source]

Check whether a given 3D point is inside the tetrahedron

Parameters

point (List[float]) – The 3D point

Return type

bool

property neighbs

Neighboring tetrahedrons (between 0 and 4)

Type

TetList, read-only

property faces

Faces of the tetrahedron

Type

TriList, read-only

property center

Barycenter of the tetrahedron

Type

Point, read-only

property verts

Vertices of the tetrahedron

Type

VertList, read-only

property comp

Compartment associated to the tetrahedron (if any)

Type

Comp or None, read-only

property Vol

Volume of the tetrahedron

Type

float, read-only

property qualityRER

Radius-edge-ratio (a quality measurement) of the tetrahedron

Type

float, read-only

toList()

Get a list holding only this element

Returns

A list with only this element

Return type

RefList

class TriReference(idx, mesh=None, **kwargs)[source]

Bases: Reference

Convenience class for accessing properties of a triangle

Parameters
  • idx (int) – Index of the triangle

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the triangle

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

property verts

Vertices of the triangle

Type

VertList, read-only

property center

Barycenter of the triangle

Type

Point, read-only

property Area

Area of the triangle

Type

float, read-only

property tetNeighbs

Neighboring tetrahedrons

Type

TetList, read-only

property triNeighbs

Neighboring triangles

Type

TriList, read-only

property patchTriNeighbs

Neighboring triangles that are also part of the same patch If the triangle is not part of a Patch, returns an empty TriList

Type

TriList, read-only

property patch

Patch associated to the triangle (if any)

Type

Patch or None, read-only

property bars

bars of the triangle

Type

VertList, read-only

property norm

Normal vector of the triangle

Type

Point, read-only

toList()

Get a list holding only this element

Returns

A list with only this element

Return type

RefList

class BarReference(idx, mesh=None, **kwargs)[source]

Bases: Reference

Convenience class for accessing properties of a bar

Parameters
  • idx (int) – Index of the bar

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the bar

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

property verts

Vertices of the bar

Type

VertList, read-only

property center
toList()

Get a list holding only this element

Returns

A list with only this element

Return type

RefList

class VertReference(idx, mesh, anonymous=False, *args, **kwargs)[source]

Bases: Reference, Point

Convenience class for accessing properties of a vertex

Can be used in the same way as a Point.

Parameters
  • idx (int) – Index of the vertex

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the vertex

Warning

When used with MPI, keep in mind that the creation of a VertReference retrieves the coordinates of the vertex and thus potentially requires synchronization across MPI ranks. When using DistMesh, there is no need for synchronization if the local=True keyword argument is used.

Note

In contrast to other references, VertReference can only be created by explicitely providing a mesh object, it cannot use the context manager syntax described in TetMesh.

toList()

Get a list holding only this element

Returns

A list with only this element

Return type

RefList

class RefList(lst=None, mesh=None, **kwargs)[source]

Bases: NamedObject

Base convenience class for geometrical element lists

Parameters
  • lst (Iterable[int] or range or Iterable[Reference] or RefList) – The list of elements

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the elements

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

Behaves like a python list but with additional functionalities.

splitToROIs(ROIPrefix, method='grid', gridSize=None)[source]

Split an element list into a square grid and define each square as a Region Of Interest

Parameters
  • ROIPrefix (str) – A prefix for the names of the ROIs, it will be suffixed with ‘x_y_z’ with x, y, and z the integer grid coordinate of a given square ROI.

  • method (str) – The method used to define the ROIs, only ‘grid’ is available for now.

  • gridSize (float) – The side length of a single square from the grid, in meters.

__getitem__(key)[source]

Access one or several element(s) from the list

Parameters

key (int or slice) – Index of the element or slice

Returns

The element(s)

Return type

Reference or RefList

Usage:

>>> lst = RefList(range(10, 20), mesh=mesh)
>>> lst[0] # Returns the 1st element from the list
Tet(10)
>>> lst[0].idx
10
>>> lst2 = lst[1:5] # Returns elements between the 2nd and the 5th (included).
>>> lst2[0]
Tet(11)

Note

The index of an element in the list is in general different from its index in the mesh.

__iter__()[source]

Get an iterator over the list elements

Usage:

>>> lst = RefList(range(3), mesh=mesh)
>>> for ref in lst:
...     print(ref.idx)
...
0
1
2
>>> [ref.idx for ref in lst]
[0, 1, 2]
append(e)[source]

Append an element to the list

Parameters

e (Reference) – The element

Note

The type of the element depends on the type of the list, one cannot add a TriReference to a TetList.

remove(e)[source]

Remove an element from the list

Parameters

e (Reference) – The element

__contains__(elem)[source]

Check whether the list contains an element

Parameters

elem (Reference) – The element

Return type

bool

Usage:

>>> ref1 = Reference(15, mesh=mesh)
>>> ref2 = Reference(25, mesh=mesh)
>>> lst = RefList(range(20), mesh=mesh)
>>> ref1 in lst
True
>>> ref2 in lst
False
__len__()[source]

Get the length of the list

Return type

int

Usage:

>>> lst = RefList(range(5), mesh=mesh)
>>> len(lst)
5
index(elem)[source]

Return the index of the element elem in this list

Parameters

elem (Reference) – The element

Returns

The index of the element in the class (not its index in the mesh)

Return type

int

Behaves like list.index().

__and__(other)[source]

Compute the intersection between two lists

Returns

All the elements that are in both lists. Keeps the order of the left operand.

Return type

RefList

Usage:

>>> lst1 = RefList(range(0, 10), mesh=mesh)
>>> lst2 = RefList(range(5, 15), mesh=mesh)
>>> lst3 = lst1 & lst2
>>> [ref.idx for ref in lst3]
[5, 6, 7, 8, 9]
__or__(other)[source]

Compute the union of two lists

Returns

All the elements of the left operand followed by all the elements of the right operand that are not in the left one.

Return type

RefList

Usage:

>>> lst1 = RefList(range(0, 10), mesh=mesh)
>>> lst2 = RefList(range(5, 15), mesh=mesh)
>>> lst3 = lst1 | lst2
>>> [ref.idx for ref in lst3]
[0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13 ,14]
__sub__(other)[source]

Compute the substraction of one list from another

Returns

All elements that are in the left operand but are absent from the right operand.

Return type

RefList

Usage:

>>> lst1 = RefList(range(0, 10), mesh=mesh)
>>> lst2 = RefList(range(5, 15), mesh=mesh)
>>> lst3 = lst1 - lst2
>>> [ref.idx for ref in lst3]
[0, 1, 2, 3, 4]
>>> lst4 = lst2 - lst1
>>> [ref.idx for ref in lst4]
[10, 11, 12, 13, 14]
__xor__(other)[source]

Compute the symetric difference of two lists

Returns

All elements in the left operand that are not in the right operand followed by all elements of the right operand that are not in the left operand.

Return type

RefList

Usage:

>>> lst1 = RefList(range(0, 10), mesh=mesh)
>>> lst2 = RefList(range(5, 15), mesh=mesh)
>>> lst3 = lst1 ^ lst2
>>> [ref.idx for ref in lst3]
[0, 1, 2, 3, 4, 10, 11, 12, 13, 14]
>>> lst4 = lst2 ^ lst1
>>> [ref.idx for ref in lst4]
[10, 11, 12, 13, 14, 0, 1, 2, 3, 4]
__add__(other)[source]

Concatenate two lists

Returns

All elements of the left operand followed by all elements of the right operand.

Return type

RefList

Usage:

>>> lst1 = RefList(range(0, 10), mesh=mesh)
>>> lst2 = RefList(range(5, 15), mesh=mesh)
>>> lst3 = lst1 + lst2
>>> [ref.idx for ref in lst3]
[0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14]
>>> lst4 = lst2 + lst1
>>> [ref.idx for ref in lst4]
[5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 0, 1, 2, 3, 4, 5, 6, 7, 8, 9]
__eq__(other)[source]

Test for the equality of two lists

Returns

True if the lists have the same elements in the same order, False otherwise

Return type

bool

Usage:

>>> lst1 = RefList(range(0, 5), mesh=mesh)
>>> lst2 = RefList(range(0, 10), mesh=mesh)
>>> lst1 == lst2
False
>>> lst1 == lst2[0:5]
True
property indices

The indices of elements in the list

Type

List[int], read-only

Only available for lists defined in distributed tetrahedral meshes:
isLocal()

Return whether the list is composed of local element indices

Returns

True if the list is local

Return type

bool

toLocal(owned=True)

Get the local version of a distributed list

If the list is local, it returns itself. Otherwise, the orginal list is filtered and only the elements that are hosted by the current MPI process are added to the local list, with their local index. The returned list is thus different for each MPI process.

Parameters

owned (bool) – Whether only owned elements should be considered (defaults to True)

Returns

A local version of the list

Return type

RefList

toGlobal(**_kwargs)

Get the global version of a distributed list

If the list is not local, it returns itself. Otherwise, the local indices of the original list are converted to global mesh indices and the list of these global indices is returned.

Returns

A global version of the list

Return type

RefList

combineWithOperator(binOp)

Combine the local versions of this list across ranks with the given operator

Parameters

binOp (Callable[[RefList, RefList], RefList]) – The binary operator used to combine the global lists

Returns

The global element list resulting from the combination of local lists

Return type

RefList

class TetList(lst=None, mesh=None, **kwargs)[source]

Bases: RefList

Convenience class for tetrahedron lists

Parameters
  • lst (Iterable[int] or range or Iterable[TetReference] or TetList) – The list of tetrahedrons

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the elements

If lst is omitted, the list will be initialized to an empty list.

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

Behaves like a RefList but with additional properties specific to tetrahedrons. Note that __getitem__() and __contains__() are modified to accept a 3D point as parameter.

__getitem__(key)[source]

Access one or several tetrahedron(s) from the list

Parameters

key (int or slice or Tuple[float, float, float] or Point) – Index of the tetrahedron, or slice, or 3D point

Returns

The tetrahedron(s)

Return type

TetReference or TetList

Usage:

>>> lst = TetList(range(10, 20), mesh=mesh)
>>> lst[0] # Returns the 1st element from the list
Tet(10)
>>> lst[0].idx
10
>>> lst2 = lst[1:5] # Returns elements between the 2nd and the 5th (included).
>>> lst2[0]
Tet(11)
>>> mesh.tets[0, 0, 0] # Returns the tetrahedron which contains the point x=0, y=0, z=0
Tet(123)

Note

The index of an element in the list is in general different from its index in the mesh. TetList is the only subclass of RefList that can take a 3D point as key of its __getitem__ special method.

__contains__(key)[source]

Check whether the list contains a tetrahedron

Parameters

key (TetReference or Tuple[float, float, float] or Point) – The element or a 3D point

Return type

bool

In addition to the arguments accepted in RefList.__contains__(), it is possible to check whether a 3D point is inside one tetrahedron in the list.

Usage:

>>> ref1 = TetReference(15, mesh=mesh)
>>> ref2 = TetReference(25, mesh=mesh)
>>> lst = TetList(range(20), mesh=mesh)
>>> ref1 in lst
True
>>> ref2 in lst
False
>>> Point(0, 0, 0) in lst
False
>>> Point(1e-6, 0, 5e-6) in lst
True
dilate(d=1)[source]

Dilates the list to tetrahedron neighbors

One cycle of dilation corresponds to adding all tetrahedrons that are neighbors of the current tetrahedron list but are not part of the list. This operation is repeated d times.

Parameters

d (int) – Topological distance to grow, defaults to 1

erode(d=1)[source]

Erodes the list, removing surface tetrahedrons

One cycle of erosion corresponds to removing all tetrahedrons that are on the surface of the list (i.e. at least one of their four neighbors is not in the list). This operation is repeated d times.

Parameters

d (int) – Topological distance to erode, defaults to 1

property surface

Surface of the tetrahedrons in the list

I.e. all triangles that are not shared by two tetrahedrons in the list.

Type

TriList, read-only

property Vol

The summed volume of all tetrahedrons in the list

Note

If the list contains duplicate elements, they will be counted only once in the total volume.

Type

float, read-only

property tris

All faces of all tetrahedrons in the list

Type

TriList, read-only

property bars

All bars of all tetrahedrons in the list

Type

BarList, read-only

property verts

All vertices of all tetrahedrons in the list

Type

VertList, read-only

append(e)

Append an element to the list

Parameters

e (Reference) – The element

Note

The type of the element depends on the type of the list, one cannot add a TriReference to a TetList.

index(elem)

Return the index of the element elem in this list

Parameters

elem (Reference) – The element

Returns

The index of the element in the class (not its index in the mesh)

Return type

int

Behaves like list.index().

remove(e)

Remove an element from the list

Parameters

e (Reference) – The element

splitToROIs(ROIPrefix, method='grid', gridSize=None)

Split an element list into a square grid and define each square as a Region Of Interest

Parameters
  • ROIPrefix (str) – A prefix for the names of the ROIs, it will be suffixed with ‘x_y_z’ with x, y, and z the integer grid coordinate of a given square ROI.

  • method (str) – The method used to define the ROIs, only ‘grid’ is available for now.

  • gridSize (float) – The side length of a single square from the grid, in meters.

class TriList(lst=None, mesh=None, **kwargs)[source]

Bases: RefList

Convenience class for triangle lists

Parameters
  • lst (Iterable[int] or range or Iterable[TriReference] or TriList) – The list of triangles

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the elements

If lst is omitted, the list will be initialized to an empty list.

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

Behaves like a RefList but with additional properties specific to triangles.

property edges

Perimeter of the triangles in the list

I.e. all bars that are not shared by two triangles in the list.

Type

BarList, read-only

property bars

All bars of all triangles in the list

Type

BarList, read-only

property verts

All vertices of all triangles in the list

Type

VertList, read-only

property Area

The summed area of all triangles in the list

Note

If the list contains duplicate elements, they will be counted only once in the total area.

Type

float, read-only

append(e)

Append an element to the list

Parameters

e (Reference) – The element

Note

The type of the element depends on the type of the list, one cannot add a TriReference to a TetList.

index(elem)

Return the index of the element elem in this list

Parameters

elem (Reference) – The element

Returns

The index of the element in the class (not its index in the mesh)

Return type

int

Behaves like list.index().

remove(e)

Remove an element from the list

Parameters

e (Reference) – The element

splitToROIs(ROIPrefix, method='grid', gridSize=None)

Split an element list into a square grid and define each square as a Region Of Interest

Parameters
  • ROIPrefix (str) – A prefix for the names of the ROIs, it will be suffixed with ‘x_y_z’ with x, y, and z the integer grid coordinate of a given square ROI.

  • method (str) – The method used to define the ROIs, only ‘grid’ is available for now.

  • gridSize (float) – The side length of a single square from the grid, in meters.

class BarList(lst=None, mesh=None, **kwargs)[source]

Bases: RefList

Convenience class for bar lists

Parameters
  • lst (Iterable[int] or range or Iterable[BarReference] or :py:class:BarList`) – The list of bars

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the elements

If lst is omitted, the list will be initialized to an empty list.

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

Behaves like a RefList but with additional properties specific to bars.

property verts

All vertices of all bars in the list

Type

VertList, read-only

append(e)

Append an element to the list

Parameters

e (Reference) – The element

Note

The type of the element depends on the type of the list, one cannot add a TriReference to a TetList.

index(elem)

Return the index of the element elem in this list

Parameters

elem (Reference) – The element

Returns

The index of the element in the class (not its index in the mesh)

Return type

int

Behaves like list.index().

remove(e)

Remove an element from the list

Parameters

e (Reference) – The element

splitToROIs(ROIPrefix, method='grid', gridSize=None)

Split an element list into a square grid and define each square as a Region Of Interest

Parameters
  • ROIPrefix (str) – A prefix for the names of the ROIs, it will be suffixed with ‘x_y_z’ with x, y, and z the integer grid coordinate of a given square ROI.

  • method (str) – The method used to define the ROIs, only ‘grid’ is available for now.

  • gridSize (float) – The side length of a single square from the grid, in meters.

class VertList(lst=None, mesh=None, **kwargs)[source]

Bases: RefList

Convenience class for vertex lists

Parameters
  • lst (Iterable[int] or range or Iterable[VertReference] or :py:class:VertList`) – The list of vertices

  • mesh (Union[TetMesh, DistMesh, None]) – Mesh object that contains the elements

If lst is omitted, the list will be initialized to an empty list.

If mesh is not provided, it will be inferred, if possible, from context managers (see TetMesh).

Behaves like a RefList.

append(e)

Append an element to the list

Parameters

e (Reference) – The element

Note

The type of the element depends on the type of the list, one cannot add a TriReference to a TetList.

index(elem)

Return the index of the element elem in this list

Parameters

elem (Reference) – The element

Returns

The index of the element in the class (not its index in the mesh)

Return type

int

Behaves like list.index().

remove(e)

Remove an element from the list

Parameters

e (Reference) – The element

splitToROIs(ROIPrefix, method='grid', gridSize=None)

Split an element list into a square grid and define each square as a Region Of Interest

Parameters
  • ROIPrefix (str) – A prefix for the names of the ROIs, it will be suffixed with ‘x_y_z’ with x, y, and z the integer grid coordinate of a given square ROI.

  • method (str) – The method used to define the ROIs, only ‘grid’ is available for now.

  • gridSize (float) – The side length of a single square from the grid, in meters.

class Point(x, y, z)[source]

Bases: ndarray

Convenience class for representing 3D points

This class inherits from numpy.ndarray and can thus be used in the same way as a numpy array. The only difference is the possibility to access invidual coordinates through the x, y, and z properties.

Parameters
  • x (float) – x coordinate

  • y (float) – y coordinate

  • z (float) – z coordinate

property x

x coordinate

Type

float, read-only

property y

y coordinate

Type

float, read-only

property z

z coordinate

Type

float, read-only

class BoundingBox(minP, maxP)[source]

Bases:

Convenience class for holding bounding box information

property min

Minimum point of the bounding box

Type

Point, read-only

property max

Maximum point of the bounding box

Type

Point, read-only

property center

Center point of the bounding box

Type

Point, read-only

class Morph(sections)[source]

Bases:

Morphological sectioning data

Contains data describing the sections of a NEURON morphology.

This class should not be instantiated directly, it should be created using one of the dedicated class methods.

classmethod Load(path)[source]

Load morphological data from a morph file

Parameters

path (str) – Path to the file

Returns

The morphological sectioning data

Return type

Morph

classmethod LoadHOC(path)[source]

Load morphological data from a NEURON hoc file

Parameters

path (str) – Path to the file

Returns

The morphological sectioning data

Return type

Morph

classmethod LoadSWC(path)[source]

Load morphological data from an swc file

Parameters

path (str) – Path to the file

Returns

The morphological sectioning data

Return type

Morph

Save(path)[source]

Save to a morph file

Parameters

path (str) – Path to the file