API for landlab.grid.network#

A class used to create and manage network models in 2D.

class NetworkModelGrid(yx_of_node, links, xy_axis_name=('x', 'y'), xy_axis_units='-', xy_of_reference=(0.0, 0.0))[source]#

Bases: NetworkGraph, GraphFields

Create a ModelGrid of just nodes and links.

Parameters:
  • yx_of_node (tuple of ndarray) – Node y and x coordinates.

  • links (array of tuple of int) – Nodes at link tail and head.

  • xy_of_reference (tuple, optional) – Coordinate value in projected space of (0., 0.) Default is (0., 0.)

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.x_of_node
array([ 0.,  0., -1.,  1.])
>>> grid.y_of_node
array([ 0.,  1.,  2.,  2.])
>>> grid.nodes_at_link
array([[0, 1],
       [2, 1],
       [1, 3]])

Define a graph of connected nodes.

Parameters:

mesh (Dataset) – xarray Dataset that defines the topology in ugrid format.

BAD_INDEX = -1#

Indicates a node is bad index.

Indicates a link is active, and can carry flux

Indicates a link has a fixed gradient value, and behaves as a boundary

Indicates a link is inactive, and cannot carry flux

BC_NODE_IS_CLOSED = 4#

Indicates a boundary node is closed

BC_NODE_IS_CORE = 0#

Indicates a node is core.

BC_NODE_IS_FIXED_GRADIENT = 2#

Indicates a boundary node has a fixed gradient.

BC_NODE_IS_FIXED_VALUE = 1#

Indicates a boundary node has a fixed value.

BC_NODE_IS_LOOPED = 3#

Indicates a boundary node is wrap-around.

VALID_LOCATIONS = ('node', 'link', 'grid')#

Grid elements on which fields can be placed.

__getitem__(name)#

Get the collection of fields from the named group.

Get array of active links.

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.active_links
array([0, 1, 2])
add_empty(name, at='node', units='-', clobber=False)#

Create and add an uninitialized array of values to the field.

Create a new array of the data field size, without initializing entries, and add it to the field as name. The units keyword gives the units of the new fields as a string. Remaining keyword arguments are the same as that for the equivalent numpy function.

This method is not valid for the group grid.

Parameters:
  • name (str) – Name of the new field to add.

  • at (str, optional) – Grid location to store values. If not given, values are assumed to be on node.

  • units (str, optional) – Optionally specify the units of the field.

  • clobber (bool, optional) – Raise an exception if adding to an already existing field.

Returns:

A reference to the newly-created array.

Return type:

numpy.ndarray

See also

numpy.empty

See for a description of optional keywords.

empty

Equivalent method that does not initialize the new array.

zeros

Equivalent method that initializes the data to 0.

add_field(name, value_array, at='node', units='-', copy=False, clobber=False)#

Add an array of values to the field.

Add an array of data values to a collection of fields and associate it with the key, name. Use the copy keyword to, optionally, add a copy of the provided array.

In the case of adding to the collection grid, the added field is a numpy scalar rather than a numpy array.

Parameters:
  • name (str) – Name of the new field to add.

  • value_array (numpy.array) – Array of values to add to the field.

  • at (str, optional) – Grid location to store values. If not given, values are assumed to be on node.

  • units (str, optional) – Optionally specify the units of the field.

  • copy (bool, optional) – If True, add a copy of the array to the field. Otherwise save add a reference to the array.

  • clobber (bool, optional) – Raise an exception if adding to an already existing field.

Returns:

The data array added to the field. Depending on the copy keyword, this could be a copy of value_array or value_array itself.

Return type:

numpy.ndarray

Raises:

ValueError – If value_array has a size different from the field.

Examples

>>> import numpy as np
>>> from landlab.field import GraphFields
>>> field = GraphFields()
>>> field.new_field_location("node", 4)
>>> values = np.ones(4, dtype=int)
>>> field.add_field("topographic__elevation", values, at="node")
array([1, 1, 1, 1])

A new field is added to the collection of fields. The saved value array is the same as the one initially created.

>>> field.at_node["topographic__elevation"] is values
True

If you want to save a copy of the array, use the copy keyword. In addition, adding values to an existing field will remove the reference to the previously saved array. The clobber=False keyword changes this behavior to raise an exception in such a case.

>>> field.add_field(
...     "topographic__elevation", values, at="node", copy=True, clobber=True
... )
array([1, 1, 1, 1])
>>> field.at_node["topographic__elevation"] is values
False
>>> field.add_field("topographic__elevation", values, at="node", clobber=False)
Traceback (most recent call last):
...
FieldError: topographic__elevation
add_full(*args, **kwds)#

Create and add an array of values, fill with fill_value.

Parameters:
  • name (str) – Name of the new field to add.

  • fill_value (scalar) – Fill value.

  • at (str, optional) – Grid location to store values. If not given, values are assumed to be on node.

  • units (str, optional) – Optionally specify the units of the field.

  • copy (bool, optional) – If True, add a copy of the array to the field. Otherwise save add a reference to the array.

  • clobber (bool, optional) – Raise an exception if adding to an already existing field.

Returns:

A reference to the newly-created array.

Return type:

numpy.ndarray

add_ones(name, at='node', units='-', clobber=False)#

Create and add an array of values, initialized to 1, to the field.

Create a new array of the data field size, filled with ones, and add it to the field as name. The units keyword gives the units of the new fields as a string. Remaining keyword arguments are the same as that for the equivalent numpy function.

This method is not valid for the group grid.

