DelaunayVoronoiDual

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

Bases: openpnm.network.GenericNetwork.GenericNetwork

Combined and interconnected Voronoi and Delaunay tessellations

A Delaunay tessellation is performed on a set of base points then the corresponding Voronoi diagram is generated. Finally, each Delaunay node is connected to it’s neighboring Voronoi vertices to create interaction between the two networks.

All pores and throats are labelled according to their network (i.e. ‘pore.delaunay’), so they can be each assigned to a different Geometry.

The dual-nature of this network is meant for modeling transport in the void and solid space simultaneously by treating one network (i.e. Delaunay) as voids and the other (i.e. Voronoi) as solid. Interaction such as heat transfer between the solid and void can be accomplished via the interconnections between the Delaunay and Voronoi nodes.

Parameters
  • points (array_like or int) – Can either be an N-by-3 array of point coordinates which will be used, or a scalar value indicating the number of points to generate

  • shape (array_like) –

    The size and shape of the domain used for generating and trimming excess points. The coordinates are treated as the outer corner of a rectangle [x, y, z] whose opposite corner lies at [0, 0, 0].

    By default, a domain size of [1, 1, 1] is used. To create a 2D network set the Z-dimension to 0.

  • name (string) – An optional name for the object to help identify it. If not given, one will be generated.

  • project (OpenPNM Project object, optional) – Each OpenPNM object must be part of a Project. If none is supplied then one will be created and this Network will be automatically assigned to it. To create a Project use openpnm.Project().

Examples

Points will be automatically generated if none are given:

>>> import openpnm as op
>>> net = op.network.DelaunayVoronoiDual(points=50, shape=[1, 1, 0])

The resulting network can be quickly visualized using opnepnm.topotools.plot_connections.

add_boundary_pores(labels=['top', 'bottom', 'front', 'back', 'left', 'right'], offset=None)[source]

Add boundary pores to the specified faces of the network

Pores are offset from the faces of the domain.

Parameters
  • labels (string or list of strings) – The labels indicating the pores defining each face where boundary pores are to be added (e.g. ‘left’ or [‘left’, ‘right’])

  • offset (scalar or array_like) – The spacing of the network (e.g. [1, 1, 1]). This must be given since it can be quite difficult to infer from the network, for instance if boundary pores have already added to other faces.

find_pore_hulls(pores=None)[source]

Finds the indices of the Voronoi nodes that define the convex hull around the given Delaunay nodes.

Parameters

pores (array_like) – The pores whose convex hull are sought. The given pores should be from the ‘delaunay’ network. If no pores are given, then the hull is found for all ‘delaunay’ pores.

Notes

This metod is not fully optimized as it scans through each pore in a for-loop, so could be slow for large networks.

find_throat_facets(throats=None)[source]

Finds the indicies of the Voronoi nodes that define the facet or ridge between the Delaunay nodes connected by the given throat.

Parameters

throats (array_like) – The throats whose facets are sought. The given throats should be from the ‘delaunay’ network. If no throats are specified, all ‘delaunay’ throats are assumed.

Notes

The method is not well optimized as it scans through each given throat inside a for-loop, so it could be slow for large networks.