pydmga Package

pydmga Package

This is a Python interface to DMGAlpha library. It has more functions and it is much simpler than the bare C++ counterpart. It uses dmga2py module to communicate with C++ library dmga.

container Module

class pydmga.container.Container(geometry)

This holds particles in simulations

add(*args)

Adds elements to this container.

May work on many different elements: * add(tuple): tuple must contain 5 elements: (id, x, y, z, r), id is integer, others are double * add(iterable): iterable must contain tuples (as above) * add(id, z, y, z, r): elements as in tuple above * add(id, z, y, z): as above, but r is assumed to be 1.0

Rises TypeError when elements are not tuples.

add_by_coords(id, x, y, z, r)

for internal use, adds one element

add_raw(id, x, y, z, r)

This adds one element without doing transformation on coordinates for boosting inserts, when we know that no transform is needed and we really need the speed...

find(id)

Returns element with given ID

Parameters:id – is integer
get(number)

Get the i-th inserted particle as a tuple (id, x, y, z)

Parameters:number – integer in range [0, size())
Returns:a particle with a given number

NOTE number is not an ID! To get element by ID use find(ID)

size()
Returns:number of particles in this container
class pydmga.container.ContainerIterator(container)

used to iterate over Container like:

for particle in container:
        pass
Parameters:container – Container for iteration
next()

diagram Module

class pydmga.diagram.Diagram(container, cache_on=False)

This represents Voronoi Diagram of spheres (Power Diagram). You can iterate over all cells in diagram by:

for cell in diagram:
        pass

where cell is of class model.Cell

get_all_cells()
Returns:all cells as list (may be heavy and time consuming!)

NOTE for generating cells use iterator (for c in diagram: pass)

get_cell(number)

returns i-th cell (if necessary, compute it, or use cache)

Parameters:i – integer in range [0, size())
Returns:i-th cell as a object of class Cell
get_cells(subset=None)

Returns cells given by subset.

Parameters:subset – iterable that holds integers in range [0, size()) or None (default: None)
  1. NOTE if subset=None is given then all cells are returned.
  2. NOTE get_cells returns a list (not generator, may be heavy and time consuming!)
  3. NOTE for generating cells use iterator (for c in diagram: pass)
size()
Returns:number of particles (cells) in this diagram
class pydmga.diagram.DiagramIterator(diagram)

To iterate over cells in diagram. Used internally in Diagram.__iter__()

next()

returns next cell or rise StopIteration if no cells

geometry Module

class pydmga.geometry.BoxGeometry(x_min, y_min, z_min, x_max, y_max, z_max, x_periodic=False, y_periodic=False, z_periodic=False)

Bases: pydmga.geometry.Geometry

This is orthogonal geometry given by bounding box (x_min, y_min, z_min, x_max, y_max, z_max). It may be periodic in any direction.

Parameters:
  • x_min,y_min,z_min,x_max,y_max,z_max – Float numbers
  • x_periodic,y_periodic,z_periodic – boolean (default: False - not periodic)
on_boundary(neighbour)
Retrun :True if neighbour id < 0, False otherwise
class pydmga.geometry.CastCylinderGeometry(radius, height, z_periodic=True, precision=32)

Bases: pydmga.geometry.Geometry

This geometry cast particles onto side of a cylinder. It may be useful in investigating properties of particular packaging of particles (haxagonal phase in DDPC or something similar) the cylinder is oriented “upwards” along z direction, and you must specify radius of the base > 0 and its height > 0. It is by default periodic in z direction. NOTE that if z_periodic is set to False, then particles with z > height or z < 0 will not be stored in the container!

Parameters:
  • radius – double > 0
  • height – double > 0
  • z_periodic – boolean (default: True - periodic in z direction)
  • precision – integer and defines the precision of the cylinder - higher the better but slower (default: 32)
adjust_area(area)

It can rigorously adjust computed area.

on_boundary(neighbour)
Retrun :True if neighbour id < 0, False otherwise
transform(id, x, y, z, r)

Casts particle onto the side of the cilinder

