landlab.graph.framed_voronoi.framed_voronoi

Implementation of the FramedVoronoiGraph and its static layout: HorizontalRectVoronoiGraph. This pattern is inspired from the developments of the HexModelGrid

Code author: sebastien lenard

class FramedVoronoiGraph[source]

Bases: DelaunayGraph

VoronoiDelaunay graph based on a fixed lattice.

Graph of an unstructured grid of Voronoi Delaunay cells and irregular patches. It is a special type of :class`~.VoronoiDelaunayGraph` in which the initial set of points is arranged in a fixed lattice (e.g. like a RasterModelGrid) named here “layout” and the core points are then moved from their initial position by a random distance, lower than a certain threshold.

Examples

>>> from landlab.graph import FramedVoronoiGraph
>>> graph = FramedVoronoiGraph((3, 3), seed=200)
>>> graph.number_of_nodes
9
>>> graph.x_of_node[2:4]
array([2., 0.])
>>> graph.y_of_node[2:4]
array([0.   , 0.749])
>>> graph.y_of_node[5]
1.251
>>> graph.number_of_links
16
>>> graph.number_of_patches
8

Create the graph.

Parameters:
  • shape (tuple of int) – Number of rows and columns of nodes.

  • xy_spacing (float or tuple of float, optional) – Node spacing along x and y coordinates. If float, same spacing x and y spacing.

  • xy_of_lower_left (tuple, optional) – Minimum x-of-node and y-of-node values. Depending on the grid, a node may not be present at this location.

  • sort (bool) – If True, nodes, links and patches are re-numbered according certain their positions. Currently not used.

  • xy_min_spacing (float or tuple of float, optional) – Final minimal spacing between nodes. Random moves of the core nodes around their position cannot be above this threshold: (xy_spacing - xy_min_spacing) / 2 If float, same minimal spacing for x and y.

  • seed (int, optional) – Seed used to generate the random x and y moves. When set, controls a pseudo-randomness of moves to ensure reproducibility. When None, seed is random and the moves of coordinates are completely random.

Returns:

A newly-created graph.

Return type:

FramedVoronoiGraph

Examples

Create a grid with 3 rows and 2 columns of nodes.

>>> from landlab.graph import FramedVoronoiGraph
>>> graph = FramedVoronoiGraph((3, 2))
>>> graph.number_of_nodes
6
__init__(shape, xy_spacing=(1.0, 1.0), xy_of_lower_left=(0.0, 0.0), sort=False, xy_min_spacing=(0.5, 0.5), seed=200)[source]

Create the graph.

Parameters:
  • shape (tuple of int) – Number of rows and columns of nodes.

  • xy_spacing (float or tuple of float, optional) – Node spacing along x and y coordinates. If float, same spacing x and y spacing.

  • xy_of_lower_left (tuple, optional) – Minimum x-of-node and y-of-node values. Depending on the grid, a node may not be present at this location.

  • sort (bool) – If True, nodes, links and patches are re-numbered according certain their positions. Currently not used.

  • xy_min_spacing (float or tuple of float, optional) – Final minimal spacing between nodes. Random moves of the core nodes around their position cannot be above this threshold: (xy_spacing - xy_min_spacing) / 2 If float, same minimal spacing for x and y.

  • seed (int, optional) – Seed used to generate the random x and y moves. When set, controls a pseudo-randomness of moves to ensure reproducibility. When None, seed is random and the moves of coordinates are completely random.

Returns:

A newly-created graph.

Return type:

FramedVoronoiGraph

Examples

Create a grid with 3 rows and 2 columns of nodes.

>>> from landlab.graph import FramedVoronoiGraph
>>> graph = FramedVoronoiGraph((3, 2))
>>> graph.number_of_nodes
6
__new__(**kwargs)
property adjacent_nodes_at_node

Get adjacent nodes.

Examples

>>> from landlab.graph import Graph

First, a simple example with no diagonals.