Parameters:
  • name (str) – Name of the new field to add.

  • at (str, optional) – Grid location to store values. If not given, values are assumed to be on node.

  • units (str, optional) – Optionally specify the units of the field.

  • clobber (bool, optional) – Raise an exception if adding to an already existing field.

Returns:

A reference to the newly-created array.

Return type:

numpy.ndarray

See also

numpy.ones

See for a description of optional keywords.

add_empty

Equivalent method that does not initialize the new array.

add_zeros

Equivalent method that initializes the data to 0.

Examples

Add a new, named field to a collection of fields.

>>> from landlab.field import GraphFields
>>> field = GraphFields()
>>> field.new_field_location("node", 4)
>>> field.add_ones("topographic__elevation", at="node")
array([ 1.,  1.,  1.,  1.])
>>> list(field.keys("node"))
['topographic__elevation']
>>> field["node"]["topographic__elevation"]
array([ 1.,  1.,  1.,  1.])
>>> field.at_node["topographic__elevation"]
array([ 1.,  1.,  1.,  1.])
add_zeros(name, at='node', units='-', clobber=False)#

Create and add an array of values, initialized to 0, to the field.

Create a new array of the data field size, filled with zeros, and add it to the field as name. The units keyword gives the units of the new fields as a string. Remaining keyword arguments are the same as that for the equivalent numpy function.

Parameters:
  • name (str) – Name of the new field to add.

  • at (str, optional) – Grid location to store values. If not given, values are assumed to be on node.

  • units (str, optional) – Optionally specify the units of the field.

  • clobber (bool, optional) – Raise an exception if adding to an already existing field.

Returns:

A reference to the newly-created array.

Return type:

array

See also

numpy.zeros

See for a description of optional keywords.

add_empty

Equivalent method that does not initialize the new array.

add_ones

Equivalent method that initializes the data to 1.

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.])
at_grid = {}#
at_node = {}#
property axis_name#

Get the name of each coordinate axis.

Returns:

The names of each axis.

Return type:

tuple of str

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.axis_name
('x', 'y')
>>> grid.axis_name = ("lon", "lat")
>>> grid.axis_name
('lon', 'lat')
>>> grid = NetworkModelGrid(
...     (y_of_node, x_of_node), nodes_at_link, xy_axis_name=("lon", "lat")
... )
>>> grid.axis_name
('lon', 'lat')
property axis_units#

Get units for each axis.

Returns:

The units (as a string) for each of a grid’s coordinates.

Return type:

tuple of str

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.axis_units
('-', '-')
>>> grid.axis_units = ("km", "km")
>>> grid.axis_units
('km', 'km')

Calculate gradients of node values at links.

Calculates the gradient in node_values at each link in the grid, returning an array of length number_of_links.

Parameters:
  • node_values (ndarray or field name (x number of nodes)) – Values at grid nodes.

  • out (ndarray, optional (x number of links)) – Buffer to hold the result.

Returns:

Gradients across active links.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> rg = RasterModelGrid((3, 4), xy_spacing=10.0)
>>> z = rg.add_zeros("topographic__elevation", at="node")
>>> z[5] = 50.0
>>> z[6] = 36.0
>>> calc_grad_at_link(rg, z)  # there are 17 links
array([ 0. ,  0. ,  0. ,  0. ,  5. ,  3.6,  0. ,  5. , -1.4, -3.6,  0. ,
       -5. , -3.6,  0. ,  0. ,  0. ,  0. ])
>>> from landlab import HexModelGrid
>>> hg = HexModelGrid((3, 3), spacing=10.0)
>>> z = hg.add_zeros("topographic__elevation", at="node", clobber=True)
>>> z[4] = 50.0
>>> z[5] = 36.0
>>> calc_grad_at_link(hg, z)  # there are 11 faces
array([ 0. ,  0. ,  0. ,  5. ,  5. ,  3.6,  3.6,  0. ,  5. , -1.4, -3.6,
        0. , -5. , -5. , -3.6, -3.6,  0. ,  0. ,  0. ])
property default_group#

Return the name of the group into which fields are put by default.

delete_field(loc, name)#

Erases an existing field.

Parameters:
  • loc (str) – Name of the group.

  • name (str) – Name of the field.

Raises:

KeyError – If the named field does not exist.

property ds#
empty(*args, **kwds)#

Uninitialized array whose size is that of the field.

Return a new array of the data field size, without initializing entries. Keyword arguments are the same as that for the equivalent numpy function.

Parameters:

group (str) – Name of the group.

See also

numpy.empty

See for a description of optional keywords.

ones

Equivalent method that initializes the data to 1.

zeros

Equivalent method that initializes the data to 0.

Examples

>>> from landlab.field import GraphFields
>>> field = GraphFields()
>>> field.new_field_location("node", 4)
>>> field.empty("node")  
array([  2.31584178e+077,  -2.68156175e+154,   9.88131292e-324,
... 2.78134232e-309]) # Uninitialized memory

Note that a new field is not added to the collection of fields.

>>> list(field.keys("node"))
[]
field_units(field, at=None)#

Get units for a field.

Returns the unit string associated with the data array in group and field.

Parameters:
  • field (str) – Name of the field withing group.

  • at (str, optional) – Name of the group.

Returns:

The units of the field.

Return type:

str

Raises:

KeyError – If either field or group does not exist.

field_values(field, at=None)#

Return the values of a field.

Given a group and a field, return a reference to the associated data array.

Parameters:
  • field (str) – Name of the field within group.

  • at (str, optional) – Name of the group.