Particles close to origin (x = 0, y = 0) will not be transformed

class pydmga.geometry.CastOrthoPlaneGeometry(from_geometry, vec)

Bases: pydmga.geometry.OrthogonalGeometry

This class defines a cast onto the given plane (xy), (xz) or (xz). It creates cast geometry around existing OrthogonalGeometry (It is something like Decorator)

Parameters:
  • from_geometry – is OrthogonalGeometry
  • vec – is one of (1,0,0), (0,1,0), (0,0,1)
on_boundary(neighbour)

may be used to filter only geometry on the casting plane

transform(id, x, y, z, r)

cast particle to a defined plane

class pydmga.geometry.CastPlaneGeometry(from_geometry, cast_vector, x_periodic=False, y_periodic=False)

Bases: pydmga.geometry.OrthogonalGeometry

NotImplementedYet due to some problems with definition of the cast

transform(id, x, y, z, r)
class pydmga.geometry.CastSphereGeometry(radius, origin=(0, 0, 0))

Bases: pydmga.geometry.Geometry

This is rather articicial but still may be useful. It defines geometry that cast all particles on the sphere with given radius > 0 and in given origin (default: (0,0,0)). Bounding box is computed as the minimal box containing given sphere. It redefines transform() to cast particles on sphere. it cannot be periodic (you do not have a way to define it properly though).

Parameters:
  • radius – is double > 0
  • origin – is tuple (default: (0,0,0))
on_boundary(neighbour)
Retrun :True if neighbour id < 0, False otherwise
transform(id, x, y, z, r)

Cast particle to the sphere taking into account particles that lies near origin (those particles are not cast - so be aware)

class pydmga.geometry.Geometry

Base (abstract) class for all geometry objects. Any geometry object must have bounding_box property that defines rectangular box that the geometry lies inside. bounding_box is a set of 6 float numbers (min_x, min_y, min_z, max_x, max_y, max_z).

adjust_area(area)

Default does nothing. It is useful for approximate geometries, where real area may be recovered in some way from approximate one.

Parameters:area – approximate area obtained from computations
Returns:adjusted area, should be exact value of the area
adjust_volume(volume)

Default does nothing. It is useful for approximate geometries, where real volume may be recovered in some way from approximate one.

Parameters:volume – approximate volume obtained from computations
Returns:adjusted volume, should be exact value of the volume
on_boundary(neighbour)

Tests if a side with a given neighbour is on boundary of the geometry (was created with some artificial wall not related to other particle) usually it is neighbour < 0.

Parameters:neighbour – a cell number, the neighbour responsible for generating some side of a cell.
transform(id, x, y, z, r)

Transformation of the input coordinates. Usually it is identity, but some geometries (eg. cast geometries) may override this to do cast onto some surface.

NOTE if writing your own Geometry you have an option to change id dynamically

class pydmga.geometry.OrthogonalGeometry(x_size, y_size, z_size, x_periodic=False, y_periodic=False, z_periodic=False)

Bases: pydmga.geometry.Geometry

Basic geometry for Molecular Dynamic simulations. It defines rectangular box as in PDB files (that is with x_size, y_size, z_size).

The bounding box for this geometry is (0.0, 0.0, 0.0, x_size, y_size, z_size)

Parameters:
  • x_size,y_size,z_size – Floats, size of the box
  • x_periodic,y_periodic,z_periodic – boolean (default: False - not periodic)
on_boundary(neighbour)
Retrun :True if neighbour id < 0, False otherwise
pydmga.geometry.random() → x in the interval [0, 1).

model Module

class pydmga.model.Cell(cell_handle)

Represents Voronoi Cell

It is created by Diagram class

Parameters:cell_handle – handle to C++ object representing the cell

it has three public attributes, which are collections of sides, edges and vertices.

It allows for constructs such as:

for side in cell.sides:
    print(side.area())
    print(side.as_coords())

print(cell.sides.size())

for edge in cell.edges:
    print(edge.as_coords())
    print(edge.id)

print cell.edges.size()

for vertex in cell.vertices:
    print(vertex.as_coords())
    print(vertex.id)