>>> node_x, node_y = [0, 0, 0, 1, 1, 1, 2, 2, 2], [0, 1, 2, 0, 1, 2, 0, 1, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links, sort=True)
>>> graph.adjacent_nodes_at_node
array([[ 1,  3, -1, -1],
       [ 2,  4,  0, -1],
       [ 5,  1, -1, -1],
       [ 4,  6,  0, -1],
       [ 5,  7,  3,  1],
       [ 8,  4,  2, -1],
       [ 7,  3, -1, -1],
       [ 8,  6,  4, -1],
       [ 7,  5, -1, -1]])

Next, we add the diagonal from node 0 to node 4.

>>> node_x, node_y = [0, 0, 0, 1, 1, 1, 2, 2, 2], [0, 1, 2, 0, 1, 2, 0, 1, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
...     (0, 4),
... )
>>> graph = Graph((node_y, node_x), links=links, sort=True)
>>> graph.adjacent_nodes_at_node
array([[ 1,  4,  3, -1, -1],
       [ 2,  4,  0, -1, -1],
       [ 5,  1, -1, -1, -1],
       [ 4,  6,  0, -1, -1],
       [ 5,  7,  3,  0,  1],
       [ 8,  4,  2, -1, -1],
       [ 7,  3, -1, -1, -1],
       [ 8,  6,  4, -1, -1],
       [ 7,  5, -1, -1, -1]])

Get the angle of each link.

Examples

>>> import numpy as np
>>> from landlab.graph import Graph
>>> node_x, node_y = ([0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1])
>>> links = ((0, 1), (1, 2), (0, 3), (1, 4), (2, 5), (3, 4), (4, 5))
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.angle_of_link * 180.0 / np.pi
array([ 0.,  0., 90., 90., 90.,  0.,  0.])
property area_of_patch

Get the area of each patch.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches)
>>> graph.area_of_patch
array([1.,  1.])
property ds
freeze()

Freeze the graph by making arrays read-only.

classmethod from_dict(meta)
classmethod from_netcdf(fname)
property frozen

Get the length of links.

Examples

>>> import numpy as np
>>> from landlab.graph import UniformRectilinearGraph
>>> graph = UniformRectilinearGraph((2, 3), spacing=(1, 2))
>>> graph.length_of_link
array([2., 2., 1., 1., 1., 2., 2.])

Return link directions into each node.

A value of 1 indicates a link points toward a given node, while a value of -1 indicates a link points away from a node.

Returns:

Link directions relative to the nodes of a grid. The shape of the matrix will be number of nodes by the maximum number of links per node. A zero indicates no link at this position.

Return type:

(n_nodes, max_links_per_node) ndarray of int

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.link_dirs_at_node
array([[-1, -1,  0,  0], [-1, -1,  1,  0], [-1,  1,  0,  0],
       [-1, -1,  1,  0], [-1, -1,  1,  1], [-1,  1,  1,  0],
       [-1,  1,  0,  0], [-1,  1,  1,  0], [ 1,  1,  0,  0]],
      dtype=int8)

Get links touching a node.

Examples

>>> from landlab.graph import Graph
>>> node_x = [0, 1, 2, 0, 1, 2, 0, 1, 2]
>>> node_y = [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.links_at_node
array([[ 0,  2, -1, -1], [ 1,  3,  0, -1], [ 4,  1, -1, -1],
       [ 5,  7,  2, -1], [ 6,  8,  5,  3], [ 9,  6,  4, -1],
       [10,  7, -1, -1], [11, 10,  8, -1], [11,  9, -1, -1]])

Get the links that define a patch.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches, sort=True)
>>> graph.links_at_patch
array([[3, 5, 2, 0],
       [4, 6, 3, 1]])
classmethod load(source)
merge(dual, node_at_cell=None, nodes_at_face=None)

Get the middle of links.

Examples

>>> import numpy as np
>>> from landlab.graph import UniformRectilinearGraph
>>> graph = UniformRectilinearGraph((2, 3), spacing=(1, 2))
>>> graph.midpoint_of_link
array([[1. , 0. ], [3. , 0. ],
       [0. , 0.5], [2. , 0.5], [4. , 0.5],
       [1. , 1. ], [3. , 1. ]])
property ndim