Returns:

The values of the field.

Return type:

array

Raises:
  • landlab.field.errors.GroupError – If group does not exist

  • landlab.field.errors.FieldError – If field does not exist

Examples

Create a group of fields called node.

>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 4)

Add a field, initialized to ones, called topographic__elevation to the node group. The field_values method returns a reference to the field’s data.

>>> _ = fields.add_ones("topographic__elevation", at="node")
>>> fields.field_values("topographic__elevation", at="node")
array([ 1.,  1.,  1.,  1.])

Raise FieldError if field does not exist in group.

>>> fields.field_values("planet_surface__temperature", at="node")
Traceback (most recent call last):
...
FieldError: planet_surface__temperature

If group does not exists, raise GroupError.

>>> fields.field_values("topographic__elevation", at="cell")
Traceback (most recent call last):
...
GroupError: cell
freeze()#

Freeze the graph by making arrays read-only.

classmethod from_dict(params)[source]#
classmethod from_file(file_like)[source]#
classmethod from_netcdf(fname)#
property frozen#
property groups#

List of group names.

Returns:

Names of field groupings.

Return type:

set

has_field(field, at=None)#

Check if a field is in a group.

Parameters:
  • field (str) – Name of the field.

  • at (str, optional) – Name of the group.

Returns:

True if the group contains the field, otherwise False.

Return type:

bool

Examples

Check if the field named topographic__elevation is contained in a group.

>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 12)
>>> _ = fields.add_ones("topographic__elevation", at="node")
>>> fields.has_field("topographic__elevation", at="node")
True
>>> fields.has_field("topographic__elevation", at="cell")
False
>>> fields = GraphFields()
>>> fields.new_field_location("node", 12)
>>> _ = fields.add_ones("topographic__elevation", at="node")
>>> fields.has_field("node", "topographic__elevation")
True
>>> fields.has_field("cell", "topographic__elevation")
False
has_group(name)#

Check if a group exists.

Parameters:

name (str) – Name of the group.

Returns:

True if the field contains group, otherwise False.

Return type:

bool

Examples

Check if the field has the groups named node or cell.

>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 12)
>>> fields.has_group("node")
True
>>> fields.has_group("cell")
False
keys(group)#

Return the field names in a group.

Parameters:

group (str) – Group name.

Returns:

Names of fields held in the given group.

Return type:

list

Examples

>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 4)
>>> list(fields.keys("node"))
[]
>>> _ = fields.add_empty("topographic__elevation", at="node")
>>> list(fields.keys("node"))
['topographic__elevation']

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]])
classmethod load(source)#

Map the largest magnitude of the links carrying flux from the node to the node.

map_downwind_node_link_max_to_node iterates across the grid and identifies the link values at each link connected to a node. It then uses the link_dirs_at_node data structure to identify links carrying flux out of the node, then maps the maximum magnitude of ‘var_name’ found on these links onto the node. If no downwind link is found, the value will be recorded as zero.

Parameters:
  • var_name (array or field name) – Values defined at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_downwind_node_link_max_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = [
...     -1.0,
...     -2.0,
...     -1.0,
...     0.0,
...     0.0,
...     0.0,
...     0.0,
...     -1.0,
...     -2.0,
...     -1.0,
...     0.0,
...     0.0,
...     0.0,
...     0.0,
...     -1.0,
...     -2.0,
...     -1.0,
... ]
>>> map_downwind_node_link_max_to_node(rmg, "grad")
array([ 1.,  2.,  1.,  0.,
        1.,  2.,  1.,  0.,
        1.,  2.,  1.,  0.])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_downwind_node_link_max_to_node(rmg, "grad", out=values_at_nodes)
>>> values_at_nodes
array([ 1.,  2.,  1.,  0.,
        1.,  2.,  1.,  0.,
        1.,  2.,  1.,  0.])
>>> rtn is values_at_nodes
True

Map the mean magnitude of the links carrying flux out of the node to the node.

map_downwind_node_link_mean_to_node iterates across the grid and identifies the link values at each link connected to a node. It then uses the link_dirs_at_node data structure to identify links carrying flux out of the node, then maps the mean magnitude of ‘var_name’ found on these links onto the node. Links with zero values are not included in the means, and zeros are returned if no upwind links are found.

Parameters:
  • var_name (array or field name) – Values defined at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_downwind_node_link_mean_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = [
...     -1.0,
...     -2.0,
...     -1.0,
...     -2.0,
...     -3.0,
...     -4.0,
...     -5.0,
...     -1.0,
...     -2.0,
...     -1.0,
...     -1.0,
...     -2.0,
...     -3.0,
...     -4.0,
...     -1.0,
...     -2.0,
...     -1.0,
... ]
>>> map_downwind_node_link_mean_to_node(rmg, "grad")
array([ 1.5,  2.5,  2.5,  5. ,
        1. ,  2. ,  2. ,  4. ,
        1. ,  2. ,  1. ,  0. ])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_downwind_node_link_mean_to_node(rmg, "grad", out=values_at_nodes)
>>> values_at_nodes
array([ 1.5,  2.5,  2.5,  5. ,
        1. ,  2. ,  2. ,  4. ,
        1. ,  2. ,  1. ,  0. ])
>>> rtn is values_at_nodes
True

Map values from a link head nodes to links.

Iterate over a grid and identify the node at the head. For each link, the value of var_name at the head node is mapped to the corresponding link.

In a RasterModelGrid, each one node has two adjacent “link heads”. This means each node value is mapped to two corresponding links.