print(cell.vertices.size())
area()
Returns:area of this cell
get_sides()

returns all sides as one list, thus may be time and space consuming the resulting list contains CellSide objects. for generating better use sides iterable collection:

for (neighbour, side) in cell.sides:
    pass

where neighbour is ??? and side is ???

get_sides_as_coords()

returns all sides as one list, thus may be time and space consuming the resulting list contains only lists of 3D coordinates of vertices of each side. for generating better use sides iterable collection:

for side in cell.sides:
    pass
volume()
Returns:volume of this cell
class pydmga.model.CellEdge(parent_cell_handle, id, current, inverse)

This class represents Edge of the Power Diagram :param parent_cell_handle: handle to C++ cell structure that contains this edge :param id: unique id of this particular edge (each edge has unique id from 0 to cell.edges.size()-1) :param current: current = (v,j), v is first vertex index (ID) inside the cell, j is edge index at v :param inverse: current = (u,i), u is second vertex index (ID) inside the cell, i is edge index at u

In inverse, i is the edge pointing backwards, i.e. from u to v

We hold inverse to allow faster retrieval of inverse edge (useful)

as_coords()
Returns:tuple of tuples: ((x_v, y_v, z_v), (x_u, y_u, z_u))

may be used as a simple representation of an interval v -> u in 3D space

as_pair()
Returns:(v, u)
as_pair_coords()
Returns:tuple of tuples: ((x_v, y_v, z_v), (x_u, y_u, z_u))

may be used as a simple representation of an interval v -> u in 3D space

Deprecated :this is replaced by as_coords
as_tuple()
Returns:tuple (id, v, j, u, i)
first()
Returns:v - id of the first vertex
first_coords()
Returns:coordinates of the first vertex as a tuple (x,y,z)
second()
Returns:u - id of the second vertex
second_coords()
Returns:coordinates of the second vertex as a tuple (x,y,z)
class pydmga.model.CellEdgesCollection(cell_handle)

This represents iterable collection of edges in a Voronoi cell

Parameters:cell_handle – handle to associated C++ cell structure

NOTE Should not be created directly. Use cell.edges to get this.

as_coords()

allows iteration:

for (vertex_1_coords, vertex_2_coords) in edges_collection.as_coords():
    pass

where vertex_1_coords and vertex_2_coords are tuples (x,y,z)

size()
Returns:number of edges in this cell
class pydmga.model.CellEdgesIterator(edges, as_coords=False)

This is for internal use to allow iteration over Edges in CellEdgesCollection

Parameters:
  • edges – reference to edges collection of class CellEdgeCollection
  • as_coords – if true, then next() will return pairs of coordinate tuples instead of CellEdge class instances, default: False

NOTE this should not be created directly. Use cell.edges to get this iterator

next()

go to the next element or throw StopIteration

if as_coords = False then CellEdge class instance is returned, otherwise tuple of coordinate tuples is returned.

class pydmga.model.CellSide(cell_handle, number, neighbour, vertices)

represents a single side of a Voronoi cell

Parameters:
  • cell_handle – handle to associated C++ cell structure
  • number – the number of the face (unique ID from 0 to cell.sides.size()-1)
  • neighbour – number of the neighbour who created this side
  • vertices – list of vertices on this side

NOTE Should not be created directly. Use cell.sides to get this.

area()
Returns:area of this side
as_coords()
Returns:this side as list of coordinates (x,y,z) of all vertices
as_list()
Returns:this side as list of vertex indices
size()

returns number of vertices on this side

class pydmga.model.CellSideIterator(cell_side)

This is for internal use to allow iteration over vertex coordinates in CellSide

next()

returns coordinates of the next vertex as a tuple (x,y,z)

class pydmga.model.CellSidesCollection(cell_handle)

This represents iterable collection of sides in a Voronoi cell

as_coords()

allows iteration:

for side in side_collection:
    pass

where side is a list of coordinate tuples (x,y,z) of consecutive vertices.

size()

returns number of sides (neighbours) in this cell

class pydmga.model.CellSidesIterator(sides, as_coords=False)