Get nodes at link head.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.node_at_link_head
array([1, 2, 3, 4, 5, 4, 5, 6, 7, 8, 7, 8])

Get nodes at link tail.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.node_at_link_tail
array([0, 1, 0, 1, 2, 3, 4, 3, 4, 5, 6, 7])
property node_layout
property node_x
property node_y
property nodes

Get identifier for each node.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1]
>>> graph = Graph((node_y, node_x))
>>> graph.nodes
array([0, 1, 2, 3, 4, 5])
property nodes_at_bottom_edge
property nodes_at_left_edge

Get nodes at either end of links.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.nodes_at_link
array([[0, 1], [1, 2],
       [0, 3], [1, 4], [2, 5],
       [3, 4], [4, 5],
       [3, 6], [4, 7], [5, 8],
       [6, 7], [7, 8]])
property nodes_at_patch

Get the nodes that define a patch.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = ([0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2])
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches, sort=True)
>>> graph.nodes_at_patch
array([[4, 3, 0, 1],
       [5, 4, 1, 2]])
property nodes_at_right_edge
property nodes_at_top_edge

Get nodes at link head.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> graph = Graph((node_y, node_x), links=links)
>>> graph.number_of_links == 12
True
property number_of_nodes

Get total number of nodes.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1]
>>> graph = Graph((node_y, node_x))
>>> graph.number_of_nodes == 6
True
property number_of_patches

Get the number of patches.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches)
>>> graph.number_of_patches == 2
True
property orientation

Get the patches on either side of each link.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = ([0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1])
>>> links = ((0, 1), (1, 2), (0, 3), (1, 4), (2, 5), (3, 4), (4, 5))
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches)
>>> graph.patches_at_link
array([[ 0, -1], [ 1, -1],
       [ 0, -1], [ 0,  1], [ 1, -1],
       [ 0, -1], [ 1, -1]])
property patches_at_node

Get the patches that touch each node.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = ([0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1])
>>> links = ((0, 1), (1, 2), (0, 3), (1, 4), (2, 5), (3, 4), (4, 5))
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches, sort=True)
>>> graph.patches_at_node
array([[ 0, -1], [ 1,  0], [ 1, -1],
       [ 0, -1], [ 0,  1], [ 1, -1]])
property perimeter_nodes

Get nodes on the convex hull of a Graph.

Examples

>>> import numpy as np
>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1]
>>> graph = Graph((node_y, node_x))
>>> np.sort(graph.perimeter_nodes)
array([0, 2, 3, 5])
property shape
sort()

Sort graph elements.

thaw()

Thaw the graph by making arrays writable.

thawed()
to_dict()
to_json()
to_netcdf(*args, **kwds)

Write graph contents to a netCDF file.

See xarray.Dataset.to_netcdf for a complete list of parameters. Below are only the most common.

Parameters:
  • path (str, optional) – Path to which to save this graph.

  • mode ({'w', 'a'}, optional) – Write (‘w’) or append (‘a’) mode. If mode=’w’, any existing file at this location will be overwritten.

  • format ({'NETCDF4', 'NETCDF4_CLASSIC', 'NETCDF3_64BIT', 'NETCDF3_CLASSIC'}, optional) –

    File format for the resulting netCDF file:

    • NETCDF4: Data is stored in an HDF5 file, using netCDF4 API features.

    • NETCDF4_CLASSIC: Data is stored in an HDF5 file, using only netCDF 3 compatible API features.

    • NETCDF3_64BIT: 64-bit offset version of the netCDF 3 file format, which fully supports 2+ GB files, but is only compatible with clients linked against netCDF version 3.6.0 or later.

    • NETCDF3_CLASSIC: The classic netCDF 3 file format. It does not handle 2+ GB files very well.

    All formats are supported by the netCDF4-python library. scipy.io.netcdf only supports the last two formats.

    The default format is NETCDF4 if you are saving a file to disk and have the netCDF4-python library available. Otherwise, xarray falls back to using scipy to write netCDF files and defaults to the NETCDF3_64BIT format (scipy does not support netCDF4).

Make arrays to store the unit vectors associated with each link.