Parameters:
  • var_name (array or field name) – Values defined at nodes.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_link_head_node_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_node["z"] = np.array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11])
>>> map_link_head_node_to_link(rmg, "z")
array([  1.,   2.,   3.,   4.,   5.,   6.,   7.,   5.,   6.,   7.,   8.,
      9.,  10.,  11.,   9.,  10.,  11.])
>>> values_at_links = rmg.empty(at="link")
>>> rtn = map_link_head_node_to_link(rmg, "z", out=values_at_links)
>>> values_at_links
array([  1.,   2.,   3.,   4.,   5.,   6.,   7.,   5.,   6.,   7.,   8.,
      9.,  10.,  11.,   9.,  10.,  11.])
>>> rtn is values_at_links
True

Map values from a link tail nodes to links.

map_link_tail_node_to_link iterates across the grid and identifies the node at the “tail”, or the “from” node for each link. For each link, the value of ‘var_name’ at the “from” node is mapped to the corresponding link.

In a RasterModelGrid, each one node has two adjacent “link tails”. This means each node value is mapped to two corresponding links.

Parameters:
  • var_name (array or field name) – Values defined at nodes.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_link_tail_node_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_node["z"] = np.array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11])
>>> map_link_tail_node_to_link(rmg, "z")
array([  0.,   1.,   2.,   0.,   1.,   2.,   3.,   4.,   5.,   6.,   4.,
         5.,   6.,   7.,   8.,   9.,  10.])
>>> values_at_links = rmg.empty(at="link")
>>> rtn = map_link_tail_node_to_link(rmg, "z", out=values_at_links)
>>> values_at_links
array([  0.,   1.,   2.,   0.,   1.,   2.,   3.,   4.,   5.,   6.,   4.,
         5.,   6.,   7.,   8.,   9.,  10.])
>>> rtn is values_at_links
True

Map (x,y) components of link data data_at_link onto nodes.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid, HexModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> link_data = np.arange(grid.number_of_links)
>>> vx, vy = grid.map_link_vector_components_to_node(link_data)
>>> vx[5:7]
array([ 7.5, 8.5])
>>> grid = HexModelGrid((3, 3))
>>> link_data = np.zeros(grid.number_of_links) + 0.5 * 3.0**0.5
>>> link_data[np.isclose(grid.angle_of_link, 0.0)] = 0.0
>>> vx, vy = grid.map_link_vector_components_to_node(link_data)
>>> vy
array([ 0.,  0.,  0.,  0.,  1.,  1.,  0.,  0.,  0.,  0.])

Map the maximum of a link’s nodes to the link.

map_max_of_link_nodes_to_link iterates across the grid and identifies the node values at both the “head” and “tail” of a given link. This function evaluates the value of ‘var_name’ at both the “to” and “from” node. The maximum value of the two node values is then mapped to the link.

Parameters:
  • var_name (array or field name) – Values defined at nodes.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_max_of_link_nodes_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field(
...     "z",
...     [
...         [0, 1, 2, 3],
...         [7, 6, 5, 4],
...         [8, 9, 10, 11],
...     ],
...     at="node",
... )
>>> map_max_of_link_nodes_to_link(rmg, "z")
array([  1.,   2.,   3.,   7.,   6.,   5.,   4.,   7.,   6.,   5.,   8.,
         9.,  10.,  11.,   9.,  10.,  11.])
>>> values_at_links = rmg.empty(at="link")
>>> rtn = map_max_of_link_nodes_to_link(rmg, "z", out=values_at_links)
>>> values_at_links
array([  1.,   2.,   3.,   7.,   6.,   5.,   4.,   7.,   6.,   5.,   8.,
         9.,  10.,  11.,   9.,  10.,  11.])
>>> rtn is values_at_links
True

Map the maximum value of a nodes’ links to the node.

map_max_of_node_links_to_node iterates across the grid and identifies the link values at each link connected to a node. This function finds the maximum value of ‘var_name’ of each set of links, and then maps this value to the node. Note no attempt is made to honor the directionality of the links.

Parameters:
  • var_name (array or field name) – Values defined at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_max_of_node_links_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = np.arange(rmg.number_of_links)
>>> map_max_of_node_links_to_node(rmg, "grad")
array([  3.,   4.,   5.,   6.,
        10.,  11.,  12.,  13.,
        14.,  15.,  16.,  16.])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_max_of_node_links_to_node(rmg, "grad", out=values_at_nodes)
>>> values_at_nodes
array([  3.,   4.,   5.,   6.,
        10.,  11.,  12.,  13.,
        14.,  15.,  16.,  16.])
>>> rtn is values_at_nodes
True

Map the mean of a link’s nodes to the link.

map_mean_of_link_nodes_to_link iterates across the grid and identifies the node values at both the “head” and “tail” of a given link. This function takes the sum of the two values of ‘var_name’ at both the “to” and “from” node. The average value of the two node values of ‘var_name’ is then mapped to the link.

Parameters:
  • var_name (array or field name) – Values defined at nodes.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_mean_of_link_nodes_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_node["z"] = np.array([0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11])
>>> map_mean_of_link_nodes_to_link(rmg, "z")
array([  0.5,   1.5,   2.5,   2. ,   3. ,   4. ,   5. ,   4.5,   5.5,
         6.5,   6. ,   7. ,   8. ,   9. ,   8.5,   9.5,  10.5])