This is for internal use to allow iteration over Sides in CellEdgesCollection it returns sides as lists of CellSide objects or as lists of coordinates.

Parameters:
  • sides – is CellSidesCollection
  • as_coords – is Boolean and controls if iterator should return coords (True) or CellSide class instances (False) (default: False)
next()

go to the next element or throw StopIteration

class pydmga.model.CellVertex(_parent_cell_handle, id, coords)

This is a class that repreents single vertex of a Voronoi structure, that is the single point in space with some id, that can be used to create connection relations (edges).

Parameters:
  • _parent_cell_handle – handle to a C++ structure of a parent cell
  • id – the number of this vertex (id) within a cell
  • coords – tuple of coordinates (x,y,z)

NOTE this class should not be created directly, you should access vertices by routines in some Diagram.

as_coords()
Returns:this object as a tuple (x, y, z)
as_tuple()
Returns:this object as a tuple (id, x, y, z)
class pydmga.model.CellVertexCollection(cell_handle)

This represents iterable collection of vertices in a Voronoi cell

this is usually accessible by Cell.vertices property.

Parameters:cell_handle – handle to C++ Cell counterpart

NOTE this should not be created directly

as_coords()

allows iteration:

for v in vertex_collection.as_coords():
    pass

where v is simple tuple of coordinates (x,y,z)

size()
Returns:number of vertices in this cell
class pydmga.model.CellVertexIterator(vertices, as_coords=False)

This is for internal use to allow iteration over Vertices in CellVertexCollection

Parameters:
  • vertices – it is a CellVertexCollection
  • as_coords – if True then the iteration will be returning tuples (x,y,z) of vertex rather than more involved CellVertex objects, default: False

NOTE using ac_coords=True will be usually faster

NOTE this should not be created directly, use CellVertexCollection from Cell to get this iterator.

next()

returns next element

shape Module

class pydmga.shape.AlphaShape(diagram_or_container)

Bases: pydmga.shape.Shape

This class represents a description of the Alpha Shape of a given Diagram

Parameters:diagram_or_container – Diagram or Container class to use for computations

NOTE if Container class is used to instantiate then the class will create Diagram for internal use

get_cell(number)
Returns:one cell of this shape

The class of shape depends on the type of the shape in Alpha Shape it is Alpha Complex (subset of...) for this cell

max_alpha_threshold()
Returns:max alpha_threshold - value of alpha at which all elements from diagram belongs to the alpha-complex for all cells

NOTE max_threshold is computed only for those cells that were actually created by for example get_cell() or get_cells() or if diagram has cache_on property set to ON (True).

size()
Returns:number of cells in this shape, the same as size in diagram
class pydmga.shape.AlphaShapeCell(shape_cell_handle, cell_handle, complex_handle)

Bases: pydmga.model.Cell

its the extension of a base Cell, it provides the same operations as Cell (for convenience) and defines an alpha shape of a single cell, i.e. the sorted list of 1, 2 or 3-D simplexes with their alhpa threshold - value of alpha at which alpha ball with r’^2 = r^2 + alpha ‘kills’ this simplex (i.e. when intersection of ball and siplex is nonempty). In our computation we use dual of the Delaunay triangulation - Power diagram, but we also allow for non-general positions of balls. The shape is accessible by shape(alpha) method, if no alpha specified then we assume alpha = max_alpha_threshold

NOTE This class should not be instantiated directly - it is returned from AlphaShape by access functions.

Parameters:
  • shape_cell_handle – tells which shape_cell C++ handle to use
  • cell_handle – tells which Voronoi cell is used (for which Voronoi Cell), this is C++ handle
  • complex_handle – tells which complex is used (should be relative to cell_handle), this is C++ handle
max_alpha_threshold()
Returns:max alpha_threshold - value of alpha at which all elements from this cell belongs to alpha complex
shape(alpha=None)

allows for iteration over all simplices in shape:

for simplex in alpha_cell.shape():
    pass

or only for those simplices that has alpha lower than given alpha threshold:

for simplex in alpha_cell.shape(10.0):
    #iterates only over simplices with alpha <= 10.0
    pass
class pydmga.shape.AlphaShapeCellIterator(alpha_shape_cell, alpha, start=0)

To iterate over all simplexes in given cell used internally in AlphaShapeCell.shape()

The iteration is done alpha-wise,m that is from the smallest alpha to the biggest.

Parameters:
  • alpha – Float, this is the upper bund on the elements we want to iterate, so only elements with element.alpha <= alpha will be present.
  • start – iteration starts here, it is sometimes useful.

example:

AlphaShapeCellIterator(shape_cell, 5.0, 10)

will create an iterator for all elements with alpha <= 5.0. Suppose there are 100 such elements. Then the iteration will start from the 10-th element.

NOTE this class should not be instantiated directly - use AlphaCellShape iteration to get it.

next()

returns next cell simplex or rise StopIteration if no cells

class pydmga.shape.AlphaSimplex(alpha_shape_cell_handle, dim, v, j, alpha_threshold, id, is_killed)

represents simple element of Power Diagram with given alpha in the coressponding alpha complex. As we construct Voronoi Diagram, we do not have literally the Alpha Complex (which is subset of Delaunay Triangulation) but we can use duality. Of course we have duality only in case of General Position of balls, but for us we only care for the time when given element of Power Diagram is touched by a growing ball, so we can treat equally elements for balls not in general position as in the case of standard elements.

We also need to note the coresspondence (PG is Power Diagram, RT is regular triangulation (Delaunay)): * VERTEX in PG <-> Tetrahedron in RT * EDGE in PG <-> Triangle in RT * SIDE in PG <-> Edge in RT

Parameters:
  • alpha_shape_cell_handle – is a handle to DMGA C++ object which holds the AlphaShapeCell structure
  • dim – is the dimension of this simplex
  • v – is a base vertex for the this simplex. For vertex it is simply this vertex, for Edge it is base vertex, for side it is starting vertex.
  • j – is edge number at v at which this simplex starts, for vertex it should be ignored, for edge it is the edge, for side it is the first edge on this side from v
  • alpha_threshold – is alpha value at which this simplex is touched by the ball with radius^2 = r^2 + alpha_threshold
  • id – is the inner id within the cell of this simplex, for vertex it is this vertex number, for edge it is id of the edge (unique) and for side it is the number of side with standard iteration
  • is_killed – is a flag that tells if this alpha_threshold is not the smallest power distance form ball center, but other alpha of other simplex. It happens when nearest point doesn’t belong to diagram part of this simplex.
  • NOTE In the structure j will be stored in edge property.
  • NOTE In the structure v will be stored in vertex property.
  • NOTE In the structure dim will be stored in dimension property.

NOTE This class should not be instantiated directly - you should use AlphaShapeCell iteration ability or accessors to get it

DIM_EDGE = 2
DIM_SIDE = 1
DIM_VERTEX = 3
as_dict()

returns this element as a dictionary:

"dimension": self.dimension, 
"vertex": self.vertex, 
"edge": self.edge, 
"alpha_threshold": self.alpha_threshold, 
"id": self.id, 
"is_killed": self.is_killed
as_dict_short()

returns this element as a dictionary:

"d": self.dimension, 
"v": self.vertex, 
"e": self.edge, 
"a": self.alpha_threshold, 
"i": self.id, 
"k": self.is_killed
as_tuple()

returns this element as a tuple:

(dimension, vertex, edge, alpha_threshold, id, is_killed)
diagram_part()
Returns:corresponding part of the diagram

Returned value is either vertex as model.CellVertex, side as model.CellSide, or edge model.CellEdge

class pydmga.shape.SASAArc(first, on_plane, on_sphere, second)

This class describes one arc of the SASA shape, that is an arc on the surface of the sphere.

:param first:first end point of the line :param on_plane: the central point of the arc lying on the Voronoi plane (plane containing some side) :param on_sphere: the central point on the surface of the sphere :param second: second end point

on_plane and on_sphere points may be used to calculate various quentities of the arc (angles, areas, etc) and to draw the arc.