For each link, the x and y components of the link’s unit vector (that is, the link’s x and y dimensions if it were shrunk to unit length but retained its orientation).

Examples

The example below is a seven-node hexagonal grid, with six nodes around the perimeter and one node (#3) in the interior. There are four horizontal links with unit vector (1,0), and 8 diagonal links with unit vector (+/-0.5, +/-sqrt(3)/2) (note: sqrt(3)/2 ~ 0.866).

>>> from landlab.graph import TriGraph
>>> graph = TriGraph((3, 2), spacing=2.0, node_layout="hex", sort=True)
>>> np.round(graph.unit_vector_at_link[:, 0], decimals=5)
array([ 1. , -0.5,  0.5, -0.5,  0.5,  1. ,  1. ,  0.5, -0.5,  0.5, -0.5,
        1. ])
>>> np.round(graph.unit_vector_at_link[:, 1], decimals=5)
array([0.     , 0.86603, 0.86603, 0.86603, 0.86603, 0.     ,
       0.     , 0.86603, 0.86603, 0.86603, 0.86603, 0.     ])
property unit_vector_at_node

Get a unit vector for each node.

Examples

>>> from landlab.graph import UniformRectilinearGraph
>>> graph = UniformRectilinearGraph((3, 3))
>>> graph.unit_vector_at_node
array([[1., 1.],
       [2., 1.],
       [1., 1.],
       [1., 2.],
       [2., 2.],
       [1., 2.],
       [1., 1.],
       [2., 1.],
       [1., 1.]])
>>> from landlab.graph import TriGraph
>>> graph = TriGraph((3, 2), spacing=2.0, node_layout="hex", sort=True)
>>> unit_vector_at_node = np.round(graph.unit_vector_at_node, decimals=5)
>>> unit_vector_at_node[:, 0]
array([2., 2., 2., 4., 2., 2., 2.])
>>> unit_vector_at_node[:, 1]
array([1.73205, 1.73205, 1.73205, 3.4641 , 1.73205, 1.73205, 1.73205])
property x_of_node

Get x-coordinate of node.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1]
>>> graph = Graph((node_y, node_x))
>>> graph.x_of_node
array([0., 1., 2., 0., 1., 2.])
property xy_of_node

Get x and y-coordinates of node.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1]
>>> graph = Graph((node_y, node_x))
>>> graph.xy_of_node[:, 0]
array([0., 1., 2., 0., 1., 2.])
>>> graph.xy_of_node[:, 1]
array([0., 0., 0., 1., 1., 1.])
property xy_of_patch

Get the centroid of each patch.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1, 2, 2, 2]
>>> links = (
...     (0, 1),
...     (1, 2),
...     (0, 3),
...     (1, 4),
...     (2, 5),
...     (3, 4),
...     (4, 5),
...     (3, 6),
...     (4, 7),
...     (5, 8),
...     (6, 7),
...     (7, 8),
... )
>>> patches = ((0, 3, 5, 2), (1, 4, 6, 3))
>>> graph = Graph((node_y, node_x), links=links, patches=patches)
>>> graph.xy_of_patch
array([[0.5, 0.5],
       [1.5, 0.5]])
property xy_spacing
property y_of_node

Get y-coordinate of node.

Examples

>>> from landlab.graph import Graph
>>> node_x, node_y = [0, 1, 2, 0, 1, 2], [0, 0, 0, 1, 1, 1]
>>> graph = Graph((node_y, node_x))
>>> graph.y_of_node
array([0., 0., 0., 1., 1., 1.])
class HorizontalRectVoronoiGraph[source]

Bases: object

The horizontal rectangular frame for the FramedVoronoiGraph.

__init__()
__new__(**kwargs)
static corner_nodes(shape)[source]
Parameters:

shape (tuple of int) – Number of rows and number of columns

Returns:

Ids of the corner nodes

Return type:

ndarray of int

Examples

>>> from landlab.graph.framed_voronoi.framed_voronoi import (
...     HorizontalRectVoronoiGraph,
... )
>>> HorizontalRectVoronoiGraph.corner_nodes((3, 4))
(11, 8, 0, 3)
static nodes_at_edge(shape)[source]
Parameters:

shape (tuple of int) – Number of rows and number of columns

Returns:

right, top, left, bottom – For each edge give the ids of the nodes present at the edge

Return type:

ndarray of int

Examples

>>> from landlab.graph.framed_voronoi.framed_voronoi import (
...     HorizontalRectVoronoiGraph,
... )
>>> HorizontalRectVoronoiGraph.nodes_at_edge((3, 3))
(array([2, 5]), array([8, 7]), array([6, 3]), array([0, 1]))
static number_of_nodes(shape)[source]
Parameters:

shape (tuple of int) – Number of rows and number of columns

Returns:

Number of nodes

Return type:

int

Examples

>>> from landlab.graph.framed_voronoi.framed_voronoi import (
...     HorizontalRectVoronoiGraph,
... )
>>> HorizontalRectVoronoiGraph.number_of_nodes((3, 2))
6
static number_of_perimeter_nodes(shape)[source]
Parameters:

shape (tuple of int) – Number of rows and number of columns

Returns:

Number of perimeter nodes

Return type:

int

Examples

>>> from landlab.graph.framed_voronoi.framed_voronoi import (
...     HorizontalRectVoronoiGraph,
... )
>>> HorizontalRectVoronoiGraph.number_of_perimeter_nodes((3, 4))
10
static perimeter_nodes(shape)[source]
Parameters:

shape (tuple of int) – Number of rows and number of columns

Returns:

Ids of the perimeter nodes

Return type:

ndarray of int

Examples

>>> from landlab.graph.framed_voronoi.framed_voronoi import (
...     HorizontalRectVoronoiGraph,
... )
>>> HorizontalRectVoronoiGraph.perimeter_nodes((3, 3))
array([2, 5, 8, 7, 6, 3, 0, 1])
static xy_of_node(shape, xy_spacing=(1.0, 1.0), xy_of_lower_left=(0.0, 0.0), xy_min_spacing=(0.5, 0.5), seed=200)[source]

The x and y coordinates of the graph’s nodes.

Calculation of the x-y coordinates is done following these steps:

  1. Generate a rectangular, regular meshgrid.

  2. Move the coordinates of the core nodes over a random distance around their initial position, within a threshold calculated from xy_spacing and xy_min_spacing.

  3. Rectify the y-coordinates of the nodes of the left and right to ensure that the leftmost node of a row has a lower y than the rightmost node. This ensures that the ids of these nodes are not modified by subsequent sorting operations on the graph and make it possible to get the perimeter nodes in simple way.

Parameters:
  • shape (tuple of int) – Number of rows and columns of nodes.

  • xy_spacing (float or tuple of float, optional) – Node spacing along x and y coordinates. If float, same spacing at x and y.

  • xy_of_lower_left (tuple, optional) – Minimum x-of-node and y-of-node values. Depending on the grid, a node may not be present at this location.

  • xy_min_spacing (float or tuple of float, optional) – Final minimal spacing between nodes. Random moves of the core nodes around their initial positions cannot be above this threshold: (xy_spacing - xy_min_spacing) / 2. If float, same minimal spacing for x and y.

  • seed (int, optional) – Seed used to generate the random x and y moves. When set, controls a pseudo-randomness of moves to ensure reproducibility. When None, seed is random and the moves of coordinates are completely random.

Returns:

x_of_node, y_of_node – The arrays of x and y coordinates.

Return type:

ndarray of float

Examples

>>> from landlab.graph.framed_voronoi.framed_voronoi import (
...     HorizontalRectVoronoiGraph,
... )
>>> x_of_node, y_of_node = HorizontalRectVoronoiGraph.xy_of_node(
...     (3, 3), seed=200
... )

Coordinates of the lower left node,

>>> x_of_node[0], y_of_node[0]
(0.0, 0.0)

x coordinates of the left and right edges,

>>> x_of_node[3], x_of_node[5]
(0.0, 2.0)

y coordinate of the middle row of the left edge,

>>> y_of_node[3]
0.749