>>> values_at_links = rmg.empty(at="link")
>>> rtn = map_mean_of_link_nodes_to_link(rmg, "z", out=values_at_links)
>>> values_at_links
array([  0.5,   1.5,   2.5,   2. ,   3. ,   4. ,   5. ,   4.5,   5.5,
         6.5,   6. ,   7. ,   8. ,   9. ,   8.5,   9.5,  10.5])
>>> rtn is values_at_links
True

Map the minimum of a link’s nodes to the link.

map_min_of_link_nodes_to_link iterates across the grid and identifies the node values at both the “head” and “tail” of a given link. This function evaluates the value of ‘var_name’ at both the “to” and “from” node. The minimum value of the two node values is then mapped to the link.

Parameters:
  • var_name (array or field name) – Values defined at nodes.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_min_of_link_nodes_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field(
...     "z",
...     [
...         [0, 1, 2, 3],
...         [7, 6, 5, 4],
...         [8, 9, 10, 11],
...     ],
...     at="node",
... )
>>> map_min_of_link_nodes_to_link(rmg, "z")
array([  0.,   1.,   2.,   0.,   1.,   2.,   3.,   6.,   5.,   4.,   7.,
         6.,   5.,   4.,   8.,   9.,  10.])
>>> values_at_links = rmg.empty(at="link")
>>> rtn = map_min_of_link_nodes_to_link(rmg, "z", out=values_at_links)
>>> values_at_links
array([  0.,   1.,   2.,   0.,   1.,   2.,   3.,   6.,   5.,   4.,   7.,
         6.,   5.,   4.,   8.,   9.,  10.])
>>> rtn is values_at_links
True

Map the minimum value of a nodes’ links to the node.

map_min_of_node_links_to_node iterates across the grid and identifies the link values at each link connected to a node. This function finds the minimum value of ‘var_name’ of each set of links, and then maps this value to the node. Note no attempt is made to honor the directionality of the links.

Parameters:
  • var_name (array or field name) – Values defined at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_min_of_node_links_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = np.arange(rmg.number_of_links)
>>> map_min_of_node_links_to_node(rmg, "grad")
array([  0.,   0.,   1.,   2.,
         3.,   4.,   5.,   6.,
        10.,  11.,  12.,  13.])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_min_of_node_links_to_node(rmg, "grad", out=values_at_nodes)
>>> values_at_nodes
array([  0.,   0.,   1.,   2.,
         3.,   4.,   5.,   6.,
        10.,  11.,  12.,  13.])
>>> rtn is values_at_nodes
True

Assign values to links using a weighted combination of node values.

Assign to each link a weighted combination of values v at nodes using the Lax-Wendroff method for upwind weighting.

c is a scalar or link vector that gives the link-parallel signed Courant number. Where c is positive, velocity is in the direction of the link; where negative, velocity is in the opposite direction.

As an example, consider 3x5 raster grid with the following values at the nodes in the central row:

0---1---2---3---4

Consider a uniform Courant value c = +0.2 at the horizontal links. The mapped link values should be:

.-0.4-.-1.4-.-2.4-.-3.4-.

Values at links when c = -0.2:

.-0.6-.-1.6-.-2.6-.-3.6-.
Parameters:
  • v ((n_nodes,) ndarray) – Values at grid nodes.

  • c (float or (n_links,) ndarray) – Courant number to use at links.

  • out ((n_links,) ndarray, optional) – If provided, place calculated values in this array. Otherwise, create a new array.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 5))
>>> grid.at_node["node_value"] = [
...     [0.0, 0.0, 0.0, 0.0, 0.0],
...     [0.0, 1.0, 2.0, 3.0, 4.0],
...     [0.0, 0.0, 0.0, 0.0, 0.0],
... ]
>>> v = grid.at_node["node_value"]
>>> c = grid.add_zeros("courant_number", at="link")
>>> c[grid.horizontal_links] = 0.2

Set values for the middle row of horizontal links.

>>> val_at_link = grid.map_node_to_link_lax_wendroff(v, c)
>>> val_at_link[9:13]
array([ 0.4,  1.4,  2.4,  3.4])
>>> val_at_link = grid.map_node_to_link_lax_wendroff(v, -c)
>>> val_at_link[9:13]
array([ 0.6,  1.6,  2.6,  3.6])

Assign values to links from upstream nodes.

Assign to each link the value v associated with whichever of its two nodes lies upstream, according to link value u.

Consider a link k with tail node t(k) and head node t(h). Nodes have value v(n). We want to assign a value to links, v’(k). The assignment is:

v'(k) = v(t(k)) where u(k) > 0,
v'(k) = v(h(k)) where u(k) <= 0

As an example, consider 3x5 raster grid with the following values at the nodes in the central row:

0---1---2---3---4

Consider a uniform velocity value u = 1 at the horizontal links. The mapped link values should be:

.-0-.-1-.-2-.-3-.

If u < 0, the link values should be:

.-1-.-2-.-3-.-4-.
Parameters:
  • v ((n_nodes,), ndarray) – Values at grid nodes.

  • u ((n_links,) ndarray) – Values at grid links.

  • out ((n_links,) ndarray, optional) – If provided, place calculated values in this array. Otherwise, create a new array.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 5))
>>> grid.at_node["node_value"] = [
...     [0.0, 0.0, 0.0, 0.0, 0.0],
...     [0.0, 1.0, 2.0, 3.0, 4.0],
...     [0.0, 0.0, 0.0, 0.0, 0.0],
... ]
>>> v = grid.at_node["node_value"]
>>> u = grid.add_zeros("advection_speed", at="link")
>>> u[grid.horizontal_links] = 1.0