class pydmga.shape.SASACell(shape_cell_handle, cell_handle, contour_handle, number)

Bases: pydmga.model.Cell

This represents a SASA contours description of a given cell. Those contours may be used to compute excluded/included area and volume. The excluded area is an area of the sphere surface outside the voronoi cell, included is the inverse. The same for the volume. Sum of all included area is a SAS Area of the system.

TODO: this should be SASShape? (naming)

border_arcs()

return: iterator over all border arcs

Border arcs are those arcs that form the border of the SAS Contour on the surface of the sphere. It is a result of intersection of the Sphere with the Voronoi Cell.

This iterator should be used when one needs to draw (visualize) the contour

excluded_area()
Returns:excluded area - the area of the sphere surface outside the cell
excluded_volume()
Returns:excluded volume - the volume of the ball outside the cell
polygon_arcs()

return: iterator over all polygon arcs

Polygon arcs are those arcs that forms the polygons inside the excluded area. Notice that border arcs form closed curve consisting of parts of circles (thus name arcs was given). Joining each point connecting two different arcs to the centers of the balls on the sphare surface with great arcs we get parts of the domes and the rest is a (sum of) spherical polygons. The border of spherical polygons consists of parts of those great arcs from the point to the center on sphere. The center on sphere is a pint on sphere where the line connecting two neighbouring particles crosses the sphere.

sas_area()
Returns:included area - the area of the sphere surface contained inside the cell
sas_volume()
Returns:included volume - the volume of the ball volume contained inside the cell
shape(mode=0)

returns the shape iterator over this SASAShape. :param mode: defines which shape to return, default is Border Arcs (mode=0), the other is polygon_arcs (mode=1)

class pydmga.shape.SASAContoursIterator(contour_iterator_handle)

To iterate over all simplexes in given cell used internally in AlphaCell.__iter__()

Parameters:contour_iterator_handle – C++ handle to an iterator
next()

returns next contour

class pydmga.shape.SASAShape(diagram_or_container)

Bases: pydmga.shape.Shape

This is a Solvent Accessible Surface Shape (SAS) for a given Diagram. It is used mainly to compute Solvent Accessible Surface Area (SASA)

SAS is a shape of the union of balls and thus is strongly connected with Alpha SHape. Alpha SHape is more like a topological description of the Shape (what is connected with what), while SAS is the description of the actual union of balls. If you are interested only in topological aspect we suggest Alpha Shapes as those are a bit faster to compute.

Parameters:diagram_or_container – Diagram or Container class to use for computations

NOTE if Container class is used to instantiate then the class will create Diagram for internal use

get_cell(number)
Returns:one cell of this shape

the class of shape depends on the type of the shape in SASA Shape it is SASACell for this cell that holds the COntour information (SASArcs)

sas_area()
Returns:area of all computed cells in this SASA Shape that is area of the intersection of sphere with its voronoi cell
sas_volume()
Returns:volume of all computed cells in this SASA Shape that is volume of the intersection of ball with its voronoi cell
size()
Returns:number of cells in this shape, the same as size in diagram
class pydmga.shape.Shape

This is abstract Shape class, it is used to define various properties of the voronoi diagram that can be inferred using distance and neighbouring information. Derived class should implement size() and get_cell(i) functions.

Shape provides standard iterator and get_cells(subset) functions.

allows iteration in the form:

for shape_cell in shape:
    pass

NOTE This is abstract class and should not be instantiated

get_cell(number)

Not implemented yet

get_cells(subset=None)

returns list of cells for given subset subset is iterable of integers or None (default: None) if None is given then all cells will be returned this returns list so may be time and space consuming for generating use iteration:

for shape_cell in shape:
    pass

where shape_cell is of class AlphaShapeCell

size()

Not implemented yet

class pydmga.shape.ShapeIterator(shape)

To iterate over cells in a Shape. Used internally in Shape.__iter__()

NOTE This is default Shape iterator and should work for all the shapes. It is necessary for the Shape to have get_cell(i) function

next()

returns next cell shape or rise StopIteration if no cells