Set values for the middle row of horizontal links.

>>> val_at_link = grid.map_node_to_link_linear_upwind(v, u)
>>> val_at_link[9:13]
array([ 0.,  1.,  2.,  3.])
>>> val_at_link = grid.map_node_to_link_linear_upwind(v, -u)
>>> val_at_link[9:13]
array([ 1.,  2.,  3.,  4.])

Map the largest magnitude of the links bringing flux into the node to the node.

map_upwind_node_link_max_to_node iterates across the grid and identifies the link values at each link connected to a node. It then uses the link_dirs_at_node data structure to identify links bringing flux into the node, then maps the maximum magnitude of ‘var_name’ found on these links onto the node. If no upwind link is found, the value will be recorded as zero.

Parameters:
  • var_name (array or field name) – Values defined at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_upwind_node_link_max_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = [
...     -1.1,
...     -1.2,
...     -1.3,
...     1.4,
...     1.5,
...     1.6,
...     -1.7,
...     -1.8,
...     -1.9,
...     2.0,
...     2.1,
...     2.2,
...     -2.3,
...     2.4,
...     2.5,
...     2.6,
...     -2.7,
... ]
>>> map_upwind_node_link_max_to_node(rmg, "grad").reshape((3, 4))
array([[ 1.4,  1.5,  1.6,  1.3],
       [ 2.1,  2.2,  2. ,  2.4],
       [ 2.5,  2.6,  2.3,  2.7]])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_upwind_node_link_max_to_node(rmg, "grad", out=values_at_nodes)
>>> values_at_nodes.reshape((3, 4))
array([[ 1.4,  1.5,  1.6,  1.3],
       [ 2.1,  2.2,  2. ,  2.4],
       [ 2.5,  2.6,  2.3,  2.7]])
>>> rtn is values_at_nodes
True

Map the mean magnitude of the links bringing flux into the node to the node.

map_upwind_node_link_mean_to_node iterates across the grid and identifies the link values at each link connected to a node. It then uses the link_dirs_at_node data structure to identify links bringing flux into the node, then maps the mean magnitude of ‘var_name’ found on these links onto the node. Links with zero values are not included in the means, and zeros are returned if no upwind links are found.

Parameters:
  • var_name (array or field name) – Values defined at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_upwind_node_link_mean_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = [
...     -1.0,
...     -2.0,
...     -1.0,
...     -2.0,
...     -3.0,
...     -4.0,
...     -5.0,
...     -1.0,
...     -2.0,
...     -1.0,
...     -1.0,
...     -2.0,
...     -3.0,
...     -4.0,
...     -1.0,
...     -2.0,
...     -1.0,
... ]
>>> map_upwind_node_link_mean_to_node(rmg, "grad")
array([ 0. ,  1. ,  2. ,  1. ,
        2. ,  2. ,  3. ,  3. ,
        1. ,  1.5,  2.5,  2.5])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_upwind_node_link_mean_to_node(rmg, "grad", out=values_at_nodes)
>>> values_at_nodes
array([ 0. ,  1. ,  2. ,  1. ,
        2. ,  2. ,  3. ,  3. ,
        1. ,  1.5,  2.5,  2.5])
>>> rtn is values_at_nodes
True

Map the the value found in one link array to a node, based on the largest magnitude value of links carrying fluxes out of the node, found in a second node array or field.

map_downwind_node_link_max_to_node iterates across the grid and identifies the link control_values at each link connected to a node. It then uses the link_dirs_at_node data structure to identify links carrying flux out of the node, then identifies the link with the maximum magnitude. The value of the second field ‘value_name’ at these links is then mapped onto the node. If no downwind link is found, the value will be recorded as zero.

Parameters:
  • control_name (array or field name) – Values defined at nodes that dictate which end of the link to draw values from.

  • value_name (array or field name) – Values defined at nodes from which values are drawn, based on control_name.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_value_at_downwind_node_link_max_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = [
...     -1.0,
...     -2.0,
...     -1.0,
...     0.0,
...     0.0,
...     0.0,
...     0.0,
...     -1.0,
...     -2.0,
...     -1.0,
...     0.0,
...     0.0,
...     0.0,
...     0.0,
...     -1.0,
...     -2.0,
...     -1.0,
... ]
>>> rmg.at_link["vals"] = np.arange(rmg.number_of_links, dtype=float)
>>> map_value_at_downwind_node_link_max_to_node(rmg, "grad", "vals")
array([  0.,   1.,   2.,   0.,
         7.,   8.,   9.,   0.,
        14.,  15.,  16.,   0.])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_value_at_downwind_node_link_max_to_node(
...     rmg, "grad", "vals", out=values_at_nodes
... )
>>> values_at_nodes
array([  0.,   1.,   2.,   0.,
         7.,   8.,   9.,   0.,
        14.,  15.,  16.,   0.])
>>> rtn is values_at_nodes
True

Map the the value found in one node array to a link, based on the maximum value found in a second node field or array.

map_value_at_max_node_to_link iterates across the grid and identifies the node values at both the “head” and “tail” of a given link. This function evaluates the value of ‘control_name’ at both the “to” and “from” node. The value of ‘value_name’ at the node with the maximum value of the two values of ‘control_name’ is then mapped to the link.

Parameters:
  • control_name (array or field name) – Name of field defined at nodes or a node array that dictates which end of the link to draw values from.

  • value_name (array or field name) – Name of field defined at nodes or node array from which values are drawn, based on control_name.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_value_at_max_node_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field(
...     "z",
...     [
...         [0, 1, 2, 3],
...         [7, 6, 5, 4],
...         [8, 9, 10, 11],
...     ],
...     at="node",
... )
>>> _ = rmg.add_field(
...     "vals_to_map",
...     [
...         [0, 10, 20, 30],
...         [70, 60, 50, 40],
...         [80, 90, 100, 110],
...     ],
...     at="node",
... )
>>> map_value_at_max_node_to_link(rmg, "z", "vals_to_map")
array([  10.,   20.,   30.,   70.,   60.,   50.,   40.,   70.,   60.,
         50.,   80.,   90.,  100.,  110.,   90.,  100.,  110.])

Map the the value found in one node array to a link, based on the minimum value found in a second node field or array.

map_value_at_min_node_to_link iterates across the grid and identifies the node values at both the “head” and “tail” of a given link. This function evaluates the value of ‘control_name’ at both the “to” and “from” node. The value of ‘value_name’ at the node with the minimum value of the two values of ‘control_name’ is then mapped to the link.

Parameters:
  • control_name (array or field name) – Name of field defined at nodes or a node array that dictates which end of the link to draw values from.

  • value_name (array or field name) – Name of field defined at nodes or node array from which values are drawn, based on control_name.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_value_at_min_node_to_link
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field(
...     "z",
...     [
...         [0, 1, 2, 3],
...         [7, 6, 5, 4],
...         [8, 9, 10, 11],
...     ],
...     at="node",
... )
>>> _ = rmg.add_field(
...     "vals_to_map",
...     [
...         [0, 10, 20, 30],
...         [70, 60, 50, 40],
...         [80, 90, 100, 110],
...     ],
...     at="node",
... )
>>> map_value_at_min_node_to_link(rmg, "z", "vals_to_map")
array([   0.,   10.,   20.,    0.,   10.,   20.,   30.,   60.,   50.,
         40.,   70.,   60.,   50.,   40.,   80.,   90.,  100.])

Map the the value found in one link array to a node, based on the largest magnitude value of links bringing fluxes into the node, found in a second node array or field.

map_upwind_node_link_max_to_node iterates across the grid and identifies the link control_values at each link connected to a node. It then uses the link_dirs_at_node data structure to identify links bringing flux into the node, then identifies the link with the maximum magnitude. The value of the second field ‘value_name’ at these links is then mapped onto the node. If no upwind link is found, the value will be recorded as zero.

Parameters:
  • control_name (array or field name) – Values defined at nodes that dictate which end of the link to draw values from.

  • value_name (array or field name) – Values defined at nodes from which values are drawn, based on control_name.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_value_at_upwind_node_link_max_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_link["grad"] = [
...     -1.0,
...     -2.0,
...     -1.0,
...     0.0,
...     0.0,
...     0.0,
...     0.0,
...     -1.0,
...     -2.0,
...     -1.0,
...     0.0,
...     0.0,
...     0.0,
...     0.0,
...     -1.0,
...     -2.0,
...     -1.0,
... ]
>>> rmg.at_link["vals"] = np.arange(rmg.number_of_links, dtype=float)
>>> map_value_at_upwind_node_link_max_to_node(rmg, "grad", "vals")
array([  0.,   0.,   1.,   2.,
         0.,   7.,   8.,   9.,
         0.,  14.,  15.,  16.])
>>> values_at_nodes = rmg.add_empty("z", at="node")
>>> rtn = map_value_at_upwind_node_link_max_to_node(
...     rmg, "grad", "vals", out=values_at_nodes
... )
>>> values_at_nodes
array([  0.,   0.,   1.,   2.,
         0.,   7.,   8.,   9.,
         0.,  14.,  15.,  16.])
>>> rtn is values_at_nodes
True

Map magnitude and sign of vectors with components (ux, uy) onto grid links.

Examples

>>> from landlab import HexModelGrid
>>> import numpy
>>> hmg = HexModelGrid((3, 2))
>>> (numpy.round(10 * map_vectors_to_links(hmg, 1.0, 0.0))).astype(int)
array([10, -5,  5, -5,  5, 10, 10,  5, -5,  5, -5, 10])

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#
new_field_location(loc, size=None)#

Add a new quantity to a field.

Create an empty group into which new fields can be added. The new group is created but no memory allocated yet. The dictionary of the new group can be through a new at_ attribute of the class instance.

Parameters:
  • loc (str) – Name of the new group to add to the field.

  • size (int, optional) – Number of elements in the new quantity. If not provided, the size is set to be the size of the first field added to the group.

Raises:

ValueError – If the field already contains the group.

Examples

Create a collection of fields and add two groups, node and cell, to it.

>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 12)
>>> fields.new_field_location("cell", 2)

The group names in the collection are retrieved with the groups attribute as a set.

>>> names = list(fields.groups)
>>> names.sort()
>>> names
['cell', 'node']

Access the new (empty) groups with the at_ attributes.

>>> fields.at_cell
FieldDataset('cell', size=2, fixed_size=True)
>>> fields.at_node
FieldDataset('node', size=12, fixed_size=True)
>>> fields.new_field_location("core_node")
>>> fields.at_core_node.size is None
True
>>> fields.at_core_node["air__temperature"] = [0, 1]
>>> fields.at_core_node.size
2

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_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])

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]])

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
ones(*args, **kwds)#

Array, initialized to 1, whose size is that of the field.

Return a new array of the data field size, filled with ones. Keyword arguments are the same as that for the equivalent numpy function.

Parameters:

group (str) – Name of the group.

See also

numpy.ones

See for a description of optional keywords.

empty

Equivalent method that does not initialize the new array.

zeros

Equivalent method that initializes the data to 0.

Examples

>>> from landlab.field import GraphFields
>>> field = GraphFields()
>>> field.new_field_location("node", 4)
>>> field.ones("node")
array([ 1.,  1.,  1.,  1.])
>>> field.ones("node", dtype=int)
array([1, 1, 1, 1])

Note that a new field is not added to the collection of fields.

>>> list(field.keys("node"))
[]
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])
reset_status_at_node()[source]#
return_array_or_field_values(field, at=None)#

Return field given a field name, or array of with the correct shape.

Given a group and a field, return a reference to the associated data array. field is either a string that is a field in the group or an array of the correct size.

This function is meant to serve like the use_field_name_or_array decorator for bound functions.

Parameters:
  • field (str or array) – Name of the field withing group.

  • at (str, optional) – Name of the group.

Returns:

The values of the field.

Return type:

numpy.ndarray

Raises:
  • landlab.field.errors.GroupError – If group does not exist

  • landlab.field.errors.FieldError – If field does not exist

Examples

Create a group of fields called node.

>>> import numpy as np
>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 4)

Add a field, initialized to ones, called topographic__elevation to the node group. The field_values method returns a reference to the field’s data.

>>> _ = fields.add_ones("topographic__elevation", at="node")
>>> fields.field_values("topographic__elevation", at="node")
array([ 1.,  1.,  1.,  1.])

Alternatively, if the second argument is an array, its size is checked and returned if correct.

>>> vals = np.array([4.0, 5.0, 7.0, 3.0])
>>> fields.return_array_or_field_values(vals, at="node")
array([ 4.,  5.,  7.,  3.])

Raise FieldError if field does not exist in group.

>>> fields.return_array_or_field_values("surface__temperature", at="node")
Traceback (most recent call last):
...
FieldError: surface__temperature

If group does not exists, Raise GroupError.

>>> fields.return_array_or_field_values("topographic__elevation", at="cell")
Traceback (most recent call last):
...
GroupError: cell

And if the array of values provided is incorrect, raise a ValueError.

>>> vals = np.array([3.0, 2.0, 1.0])
>>> fields.return_array_or_field_values(vals, at="node")
Traceback (most recent call last):
...
ValueError: Array has incorrect size.
size(group)#

Return the size of the arrays stored in a group.

Parameters:

group (str) – Group name.

Returns:

Array size.

Return type:

int

Examples

>>> from landlab.field import GraphFields
>>> fields = GraphFields()
>>> fields.new_field_location("node", 4)
>>> fields.size("node")
4
sort()#

Sort graph elements.

Get array of the status of all links.

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.status_at_link
array([0, 0, 0], dtype=uint8)
property status_at_node#

Get array of the boundary status for each node.

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.status_at_node
array([0, 0, 0, 0], dtype=uint8)
>>> grid.status_at_link
array([0, 0, 0], dtype=uint8)

Now we change the status at node 0 to a closed boundary. This will result in changing the link status as well.

>>> grid.status_at_node = [
...     grid.BC_NODE_IS_CLOSED,
...     grid.BC_NODE_IS_CORE,
...     grid.BC_NODE_IS_CORE,
...     grid.BC_NODE_IS_CORE,
... ]
>>> grid.status_at_node
array([4, 0, 0, 0], dtype=uint8)
>>> grid.status_at_link
array([4, 0, 0], dtype=uint8)
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])

Get array of the x-coordinates of link midpoints.

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.x_of_link
array([ 0. , -0.5,  0.5])
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_reference#

Return the coordinates (x, y) of the reference point.

For RasterModelGrid and HexModelGrid the reference point is the minimum of x_of_node and of y_of_node. By default it is (0, 0). For VoronoiDelaunayGrid the reference point is (0, 0). For RadialModelGrid it is the (x, y) of the center point.

The intention of these coordinates is to provide a method to store the large float values of projected coordinates.

Example

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid(
...     (y_of_node, x_of_node), nodes_at_link, xy_of_reference=(12345, 678910)
... )
>>> grid.xy_of_reference
(12345, 678910)
>>> grid.xy_of_reference = (98765, 43210)
>>> grid.xy_of_reference
(98765, 43210)

Get array of the y-coordinates of link midpoints.

Examples

>>> from landlab import NetworkModelGrid
>>> y_of_node = (0, 1, 2, 2)
>>> x_of_node = (0, 0, -1, 1)
>>> nodes_at_link = ((1, 0), (2, 1), (3, 1))
>>> grid = NetworkModelGrid((y_of_node, x_of_node), nodes_at_link)
>>> grid.y_of_link
array([ 0.5,  1.5,  1.5])
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.])
zeros(*args, **kwds)#

Array, initialized to 0, whose size is that of the field.

Parameters:
  • group (str) – Name of the group.

  • size (Return a new array of the data field) –

  • Keyword (filled with zeros.) –

  • function. (arguments are the same as that for the equivalent numpy) –

  • grid. (This method is not valid for the group) –

See also

numpy.zeros

See for a description of optional keywords.

empty

Equivalent method that does not initialize the new array.

ones

Equivalent method that initializes the data to 1.

Examples

>>> from landlab.field import GraphFields
>>> field = GraphFields()
>>> field.new_field_location("node", 4)
>>> field.zeros("node")
array([ 0.,  0.,  0.,  0.])

Note that a new field is not added to the collection of fields.

>>> list(field.keys("node"))
[]