API for landlab.grid.raster#

A class used to create and manage regular raster grids for 2D numerical models in Landlab.

Do NOT add new documentation here. Grid documentation is now built in a semi- automated fashion. To modify the text seen on the web, edit the files docs/text_for_[gridfile].py.txt.

class RasterModelGrid(shape, xy_spacing=1.0, xy_of_lower_left=(0.0, 0.0), xy_of_reference=(0.0, 0.0), xy_axis_name=('x', 'y'), xy_axis_units='-', bc=None)[source]#

Bases: DiagonalsMixIn, DualUniformRectilinearGraph, ModelGrid

A 2D uniform rectilinear grid.

Examples

Create a uniform rectilinear grid that has 4 rows and 5 columns of nodes. Nodes along the edges will be open. That is, links connecting these nodes to core nodes are active.

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5))
>>> rmg.number_of_node_rows, rmg.number_of_node_columns
(4, 5)
>>> rmg.number_of_active_links
17

Set the nodes along the top edge of the grid to be closed boundaries. This means that any links touching these nodes will be inactive.

>>> rmg = RasterModelGrid((4, 5), bc={"top": "closed"})
>>> rmg.number_of_node_rows, rmg.number_of_node_columns
(4, 5)
>>> rmg.number_of_active_links
14

A RasterModelGrid can have different node spacings in the x and y directions.

>>> grid = RasterModelGrid((4, 5), xy_spacing=(2, 1))
>>> grid.dx, grid.dy
(2.0, 1.0)
>>> grid.node_y.reshape(grid.shape)
array([[ 0.,  0.,  0.,  0.,  0.],
       [ 1.,  1.,  1.,  1.,  1.],
       [ 2.,  2.,  2.,  2.,  2.],
       [ 3.,  3.,  3.,  3.,  3.]])
>>> grid.node_x.reshape(grid.shape)
array([[ 0.,  2.,  4.,  6.,  8.],
       [ 0.,  2.,  4.,  6.,  8.],
       [ 0.,  2.,  4.,  6.,  8.],
       [ 0.,  2.,  4.,  6.,  8.]])

Create a 2D grid with equal spacing.

Optionally takes numbers of rows and columns and cell size as inputs. If this are given, calls initialize() to set up the grid. At the moment, num_rows and num_cols MUST be specified. Both must be >=3 to allow correct automated setup of boundary conditions.

Parameters:
  • shape (tuple of int) – Shape of the grid in nodes as (nrows, ncols).

  • xy_spacing (tuple or float, optional) – dx and dy spacing. Either provided as a float or a (dx, dy) tuple.

  • xy_of_lower_left (tuple, optional) – (x, y) coordinates of the lower left corner.

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

  • xy_axis_name (tuple of str) – Name to use for each axis.

  • xy_axis_units (tuple of str, or str) – Units for coordinates of each axis.

  • bc (dict, optional) – Edge boundary conditions.

Returns:

A newly-created grid.

Return type:

RasterModelGrid

Notes

The option for NOT giving rows, cols, and dx no longer works, because the field init requires num_active_cells, etc., to be defined. Either we force users to give arguments on instantiation, or set it up such that one can create a zero-node grid.

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', 'patch', 'corner', 'face', 'cell', 'grid')#

Grid elements on which fields can be placed.

__getitem__(name)#

Get the collection of fields from the named group.

__getstate__()[source]#

Get state for pickling.

__init__(shape, xy_spacing=1.0, xy_of_lower_left=(0.0, 0.0), xy_of_reference=(0.0, 0.0), xy_axis_name=('x', 'y'), xy_axis_units='-', bc=None)[source]#

Create a 2D grid with equal spacing.

Optionally takes numbers of rows and columns and cell size as inputs. If this are given, calls initialize() to set up the grid. At the moment, num_rows and num_cols MUST be specified. Both must be >=3 to allow correct automated setup of boundary conditions.

Parameters:
  • shape (tuple of int) – Shape of the grid in nodes as (nrows, ncols).

  • xy_spacing (tuple or float, optional) – dx and dy spacing. Either provided as a float or a (dx, dy) tuple.

  • xy_of_lower_left (tuple, optional) – (x, y) coordinates of the lower left corner.

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

  • xy_axis_name (tuple of str) – Name to use for each axis.

  • xy_axis_units (tuple of str, or str) – Units for coordinates of each axis.

  • bc (dict, optional) – Edge boundary conditions.

Returns:

A newly-created grid.

Return type:

RasterModelGrid

Notes

The option for NOT giving rows, cols, and dx no longer works, because the field init requires num_active_cells, etc., to be defined. Either we force users to give arguments on instantiation, or set it up such that one can create a zero-node grid.

__setstate__(state_dict)[source]#

Set state for of RasterModelGrid from pickled state_dict.

property active_adjacent_corners_at_corner#

Adjacent corners for each grid corner.

property active_adjacent_nodes_at_node#

Adjacent nodes for each grid node.

For each grid node, get the adjacent nodes ordered counterclockwise starting from the positive x axis.

Examples

>>> from landlab import RasterModelGrid, HexModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.active_adjacent_nodes_at_node[(-1, 6, 2),]
array([[-1, -1, -1, -1],
       [ 7, 11,  5,  1],
       [-1,  7, -1, -1]])

Setting a node to closed causes all links touching it to be inactive.

>>> grid.status_at_node[6] = grid.BC_NODE_IS_CLOSED
>>> grid.active_adjacent_nodes_at_node[(-1, 6, 2),]
array([[-1, -1, -1, -1],
       [-1, -1, -1, -1],
       [-1,  7, -1, -1]])
>>> grid.active_adjacent_nodes_at_node[7]
array([ 8, 12, -1,  2])
>>> grid.active_adjacent_nodes_at_node[2]
array([-1,  7, -1, -1])
>>> grid = HexModelGrid((3, 2))
>>> grid.status_at_node[0] = grid.BC_NODE_IS_CLOSED
>>> grid.active_adjacent_nodes_at_node
array([[-1, -1, -1, -1, -1, -1],
       [-1,  3, -1, -1, -1, -1],
       [ 3, -1, -1, -1, -1, -1],
       [ 4,  6,  5,  2, -1,  1],
       [-1,  3, -1, -1, -1, -1],
       [-1, -1,  3, -1, -1, -1],
       [-1,  3, -1, -1, -1, -1]])
property active_d8#
property active_d8_dirs_at_corner#
property active_d8_dirs_at_node#
property active_diagonal_dirs_at_corner#
property active_diagonal_dirs_at_node#
property active_diagonals#
property active_face_dirs_at_corner#

Return face directions into each corner.

property active_faces#

Get array of active faces.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.active_faces
array([0, 1, 2, 3, 4, 5, 6])
>>> grid.status_at_node[6] = grid.BC_NODE_IS_CLOSED
>>> grid.active_faces
array([0, 2, 5])

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. Note that inactive links have a value of 0, but active and fixed links are both reported normally.

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 import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
.--->.--->.--->.
^    ^    ^    ^
|    |    |    |
.--->5--->.--->.
^    ^    ^    ^
|    |    |    |
.--->.--->.--->.
>>> grid.active_link_dirs_at_node[5]
array([-1, -1,  1,  1], dtype=int8)

If we set the nodes along the left edge to be closed, the links attached to those nodes become inactive and so their directions are reported as 0.

>>> grid.status_at_node[grid.nodes_at_left_edge] = grid.BC_NODE_IS_CLOSED
>>> grid.active_link_dirs_at_node
array([[ 0,  0,  0,  0], [ 0, -1,  0,  0], [ 0, -1,  0,  0], [ 0,  0,  0,  0],
       [ 0,  0,  0,  0], [-1, -1,  0,  1], [-1, -1,  1,  1], [ 0,  0,  1,  0],
       [ 0,  0,  0,  0], [ 0,  0,  0,  1], [ 0,  0,  0,  1], [ 0,  0,  0,  0]],
       dtype=int8)

Get array of active links.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.active_links
array([ 4,  5,  7,  8,  9, 11, 12])
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_corners_at_corner#

Get adjacent corners.

property adjacent_faces_at_face#
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]])
property all_corner_azimuths_map#

Get azimuths from every corner to every other corner.

property all_corner_distances_map#

Get distances from every corner to every other corner.

property all_node_azimuths_map#

Get azimuths from every node to every other node.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> angles = grid.all_node_azimuths_map

The shape of the array is number_of_nodes by number_of_nodes and azimuth from a node to itself is zero.

>>> angles.shape == (grid.number_of_nodes, grid.number_of_nodes)
True
>>> angles.diagonal()
array([ 0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.])

Angles are measured in radians and increase clockwise starting at north.

>>> angles *= 180.0 / np.pi
>>> angles[0, :4]
array([  0.,  90.,  90.,  90.])
>>> angles[0, ::4]
array([ 0.,  0.,  0.])
>>> angles[0, ::5]
array([  0.,  45.,  45.])
property all_node_distances_map#

Get distances from every node to every other node.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> distances = grid.all_node_distances_map

The shape of the array is number_of_nodes by number_of_nodes and distance from a node to itself is zero.

>>> distances.shape == (grid.number_of_nodes, grid.number_of_nodes)
True
>>> distances.diagonal()
array([ 0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.])

The distances from the first node to all nodes in its row and all the nodes in its column.

>>> distances[0, :4]
array([ 0.,  1.,  2.,  3.])
>>> distances[0, ::4]
array([ 0.,  1.,  2.])
property angle_of_face#

Get the angle of each face.

See also

angle_of_link

property angle_of_face_about_head#

Find and return the angle of a face about the corner at the face head.

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

Find and return the angle of a link about the node at the link head.

Because links have direction, their angle can be specified as an angle about either the node at the link head, or the node at the link tail. The default behaviour of angle_of_link is to return the angle about the link tail, but this method gives the angle about the link head.

Examples

>>> from landlab import HexModelGrid
>>> import numpy as np
>>> grid = HexModelGrid((3, 2), node_layout="hex")
>>> np.round(grid.angle_of_link[:3] / np.pi * 3.0)
array([ 0., 2.,  1.])
>>> np.round(grid.angle_of_link_about_head[:3] / np.pi * 3.0)  # 60 deg segments
array([ 3.,  5.,  4.])
property area_of_cell#

Get the area of each cell.

See also

area_of_patch

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.])
as_dataarray(name, at=None, time=None)#

Create an xarray DataArray representation of a grid field.

Parameters:
  • name (str) – Name of a field. This can either be a canonical field name (of the form “at_<element>:<field_name>”, or just the field name. In the latter case, use the at keyword to specify where the field is defined.

  • at (str, optional) – The grid elements on which the field is defined. Use this only if name is not a canonical field name that already contains the grid element.

Returns:

The field represented as a newly-created xarray DataArray.

Return type:

DataArray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> _ = grid.add_full("elevation", 3.0, at="node")
>>> grid.as_dataarray("at_node:elevation")
<xarray.DataArray 'at_node:elevation' (node: 12)>
array([ 3.,  3.,  3.,  3.,  3.,  3.,  3.,  3.,  3.,  3.,  3.,  3.])
Dimensions without coordinates: node
>>> all(
...     grid.as_dataarray("at_node:elevation")
...     == grid.as_dataarray("elevation", at="node")
... )
True
as_dataset(include='*', exclude=None, time=None)[source]#

Create an xarray Dataset representation of a grid.

This method creates a new xarray Dataset object that contains the grid’s data fields. A particular grid type (e.g. a RasterModelGrid) should define its own as_dataset method to represent that particular grid type as a Dataset and then call this method to add the data fields.

Parameters:
  • include (str or iterable or str) – Glob-style patterns of fields to include in the dataset.

  • exclude (str or iterable or str) – Glob-style patterns of fields to exclude from the dataset.

Returns:

An xarray Dataset representation of a ModelGrid.

Return type:

Dataset

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))

Add some fields to the grid. Notice that we have defined a field named “elevation” at both nodes and links.

>>> _ = grid.add_full("elevation", 3.0, at="node")
>>> _ = grid.add_full("elevation", 4.0, at="link")
>>> _ = grid.add_full("temperature", 5.0, at="node")
>>> ds = grid.as_dataset()
>>> sorted(ds.dims.items())
[('dim', 2), ('link', 17), ('node', 12)]
>>> sorted([var for var in ds.data_vars if var.startswith("at_")])
['at_link:elevation', 'at_node:elevation', 'at_node:temperature']
>>> grid.event_layers.add(1.0, rho=0.5)
>>> ds = grid.as_dataset()
>>> sorted(ds.dims.items())
[('cell', 2), ('dim', 2), ('layer', 1), ('link', 17), ('node', 12)]
>>> sorted([var for var in ds.data_vars if var.startswith("at_")])
['at_layer_cell:rho', 'at_layer_cell:thickness', 'at_link:elevation',
 'at_node:elevation', 'at_node:temperature']
at_cell = {}#
at_corner = {}#
at_face = {}#
at_grid = {}#
property at_layer#

EventLayers for each cell.

at_node = {}#
at_patch = {}#
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 RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.axis_name
('x', 'y')
>>> grid.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 RasterModelGrid
>>> mg = RasterModelGrid((4, 5), xy_spacing=(3.0, 2.0))
>>> mg.axis_units
('-', '-')
>>> mg.axis_units = ("degrees_north", "degrees_east")
>>> mg.axis_units
('degrees_north', 'degrees_east')
>>> mg.axis_units = "m"
>>> mg.axis_units
('m', 'm')
property boundary_corners#

Get array of boundary corners.

See also

boundary_nodes

property boundary_nodes#

Get array of boundary nodes.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.boundary_nodes
array([ 0,  1,  2,  3,  4,  5,  9, 10, 14, 15, 16, 17, 18, 19])
calc_aspect_at_cell_subtriangles(elevs='topographic__elevation', subtriangle_unit_normals=None, unit='degrees')#

Get tuple of arrays of aspect of each of the eight cell subtriangles.

Aspect is returned as radians clockwise of north, unless input parameter units is set to ‘degrees’.

If subtriangle_unit_normals is provided the aspect will be calculated from these data.

If it is not, it will be derived from elevation data at the nodes, which can either be a string referring to a grid field (default: ‘topographic__elevation’), or an nnodes-long numpy array of the values themselves.

Parameters:
  • elevs (str or array (optional)) – Node field name or node array of elevations. If subtriangle_unit_normals is not provided, must be set, but unused otherwise.

  • subtriangle_unit_normals (tuple of 8 (ncels, 3) arrays (optional)) – The unit normal vectors for the eight subtriangles of each cell, if already known. Order is from north of east, counter clockwise to south of east (East North East, North North East, North North West, West North West, West South West, South South West, South South East, East South East).

  • unit ({'degrees', 'radians'}) – Controls the unit that the aspect is returned as.

Returns:

each a length num-cells array Len-8 tuple of the aspect of each of the eight cell subtriangles. Aspect is returned as angle clockwise of north. Units are given as radians unless input parameter units is set to ‘degrees’. Order is from north of east, counter clockwise to south of east (East North East, North North East, North North West, West North West, West South West, South South West, South South East, East South East).

Return type:

(a_ENE, a_NNE, a_NNW, a_WNW, a_WSW, a_SSW, a_SSE, a_ESE)

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> z = np.array([1.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0])
>>> eight_tris = mg.calc_unit_normals_at_cell_subtriangles(z)
>>> A = mg.calc_aspect_at_cell_subtriangles(z, eight_tris)
>>> A0 = mg.calc_aspect_at_cell_subtriangles(z)
>>> np.allclose(A, A0)
True
>>> type(A) is tuple
True
>>> len(A)
8
>>> len(A[0]) == mg.number_of_cells
True
>>> A0
(array([ 180.]), array([ 270.]), array([ 90.]), array([ 180.]),
 array([ 0.]), array([ 90.]), array([ 270.]), array([ 0.]))
calc_aspect_at_node(slope_component_tuple=None, elevs='topographic__elevation', unit='degrees', ignore_closed_nodes=True)#

Get array of aspect of a surface.

Calculates at returns the aspect of a surface. Aspect is returned as radians clockwise of north, unless input parameter units is set to ‘degrees’.

If slope_component_tuple is provided, i.e., (slope_x, slope_y), the aspect will be calculated from these data.

If it is not, it will be derived from elevation data at the nodes, which can either be a string referring to a grid field (default: ‘topographic__elevation’), or an nnodes-long numpy array of the values themselves.

If ignore_closed_nodes is False, all proximal elevation values will be used in the calculation. If True, only unclosed nodes are used.

Parameters:
  • slope_component_tuple ((slope_x_array, slope_y_array) (optional)) – Tuple of components of slope in the x and y directions, defined on nodes, if already known. If not, provide elevs.

  • elevs (str or array (optional)) – Node field name or node array of elevations. If slope_component_tuple is not provided, must be set, but unused otherwise.

  • unit ({'degrees', 'radians'}) – Controls the unit that the aspect is returned as.

  • ignore_closed_nodes (bool) – If True, do not incorporate values at closed nodes into the calc.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 4))
>>> z = mg.node_x**2 + mg.node_y**2
>>> mg.calc_aspect_at_node(elevs=z)
array([ 225.        ,  240.16585039,  255.2796318 ,  258.69006753,
        209.83414961,  225.        ,  243.54632481,  248.77808974,
        194.7203682 ,  206.45367519,  225.        ,  231.94498651,
        191.30993247,  201.22191026,  218.05501349,  225.        ])
>>> z = z.max() - z
>>> mg.calc_aspect_at_node(elevs=z)
array([ 45.        ,  60.16585039,  75.2796318 ,  78.69006753,
        29.83414961,  45.        ,  63.54632481,  68.77808974,
        14.7203682 ,  26.45367519,  45.        ,  51.94498651,
        11.30993247,  21.22191026,  38.05501349,  45.        ])
>>> mg = RasterModelGrid((4, 4), xy_spacing=(3.0, 2.0))
>>> z = mg.node_x**2 + mg.node_y**2
>>> mg.calc_aspect_at_node(elevs=z)
array([ 236.30993247,  247.52001262,  259.97326008,  262.40535663,
        220.75264634,  234.41577266,  251.13402374,  255.29210302,
        201.54258265,  215.47930877,  235.73541937,  242.24162456,
        196.69924423,  209.43534223,  229.19345757,  236.30993247])

Note that a small amount of asymmetry arises at the grid edges due to the “missing” nodes beyond the edge of the grid.

calc_diff_at_d8(node_values, out=None)#

Calculate differences of node values over links and diagonals.

Calculates the difference in quantity node_values at each link in the grid.

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

  • out (ndarray, optional) – Buffer to hold the result.

Returns:

Differences across links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4), xy_spacing=(4, 3))
>>> z = [
...     [60.0, 60.0, 60.0, 60.0],
...     [60.0, 60.0, 0.0, 0.0],
...     [60.0, 0.0, 0.0, 0.0],
... ]
>>> grid.calc_diff_at_d8(z)
array([  0.,   0.,   0.,   0.,   0., -60., -60.,   0., -60.,   0.,   0.,
       -60.,   0.,   0., -60.,   0.,   0.,   0.,   0., -60.,   0., -60.,
       -60., -60.,   0., -60.,   0.,   0.,   0.])
calc_diff_at_diagonal(node_values, out=None)#

Calculate differences of node values over diagonals.

Calculates the difference in quantity node_values at each link in the grid.

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

  • out (ndarray, optional) – Buffer to hold the result.

Returns:

Differences across links.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4), xy_spacing=(4, 3))
>>> z = [
...     [5.0, 5.0, 5.0, 5.0],
...     [5.0, 5.0, 0.0, 0.0],
...     [5.0, 0.0, 0.0, 0.0],
... ]
>>> grid.calc_diff_at_diagonal(z)
array([ 0.,  0., -5.,  0., -5., -5., -5.,  0., -5.,  0.,  0.,  0.])

Calculate differences in node_values at links.

Parameters:
  • value_at_node (array_like or field name) – Values at nodes.

  • out (ndarray, optional) – Buffer to hold result. If None, create a new array.

Returns:

Differences of the nodes values for each link.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3))
>>> node_values = [
...     [0.0, 0.0, 0.0],
...     [1.0, 3.0, 1.0],
...     [2.0, 2.0, 2.0],
... ]
>>> grid.calc_diff_at_link(node_values)
array([ 0.,  0.,  1.,  3.,  1.,  2., -2.,  1., -1.,  1.,  0.,  0.])
>>> out = np.empty(grid.number_of_links, dtype=float)
>>> rtn = grid.calc_diff_at_link(node_values, out=out)
>>> rtn is out
True
>>> out
array([ 0.,  0.,  1.,  3.,  1.,  2., -2.,  1., -1.,  1.,  0.,  0.])
>>> grid = RasterModelGrid((3, 3), xy_spacing=(2, 1))
>>> grid.calc_diff_at_link(node_values)
array([ 0.,  0.,  1.,  3.,  1.,  2., -2.,  1., -1.,  1.,  0.,  0.])
>>> _ = grid.add_field("elevation", node_values, at="node")
>>> grid.calc_diff_at_link("elevation")
array([ 0.,  0.,  1.,  3.,  1.,  2., -2.,  1., -1.,  1.,  0.,  0.])
calc_distances_of_nodes_to_point(coord, get_az=None, node_subset=None, out_distance=None, out_azimuth=None)#

Get distances for nodes to a given point.

Returns an array of distances for each node to a provided point. If “get_az” is set to ‘angles’, returns both the distance array and an array of azimuths from up/north. If it is set to ‘displacements’, it returns the azimuths as a 2xnnodes array of x and y displacements. If it is not set, returns just the distance array.

If “node_subset” is set as an ID, or list/array/etc of IDs method returns just the distance (and optionally azimuth) for that node. Point is provided as a tuple (x,y).

If out_distance (& out_azimuth) are provided, these arrays are used to store the outputs. This is recommended for memory management reasons if you are working with node subsets.

Note

Angles are returned in radians but measured clockwise from north.

Parameters:
  • coord (tuple of float) – Coodinates of point as (x, y).

  • get_az ({None, 'angles', 'displacements'}, optional) – Optionally calculate azimuths as either angles or displacements. The calculated values will be returned along with the distances as the second item of a tuple.

  • node_subset (array_like, optional) – Calculate distances on a subset of grid nodes. The default is to calculate distances from the provided points to all nodes.

  • out_distance (array_like, optional) – If provided, put the calculated distances here. Otherwise, create a new array.

  • out_azimuth (array_like, optional) – If provided, put the calculated distances here. Otherwise, create a new array.

Returns:

If get_az is None return the array of distances. Otherwise, return a tuple of distances and azimuths.

Return type:

ndarray or tuple of ndarray

Notes

Once you start working with node subsets in Landlab, which can change size between loops, it’s quite possible for Python’s internal memory management to crap out after large numbers of loops (~>10k). This is to do with the way it block allocates memory for arrays of differing lengths, then cannot free this memory effectively. The solution - as implemented here - is to pre-allocate all arrays as nnodes long, then only work with the first [len_subset] entries by slicing (in a pseudo-C-style). Care has to be taken not to “accidentally” allow Python to allocate a new array you don’t have control over. Then, to maintain efficient memory allocation, we create some “dummy” nnode-long arrays to store intermediate parts of the solution in.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))

Calculate distances from point at (2., 1.) to a subset of nodes on the grid.

>>> grid.calc_distances_of_nodes_to_point((2, 1), node_subset=(2, 6, 7, 8, 12))
array([ 1.,  1.,  0.,  1.,  1.])

Calculate distances from a point to all nodes on the grid.

>>> dist = grid.calc_distances_of_nodes_to_point((2, 1))
>>> dist.shape == (grid.number_of_nodes,)
True
>>> dist.take((2, 6, 7, 8, 12))
array([ 1.,  1.,  0.,  1.,  1.])

Put the distances into a buffer.

>>> out = np.empty(grid.number_of_nodes, dtype=float)
>>> dist = grid.calc_distances_of_nodes_to_point((2, 1), out_distance=out)
>>> out is dist
True
>>> out.take((2, 6, 7, 8, 12))
array([ 1.,  1.,  0.,  1.,  1.])

Calculate azimuths along with distances. The azimuths are calculated in radians but measured clockwise from north.

>>> (_, azim) = grid.calc_distances_of_nodes_to_point((2, 1), get_az="angles")
>>> azim.take((2, 6, 7, 8, 12)) * 180.0 / np.pi
array([ 180.,  270.,    0.,   90.,    0.])
>>> (_, azim) = grid.calc_distances_of_nodes_to_point(
...     (2, 1), get_az="angles", node_subset=(1, 3, 11, 13)
... )
>>> azim * 180.0 / np.pi
array([ 225.,  135.,  315.,   45.])

When calculating displacements, the first row contains displacements in x and the second displacements in y.

>>> (_, azim) = grid.calc_distances_of_nodes_to_point(
...     (2, 1), get_az="displacements", node_subset=(2, 6, 7, 8, 12)
... )
>>> azim
array([[ 0., -1.,  0.,  1.,  0.],
       [-1.,  0.,  0.,  0.,  1.]])
calc_flux_div_at_cell(unit_flux, out=None)#

Calculate divergence of link-based fluxes at cells.

This function is very similar to the function calc_flux_div_at_node.

Given a flux per unit width across each cell face in the grid, calculate the net outflux (or influx, if negative) divided by cell area, at each cell.

Parameters:

unit_flux_at_links_across_faces (ndarray or field name) – Flux per unit width along links at faces (x number of faces) or link field.

Returns:

Flux divergence at cells.

Return type:

ndarray (x number of cells)

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.grid.divergence import calc_flux_div_at_cell
>>> rg = RasterModelGrid((3, 4), xy_spacing=10.0)
>>> import numpy as np
>>> z = rg.add_zeros("topographic__elevation", at="node")
>>> z[5] = 50.0
>>> z[6] = 36.0
>>> lg = rg.calc_grad_at_link(z)  # there are 17 links
>>> lg
array([ 0. ,  0. ,  0. ,  0. ,  5. ,  3.6,  0. ,  5. , -1.4, -3.6,  0. ,
       -5. , -3.6,  0. ,  0. ,  0. ,  0. ])
>>> fg = lg[rg.link_at_face]  # there are 7 faces
>>> fg
array([ 5. ,  3.6,  5. , -1.4, -3.6, -5. , -3.6])
>>> calc_flux_div_at_cell(rg, -fg)
array([ 1.64,  0.94])
>>> rg.set_status_at_node_on_edges(right=rg.BC_NODE_IS_CLOSED)
>>> rg.set_status_at_node_on_edges(top=rg.BC_NODE_IS_CLOSED)
>>> unit_flux_at_faces = np.zeros(rg.number_of_faces)
>>> unit_flux_at_faces[rg.active_faces] = -fg[rg.active_faces]
>>> calc_flux_div_at_cell(rg, unit_flux_at_faces)
array([ 1.14,  0.22])
>>> _ = rg.add_field("neg_grad_at_link", -lg, at="link")
>>> calc_flux_div_at_cell(rg, "neg_grad_at_link")
array([ 1.64,  0.94])

Notes

Performs a numerical flux divergence operation at cells.

calc_flux_div_at_node(unit_flux, out=None)#

Calculate divergence of link-based fluxes at nodes.

Given a flux per unit width across each face in the grid, calculate the net outflux (or influx, if negative) divided by cell area, at each node (zero or “out” value for nodes without cells).

Parameters:

unit_flux (ndarray or field name) – Flux per unit width along links (x number of links).

Returns:

Flux divergence at nodes.

Return type:

ndarray (x number of nodes)

Examples

>>> from landlab import RasterModelGrid
>>> from landlab.grid.divergence import calc_flux_div_at_node
>>> rg = RasterModelGrid((3, 4), xy_spacing=10.0)
>>> z = rg.add_zeros("topographic__elevation", at="node")
>>> z[5] = 50.0
>>> z[6] = 36.0
>>> lg = rg.calc_grad_at_link(z)  # there are 17 links
>>> lg
array([ 0. ,  0. ,  0. ,  0. ,  5. ,  3.6,  0. ,  5. , -1.4, -3.6,  0. ,
       -5. , -3.6,  0. ,  0. ,  0. ,  0. ])
>>> calc_flux_div_at_node(rg, -lg)
array([ 0.  ,  0.  ,  0.  ,  0.  ,  0.  ,  1.64,  0.94,  0.  ,  0.  ,
        0.  ,  0.  ,  0.  ])
>>> rg.set_status_at_node_on_edges(right=rg.BC_NODE_IS_CLOSED)
>>> rg.set_status_at_node_on_edges(top=rg.BC_NODE_IS_CLOSED)
>>> unit_flux_at_links = np.zeros(rg.number_of_links)
>>> unit_flux_at_links[rg.active_links] = -lg[rg.active_links]
>>> calc_flux_div_at_node(rg, unit_flux_at_links)
array([ 0.  ,  0.  ,  0.  ,  0.  ,  0.  ,  1.14,  0.22,  0.  ,  0.  ,
        0.  ,  0.  ,  0.  ])
>>> _ = rg.add_field("neg_grad_at_link", -lg, at="link")
>>> calc_flux_div_at_node(rg, "neg_grad_at_link")
array([ 0.  ,  0.  ,  0.  ,  0.  ,  0.  ,  1.64,  0.94,  0.  ,  0.  ,
        0.  ,  0.  ,  0.  ])

Notes

Performs a numerical flux divergence operation on nodes.

calc_grad_across_cell_corners(node_values, *args, **kwds)#

grid.calc_grad_across_cell_corners(node_values, [cell_ids], out=None)

Get gradients to diagonally opposite nodes.

Calculate gradient of the value field provided by node_values to the values at diagonally opposite nodes. The returned gradients are ordered as upper-right, upper-left, lower-left and lower-right.

Parameters:
  • node_values (array_like or field name) – Quantity to take the gradient of defined at each node.

  • cell_ids (array_like, optional) – If provided, cell ids to measure gradients. Otherwise, find gradients for all cells.

  • out (array_like, optional) – Alternative output array in which to place the result. Must be of the same shape and buffer length as the expected output.

Returns:

Gradients to each diagonal node.

Return type:

(N, 4) ndarray

Examples

Create a grid with two cells.

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> x = np.array([1.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 1.0, 3.0, 3.0, 3.0, 3.0])

A decrease in quantity to a diagonal node is a negative gradient.

>>> from math import sqrt
>>> grid.calc_grad_across_cell_corners(x) * sqrt(2.0)
array([[ 3.,  3.,  1.,  0.],
       [ 2.,  2., -1.,  0.]])
>>> grid = RasterModelGrid((3, 4), xy_spacing=(4, 3))
>>> grid.calc_grad_across_cell_corners(x)
array([[ 0.6,  0.6,  0.2,  0. ],
       [ 0.4,  0.4, -0.2,  0. ]])
calc_grad_across_cell_faces(node_values, *args, **kwds)#

grid.calc_grad_across_cell_faces(node_values, [cell_ids], out=None)

Get gradients across the faces of a cell.

Calculate gradient of the value field provided by node_values across each of the faces of the cells of a grid. The returned gradients are ordered as right, top, left, and bottom.

Note that the returned gradients are masked to exclude neighbor nodes which are closed. Beneath the mask is the value -1.

Parameters:
  • node_values (array_like or field name) – Quantity to take the gradient of defined at each node.

  • cell_ids (array_like, optional) – If provided, cell ids to measure gradients. Otherwise, find gradients for all cells.

  • out (array_like, optional) – Alternative output array in which to place the result. Must be of the same shape and buffer length as the expected output.

Returns:

Gradients for each face of the cell.

Return type:

(N, 4) Masked ndarray

Examples

Create a grid with two cells.

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> x = np.array([0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 1.0, 3.0, 3.0, 3.0, 3.0])

A decrease in quantity across a face is a negative gradient.

>>> grid.calc_grad_across_cell_faces(x)
masked_array(data =
 [[ 1.  3.  0.  0.]
 [ 0.  2. -1. -1.]],
             mask =
 False,
       fill_value = 1e+20)
>>> grid = RasterModelGrid((3, 4), xy_spacing=(1, 2))
>>> grid.calc_grad_across_cell_faces(x)
masked_array(data =
 [[ 1.   1.5  0.   0. ]
 [ 0.   1.  -1.  -0.5]],
              mask =
 False,
       fill_value = 1e+20)

grid.calc_grad_along_node_links(node_values, [cell_ids], out=None)

Get gradients along links touching a node.

Calculate gradient of the value field provided by node_values across each of the faces of the nodes of a grid. The returned gradients are ordered as right, top, left, and bottom. All returned values follow our standard sign convention, where a link pointing N or E and increasing in value is positive, a link pointing S or W and increasing in value is negative.

Note that the returned gradients are masked to exclude neighbor nodes which are closed. Beneath the mask is the value -1.

Parameters:
  • node_values (array_like or field name) – Quantity to take the gradient of defined at each node.

  • node_ids (array_like, optional) – If provided, node ids to measure gradients. Otherwise, find gradients for all nodes.

  • out (array_like, optional) – Alternative output array in which to place the result. Must be of the same shape and buffer length as the expected output.

Returns:

Gradients for each link of the node. Ordering is E,N,W,S.

Return type:

(N, 4) Masked ndarray

Examples

Create a grid with nine nodes.

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3))
>>> x = np.array([0.0, 0.0, 0.0, 0.0, 1.0, 2.0, 2.0, 2.0, 2.0])

A decrease in quantity across a face is a negative gradient.

>>> grid.calc_grad_along_node_links(x)
masked_array(data =
 [[-- -- -- --]
 [-- 1.0 -- --]
 [-- -- -- --]
 [1.0 -- -- --]
 [1.0 1.0 1.0 1.0]
 [-- -- 1.0 --]
 [-- -- -- --]
 [-- -- -- 1.0]
 [-- -- -- --]],
             mask =
 [[ True  True  True  True]
 [ True False  True  True]
 [ True  True  True  True]
 [False  True  True  True]
 [False False False False]
 [ True  True False  True]
 [ True  True  True  True]
 [ True  True  True False]
 [ True  True  True  True]],
       fill_value = 1e+20)
>>> grid = RasterModelGrid((3, 3), xy_spacing=(4, 2))
>>> grid.calc_grad_along_node_links(x)
masked_array(data =
 [[-- -- -- --]
 [-- 0.5 -- --]
 [-- -- -- --]
 [0.25 -- -- --]
 [0.25 0.5 0.25 0.5]
 [-- -- 0.25 --]
 [-- -- -- --]
 [-- -- -- 0.5]
 [-- -- -- --]],
             mask =
 [[ True  True  True  True]
 [ True False  True  True]
 [ True  True  True  True]
 [False  True  True  True]
 [False False False False]
 [ True  True False  True]
 [ True  True  True  True]
 [ True  True  True False]
 [ True  True  True  True]],
       fill_value = 1e+20)
calc_grad_at_d8(node_values, out=None)#

Calculate gradients over all diagonals and links.

Parameters:
  • node_values (array_like or field name) – Values at nodes.

  • out (ndarray, optional) – Buffer to hold result. If None, create a new array.

Examples

>>> from landlab import RasterModelGrid
>>> import numpy as np
>>> grid = RasterModelGrid((3, 4), xy_spacing=(4, 3))
>>> z = [
...     [60.0, 60.0, 60.0, 60.0],
...     [60.0, 60.0, 0.0, 0.0],
...     [60.0, 0.0, 0.0, 0.0],
... ]
>>> grid.calc_grad_at_d8(z)
array([  0.,   0.,   0.,   0.,   0., -20., -20.,   0., -15.,   0.,   0.,
       -20.,   0.,   0., -15.,   0.,   0.,   0.,   0., -12.,   0., -12.,
       -12., -12.,   0., -12.,   0.,   0.,   0.])
calc_grad_at_diagonal(node_values, out=None)#

Calculate gradients over all diagonals.

Parameters:
  • node_values (array_like or field name) – Values at nodes.

  • out (ndarray, optional) – Buffer to hold result. If None, create a new array.

Examples

>>> from landlab import RasterModelGrid
>>> import numpy as np
>>> grid = RasterModelGrid((3, 4), xy_spacing=(4, 3))
>>> z = [
...     [5.0, 5.0, 5.0, 5.0],
...     [5.0, 5.0, 0.0, 0.0],
...     [5.0, 0.0, 0.0, 0.0],
... ]
>>> grid.calc_grad_at_diagonal(z)
array([ 0.,  0., -1.,  0., -1., -1., -1.,  0., -1.,  0.,  0.,  0.])

Calculate gradients in node_values at links.

Parameters:
  • value_at_node (array_like or field name) – Values at nodes.

  • out (ndarray, optional) – Buffer to hold result. If None, create a new array.

Returns:

Gradients of the nodes values for each link.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3))
>>> node_values = [0.0, 0.0, 0.0, 1.0, 3.0, 1.0, 2.0, 2.0, 2.0]
>>> grid.calc_grad_at_link(node_values)
array([ 0.,  0.,  1.,  3.,  1.,  2., -2.,  1., -1.,  1.,  0.,  0.])
>>> out = np.empty(grid.number_of_links, dtype=float)
>>> rtn = grid.calc_grad_at_link(node_values, out=out)
>>> rtn is out
True
>>> out
array([ 0.,  0.,  1.,  3.,  1.,  2., -2.,  1., -1.,  1.,  0.,  0.])
>>> grid = RasterModelGrid((3, 3), xy_spacing=(2, 1))
>>> grid.calc_grad_at_link(node_values)
array([ 0.,  0.,  1.,  3.,  1.,  1., -1.,  1., -1.,  1.,  0.,  0.])
>>> _ = grid.add_field("elevation", node_values, at="node")
>>> grid.calc_grad_at_link("elevation")
array([ 0.,  0.,  1.,  3.,  1.,  1., -1.,  1., -1.,  1.,  0.,  0.])
calc_grad_at_patch(elevs='topographic__elevation', ignore_closed_nodes=True, subtriangle_unit_normals=None, slope_magnitude=None)#

Calculate the components of the gradient of each raster patch.

Returns the mean gradient of the four possible patch subtriangles, in radians.

If ignore_closed_nodes is True, closed nodes do not affect gradient calculations. If more than one closed node is present in a patch, the patch gradients in both x and y directions are set to zero.

Parameters:
  • elevs (str or ndarray, optional) – Field name or array of node values.

  • ignore_closed_nodes (bool) – If True, do not incorporate values at closed nodes into the calc.

  • subtriangle_unit_normals (tuple of 4 (npatches, 3) arrays (optional)) – The unit normal vectors for the four subtriangles of each patch, if already known. Order is TR, TL, BL, BR.

  • slope_magnitude (array with size num_patches (optional)) – The mean slope of each patch, if already known. Units must be the same as provided here!

Returns:

gradient_tuple – Len-2 tuple of arrays giving components of gradient in the x and y directions, in the units of radians.

Return type:

(x_component_at_patch, y_component_at_patch)

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> z = mg.node_y
>>> (x_grad, y_grad) = mg.calc_grad_at_patch(elevs=z)
>>> np.allclose(y_grad, np.pi / 4.0)
True
>>> np.allclose(x_grad, 0.0)
True
>>> z = mg.node_x.copy()
>>> mg.set_closed_boundaries_at_grid_edges(True, True, True, True)
>>> mg.status_at_node[11] = mg.BC_NODE_IS_CLOSED
>>> mg.status_at_node[[9, 2]] = mg.BC_NODE_IS_FIXED_VALUE
>>> z[11] = 100.0  # this should get ignored now
>>> z[9] = 2.0  # this should be felt by patch 7 only
>>> z[2] = 1.0  # should be felt by patches 1 and 2
>>> xgrad, ygrad = mg.calc_grad_at_patch(elevs=z, ignore_closed_nodes=True)
>>> (xgrad.reshape((3, 4)) * 4.0 / np.pi)[1, 1:]
array([ 1.,  1., -1.])
>>> np.allclose(ygrad[1:3], xgrad[1:3])
True
calc_hillshade_at_node(alt=45.0, az=315.0, slp=None, asp=None, unit='degrees', elevs='topographic__elevation')#

Get array of hillshade.

Code author: Katy Barnhart <katherine.barnhart@colorado.edu>

Parameters:
  • alt (float) – Sun altitude (from horizon) - defaults to 45 degrees

  • az (float) – Sun azimuth (CW from north) - defaults to 315 degrees

  • slp (float) – slope of cells at surface - optional

  • asp (float) – aspect of cells at surface (from north) - optional (with slp)

  • unit (string) –

    ‘degrees’ (default) or ‘radians’ - only needed if slp and asp

    are not provided

  • specified (If slp and asp are both not) –

  • as ('elevs' must be provided) –

  • an (a grid field name (defaults to 'topographic__elevation') or) –

  • case (nnodes-long array of elevation values. In this) –

  • will (the method) –

  • hillshade (calculate local slopes and aspects internally as part of the) –

  • production.

Returns:

Hillshade at each cell.

Return type:

ndarray of float

Notes

code taken from GeospatialPython.com example from December 14th, 2014 DEJH found what looked like minor sign problems, and adjusted to follow the ArcGIS algorithm <http://help.arcgis.com/en/arcgisdesktop/10.0/help/index.html#/How_Hillshade_works/009z000000z2000000/>.

Remember when plotting that bright areas have high values. cmap=’Greys’ will give an apparently inverted color scheme. cmap=’gray’ has white associated with the high values, so is recommended for plotting.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((5, 5), xy_spacing=1.0)
>>> z = mg.x_of_node * np.tan(60.0 * np.pi / 180.0)
>>> mg.calc_hillshade_at_node(elevs=z, alt=30.0, az=210.0)
array([ 0.625,  0.625,  0.625,  0.625,  0.625,  0.625,  0.625,  0.625,
        0.625,  0.625,  0.625,  0.625,  0.625,  0.625,  0.625,  0.625,
        0.625,  0.625,  0.625,  0.625,  0.625,  0.625,  0.625,  0.625,
        0.625])
calc_net_flux_at_node(unit_flux_at_links, out=None)#

Calculate net link fluxes at nodes.

Given a flux per unit width along each link in the grid, calculate the net outflux (or influx, if negative) at each node. Fluxes are treated as zero for links that have no faces, and net fluxes are treated as zero for nodes that have no cell.

Parameters:
  • unit_flux_at_links (ndarray or field name) – Flux per unit width associated with links.

  • out (ndarray, optional) – Buffer to hold the result.

Returns:

Net flux at nodes.

Return type:

ndarray (x number of cells)

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
>>> lg = rg.calc_grad_at_link(z)  # there are 17 links
>>> lg
array([ 0. ,  0. ,  0. ,  0. ,  5. ,  3.6,  0. ,  5. , -1.4, -3.6,  0. ,
       -5. , -3.6,  0. ,  0. ,  0. ,  0. ])
>>> calc_net_flux_at_node(rg, -lg)
array([   0.,    0.,    0.,    0.,    0.,  164.,   94.,    0.,    0.,
          0.,    0.,    0.])
>>> rg.set_status_at_node_on_edges(right=rg.BC_NODE_IS_CLOSED)
>>> rg.set_status_at_node_on_edges(top=rg.BC_NODE_IS_CLOSED)
>>> unit_flux_at_links = np.zeros(rg.number_of_links)
>>> unit_flux_at_links[rg.active_links] = -lg[rg.active_links]
>>> nlfn = calc_net_flux_at_node(rg, unit_flux_at_links)
>>> np.round(nlfn)
array([   0.,    0.,    0.,    0.,    0.,  114.,   22.,    0.,    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
>>> lg = hg.calc_grad_at_link(z)  # there are ? links
>>> lg
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. ])
>>> nlfn = calc_net_flux_at_node(hg, -lg)
>>> np.round(nlfn)
array([   0.,    0.,    0.,    0.,  152.,   96.,    0.,    0.,    0.,    0.])

Notes

This is essentially a line integral for the fluxes along the boundaries of each cell. Hence, the resulting output has dimensions of total flux (so, if the unit flux happens to be mass per time per face width, the output will be in mass per unit time). Because a line integral is undefined where there are no cells (i.e., perimeter nodes), the result is given as zeros for these nodes. The current algorithm uses fancy indexing (calling _calc_net_face_flux_at_cells) and could probably be made faster.

calc_slope_at_cell_subtriangles(elevs='topographic__elevation', subtriangle_unit_normals=None)#

Calculate the slope (positive magnitude of gradient) at each of the eight cell subtriangles.

Parameters:
  • elevs (str or ndarray, optional) – Field name or array of node values.

  • subtriangle_unit_normals (tuple of 8 (ncells, 3) arrays (optional)) – The unit normal vectors for the eight subtriangles of each cell, if already known. Order is from north of east, counter clockwise to south of east (East North East, North North East, North North West, West North West, West South West, South South West, South South East, East South East).

Returns:

each a length num-cells array Len-8 tuple of the slopes (positive gradient magnitude) of each of the eight cell subtriangles, in radians. Order is from north of east, counter clockwise to south of east (East North East, North North East, North North West, West North West, West South West, South South West, South South East, East South East).

Return type:

(s_ENE, s_NNE, s_NNW, s_WNW, s_WSW, s_SSW, s_SSE, s_ESE)

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> z = np.array(
...     [np.sqrt(3.0), 0.0, 4.0 / 3.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0 / np.sqrt(3.0)]
... )
>>> eight_tris = mg.calc_unit_normals_at_cell_subtriangles(z)
>>> S = mg.calc_slope_at_cell_subtriangles(z, eight_tris)
>>> S0 = mg.calc_slope_at_cell_subtriangles(z)
>>> np.allclose(S, S0)
True
>>> type(S) is tuple
True
>>> len(S)
8
>>> len(S[0]) == mg.number_of_cells
True
>>> np.allclose(S[0], S[1])
True
>>> np.allclose(S[2], S[3])
True
>>> np.allclose(S[4], S[5])
True
>>> np.allclose(S[6], S[7])
True
>>> np.allclose(np.rad2deg(S[0])[0], 30.0)
True
>>> np.allclose(np.rad2deg(S[2])[0], 45.0)
True
>>> np.allclose(np.rad2deg(S[4])[0], 60.0)
True
>>> np.allclose(np.cos(S[6])[0], 3.0 / 5.0)
True
calc_slope_at_node(elevs='topographic__elevation', method='patch_mean', ignore_closed_nodes=True, return_components=False)#

Array of slopes at nodes, averaged over neighboring patches.

Produces a value for node slope (i.e., mean gradient magnitude) at each node in a manner analogous to a GIS-style slope map. If method==’patch_mean’, it averages the gradient on each of the patches surrounding the node; if method==’Horn’, it returns the resolved slope direction. Directional information can still be returned through use of the return_components keyword. All values are returned in radians, including the components; take the tan to recover the rise/run.

Note that under these definitions, it is not always true that:

mag, cmp = mg.calc_slope_at_node(z)
mag**2 == cmp[0]**2 + cmp[1]**2  # only if method=='Horn'

If ignore_closed_nodes is False, all proximal elevation values will be used in the calculation. If True, only unclosed nodes are used.

This is a verion of this code specialized for a raster. It subdivides the four square patches around each node into subtriangles, in order to ensure more correct solutions that incorporate equally weighted information from all surrounding nodes on rough surfaces.

Parameters:
  • elevs (str or ndarray, optional) – Field name or array of node values.

  • method ({'patch_mean', 'Horn'}) – Controls the slope algorithm. Current options are ‘patch_mean’, which takes the mean slope of each pf the four neighboring square patches, and ‘Horn’, which is the standard ArcGIS slope algorithm. These produce very similar solutions; the Horn method gives a vector mean and the patch_mean gives a scalar mean.

  • ignore_closed_nodes (bool) – If True, do not incorporate values at closed nodes into the calc.

  • return_components (bool) – If True, return a tuple, (array_of_magnitude, (array_of_slope_x_radians, array_of_slope_y_radians)). If false, return an array of floats of the slope magnitude.

Returns:

If return_components, returns (array_of_magnitude, (array_of_slope_x_radians, array_of_slope_y_radians)). If not return_components, returns an array of slope magnitudes.

Return type:

float array or length-2 tuple of float arrays

Examples

>>> import numpy as np
>>> from landlab import RadialModelGrid, RasterModelGrid
>>> mg = RasterModelGrid((5, 5))
>>> z = mg.node_x
>>> slopes = mg.calc_slope_at_node(elevs=z)
>>> np.allclose(slopes, np.pi / 4.0)
True
>>> mg = RasterModelGrid((4, 5), xy_spacing=2.0)
>>> z = -mg.node_y
>>> slope_mag, cmp = mg.calc_slope_at_node(elevs=z, return_components=True)
>>> np.allclose(slope_mag, np.pi / 4.0)
True
>>> np.allclose(cmp[0], 0.0)
True
>>> np.allclose(cmp[1], -np.pi / 4.0)
True
>>> mg = RasterModelGrid((4, 4))
>>> z = mg.node_x**2 + mg.node_y**2
>>> slopes, cmp = mg.calc_slope_at_node(z, return_components=True)
>>> slopes
array([ 0.95531662,  1.10991779,  1.32082849,  1.37713803,  1.10991779,
        1.20591837,  1.3454815 ,  1.38904403,  1.32082849,  1.3454815 ,
        1.39288142,  1.41562833,  1.37713803,  1.38904403,  1.41562833,
        1.43030663])

Check radial symmetry.

>>> np.allclose(cmp[0].reshape((4, 4))[:, 0], cmp[1].reshape((4, 4))[0, :])
True
calc_slope_at_patch(elevs='topographic__elevation', ignore_closed_nodes=True, subtriangle_unit_normals=None)#

Calculate the slope (positive magnitude of gradient) at raster patches.

Returns the mean of the slopes of the four possible patch subtriangles.

If ignore_closed_nodes is True, closed nodes do not affect slope calculations. If more than one closed node is present in a patch, the patch slope is set to zero.

Parameters:
  • elevs (str or ndarray, optional) – Field name or array of node values.

  • ignore_closed_nodes (bool) – If True, do not incorporate values at closed nodes into the calc.

  • subtriangle_unit_normals (tuple of 4 (npatches, 3) arrays (optional)) – The unit normal vectors for the four subtriangles of each patch, if already known. Order is TR, TL, BL, BR.

Returns:

slopes_at_patch – The slope (positive gradient magnitude) of each patch, in radians.

Return type:

n_patches-long array

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> z = mg.node_x
>>> S = mg.calc_slope_at_patch(elevs=z)
>>> S.size == mg.number_of_patches
True
>>> np.allclose(S, np.pi / 4.0)
True
>>> z = mg.node_y**2
>>> mg.calc_slope_at_patch(elevs=z).reshape((3, 4))
array([[ 0.78539816,  0.78539816,  0.78539816,  0.78539816],
       [ 1.24904577,  1.24904577,  1.24904577,  1.24904577],
       [ 1.37340077,  1.37340077,  1.37340077,  1.37340077]])
>>> z = mg.node_x.copy()
>>> mg.set_closed_boundaries_at_grid_edges(True, True, True, True)
>>> mg.status_at_node[11] = mg.BC_NODE_IS_CLOSED
>>> mg.status_at_node[9] = mg.BC_NODE_IS_FIXED_VALUE
>>> z[11] = 100.0  # this should get ignored now
>>> z[9] = 2.0  # this should be felt by patch 7 only
>>> mg.calc_slope_at_patch(elevs=z, ignore_closed_nodes=True).reshape(
...     (3, 4)
... ) * 4.0 / np.pi
array([[ 0.,  0.,  0.,  0.],
       [ 0.,  1.,  1.,  1.],
       [ 0.,  0.,  0.,  0.]])
calc_unit_normal_at_patch(elevs='topographic__elevation')[source]#

Calculate and return the unit normal vector <a, b, c> to a patch.

This method is not defined on a raster, as there is no unique unit normal for a square patch. Use _calc_unit_normals_to_patch_subtriangles instead.

calc_unit_normals_at_cell_subtriangles(elevs='topographic__elevation')#

Calculate unit normals on a cell.

Calculate the eight unit normal vectors <a, b, c> to the eight subtriangles of a four-cornered (raster) cell.

Parameters:

elevs (str or ndarray, optional) – Field name or array of node values.

Returns:

each a num-cells x length-3 array Len-8 tuple of the eight unit normal vectors <a, b, c> for the eight subtriangles in the cell. Order is from north of east, counter clockwise to south of east (East North East, North North East, North North West, West North West, West South West, South South West, South South East, East South East).

Return type:

(n_ENE, n_NNE, n_NNW, n_WNW, n_WSW, n_SSW, n_SSE, n_ESE)

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> z = mg.node_x**2
>>> eight_tris = mg.calc_unit_normals_at_cell_subtriangles(z)
>>> type(eight_tris) is tuple
True
>>> len(eight_tris)
8
>>> eight_tris[0].shape == (mg.number_of_cells, 3)
True
>>> eight_tris
(array([[-0.9486833 ,  0.        ,  0.31622777]]),
 array([[-0.9486833 ,  0.        ,  0.31622777]]),
 array([[-0.70710678,  0.        ,  0.70710678]]),
 array([[-0.70710678,  0.        ,  0.70710678]]),
 array([[-0.70710678,  0.        ,  0.70710678]]),
 array([[-0.70710678,  0.        ,  0.70710678]]),
 array([[-0.9486833 ,  0.        ,  0.31622777]]),
 array([[-0.9486833 ,  0.        ,  0.31622777]]))
calc_unit_normals_at_patch_subtriangles(elevs='topographic__elevation')#

Calculate unit normals on a patch.

Calculate the four unit normal vectors <a, b, c> to the four possible subtriangles of a four-cornered (raster) patch.

Parameters:

elevs (str or ndarray, optional) – Field name or array of node values.

Returns:

(n_TR, n_TL, n_BL, n_BR) – Len-4 tuple of the four unit normal vectors <a, b, c> for the four possible subtriangles in the patch. Order is (topright, topleft, bottomleft, bottomright).

Return type:

each a num-patches x length-3 array

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> z = mg.node_x**2
>>> four_tris = mg.calc_unit_normals_at_patch_subtriangles(z)
>>> type(four_tris) is tuple
True
>>> len(four_tris)
4
>>> np.allclose(four_tris[0], four_tris[1])
True
>>> np.allclose(four_tris[2], four_tris[3])
True
>>> np.allclose(four_tris[0], four_tris[2])
True
>>> np.allclose(np.square(four_tris[0]).sum(axis=1), 1.0)
True
>>> four_tris[0]
array([[-0.70710678,  0.        ,  0.70710678],
       [-0.9486833 ,  0.        ,  0.31622777],
       [-0.98058068,  0.        ,  0.19611614],
       [-0.98994949,  0.        ,  0.14142136],
       [-0.70710678,  0.        ,  0.70710678],
       [-0.9486833 ,  0.        ,  0.31622777],
       [-0.98058068,  0.        ,  0.19611614],
       [-0.98994949,  0.        ,  0.14142136],
       [-0.70710678,  0.        ,  0.70710678],
       [-0.9486833 ,  0.        ,  0.31622777],
       [-0.98058068,  0.        ,  0.19611614],
       [-0.98994949,  0.        ,  0.14142136]])
calculate_slope_aspect_at_nodes_burrough(ids=None, vals='Elevation')[source]#

Calculate topographic slope.

Calculates the local topographic slope (i.e., the down-dip slope, and presented as positive), and the aspect (dip direction in degrees clockwise from north), at the given nodes, ids. All ids must be of core nodes. This method uses Burrough’s 1998 Pg. 190 method similar to the methods used by ArcMap to calculate slope and aspect.

If ids is not provided, the slope will be returned for nodes at all cells.

vals is either the name of an existing grid field from which to draw topographic data, or an array of values to use. If an array of values is passed, it must be nnodes long. If vals is not provided, this method will default to trying to use the field ‘Elevation’.

Returns:

(slope, aspect)slope, a len(ids) array of slopes at each node provided. aspect, a len(ids) array of aspects at each node provided.

Return type:

tuple of float

Examples

>>> import pytest
>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4), xy_spacing=(4, 4))
>>> z = np.array([0.0, 0.0, 0.0, 0.0, 3.0, 3.0, 3.0, 3, 6.0, 6.0, 6.0, 6.0])
>>> slope, aspect = grid.calculate_slope_aspect_at_nodes_burrough(vals=z)
>>> np.tan(slope)
array([ 0.75,  0.75])
>>> np.degrees(aspect)
array([ 180.,  180.])

We recommend using the following functions instead of this one: - calc_slope_at_node - calc_aspect_at_node Notice that calc_slope_at_node and :py:meth:`~landlab.grid.RasterModelGrid.calc_aspect_at_node return values for all nodes, not just core nodes. In addition, :py:meth:`~landlab.grid.RasterModelGrid.calc_aspect_at_node returns compass-style angles in degrees.

>>> np.tan(grid.calc_slope_at_node(elevs=z)[grid.core_nodes])
array([ 0.75,  0.75])
>>> grid.calc_aspect_at_node(elevs=z)[grid.core_nodes]
array([ 180.,  180.])
property cell_area_at_node#

Cell areas in a nnodes-long array.

Zeros are entered at all perimeter nodes, which lack cells.

Returns:

Cell areas as an n_nodes-long array.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5), xy_spacing=(3, 4))
>>> grid.status_at_node[7] = grid.BC_NODE_IS_CLOSED
>>> grid.cell_area_at_node
array([  0.,   0.,   0.,   0.,   0.,
         0.,  12.,  12.,  12.,   0.,
         0.,  12.,  12.,  12.,   0.,
         0.,   0.,   0.,   0.,   0.])
property cell_at_node#
property cell_grid_shape#

Get the shape of the cellular grid (grid with only cells).

Returns:

shape – The shape of the cellular grid as number of cell rows and cell columns.

Return type:

tuple of ints

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.cell_grid_shape
(1, 2)
cell_vector_to_raster(u, flip_vertically=False)[source]#

Unravel a 1D array.

Converts cell vector u to a 2D array and returns it, so that it can be plotted, output, etc.

If the optional argument flip_vertically=True, the function returns an array that has the rows in reverse order, for use in plot commands (such as the image display functions) that put the (0,0) axis at the top left instead of the bottom left.

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5))
>>> u = rmg.zeros(centering="cell")
>>> u = u + range(0, len(u))
>>> u
array([ 0.,  1.,  2.,  3.,  4.,  5.])
>>> ur = rmg.cell_vector_to_raster(u)
>>> ur
array([[ 0.,  1.,  2.],
       [ 3.,  4.,  5.]])
>>> ur = rmg.cell_vector_to_raster(u, flip_vertically=True)
>>> ur
array([[ 3.,  4.,  5.],
       [ 0.,  1.,  2.]])
property cells_at_corner#
property cells_at_corners_of_grid#

Get array of cells in cellular grid (grid with only cells) corners.

Return the IDs to the corner cells of the cellular grid, sorted by ID.

Returns:

Array of corner node IDs.

Return type:

(4, ) ndarray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.cells_at_corners_of_grid
array([0, 2, 3, 5])
property cells_at_face#
property cells_present_at_corner#

A boolean array, False where a cell has a closed corner or is

property cells_present_at_face#

A boolean array, False where a cell has a closed corner or is

property closed_boundary_corners#

Get array of closed boundary corners.

property closed_boundary_nodes#

Get array of closed boundary nodes.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.closed_boundary_nodes
array([15, 16, 17, 18, 19])
property core_cells#

Get array of core cells.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))

Initially all of the cells are “core”.

>>> grid.core_cells
array([0, 1, 2,
       3, 4, 5])

Setting a node to closed causes its cell to no longer be core.

>>> grid.status_at_node[8] = grid.BC_NODE_IS_CLOSED
>>> grid.core_cells
array([0, 1, 3, 4, 5])
property core_corners#

Get array of core corners.

See also

core_nodes

property core_nodes#

Get array of core nodes.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.core_nodes
array([ 6,  7,  8, 11, 12, 13])
property core_patches#

Get array of core patches.

See also

core_cells

property corner_at_core_patch#

Get array of corners associated with core patches.

property corner_at_face_head#

Get corners at face head.

property corner_at_face_tail#

Get corners at face tail.

property corner_at_patch#
property corner_nodes#
property corner_x#
property corner_y#
property corners#

A shaped array of corner ids.

See also

nodes

property corners_at_bottom_edge#
property corners_at_cell#

Get the corners that define a cell.

See also

nodes_at_patch

property corners_at_d8#
property corners_at_diagonal#

Nodes at diagonal tail and head.

property corners_at_face#

Get corners at either end of faces.

See also

nodes_at_link

property corners_at_left_edge#
property corners_at_nodes_of_grid#

Nodes at nodes of grid.

property corners_at_right_edge#
property corners_at_top_edge#
property d8_adjacent_corners_at_corner#
property d8_adjacent_nodes_at_node#
property d8_dirs_at_corner#
property d8_dirs_at_node#
property d8_status_at_corner#
property d8_status_at_node#
property d8s_at_corner#

Links and diagonals attached to corners.

See also

d8s_at_node

property d8s_at_node#

Links and diagonals attached to nodes.

Returns:

Links and diagonals at each node.

Return type:

ndarray of int, shape (n_nodes, 8)

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.d8s_at_node.shape == (grid.number_of_nodes, 8)
True
>>> grid.d8s_at_node
array([[ 0,  3, -1, -1, 17, -1, -1, -1],
       [ 1,  4,  0, -1, 19, 18, -1, -1],
       [ 2,  5,  1, -1, 21, 20, -1, -1],
       [-1,  6,  2, -1, -1, 22, -1, -1],
       [ 7, 10, -1,  3, 23, -1, -1, 18],
       [ 8, 11,  7,  4, 25, 24, 17, 20],
       [ 9, 12,  8,  5, 27, 26, 19, 22],
       [-1, 13,  9,  6, -1, 28, 21, -1],
       [14, -1, -1, 10, -1, -1, -1, 24],
       [15, -1, 14, 11, -1, -1, 23, 26],
       [16, -1, 15, 12, -1, -1, 25, 28],
       [-1, -1, 16, 13, -1, -1, 27, -1]])
>>> np.all(grid.d8s_at_node[:, :4] == grid.links_at_node)
True
>>> diagonals_at_node = grid.d8s_at_node[:, 4:] - grid.number_of_links
>>> diagonals_at_node[grid.d8s_at_node[:, 4:] == -1] = -1
>>> np.all(diagonals_at_node == grid.diagonals_at_node)
True
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 diagonal_adjacent_corners_at_corner#

Get adjacent corners along diagonals.

property diagonal_adjacent_nodes_at_node#

Get adjacent nodes along diagonals.

Order is the landlab standard, counterclockwise starting from east.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.diagonal_adjacent_nodes_at_node
array([[ 4, -1, -1, -1], [ 5,  3, -1, -1], [-1,  4, -1, -1],
       [ 7, -1, -1,  1], [ 8,  6,  0,  2], [-1,  7,  1, -1],
       [10, -1, -1,  4], [11,  9,  3,  5], [-1, 10,  4, -1],
       [-1, -1, -1,  7], [-1, -1,  6,  8], [-1, -1,  7, -1]])
property diagonal_dirs_at_corner#

Directions of diagonals attached to corners.

property diagonal_dirs_at_node#

Directions of diagonals attached to nodes.

Array of diagonal directions relative to each node. If the diagonal is directed away from the node the direction is -1, if the diagonal is directed toward the node its direction is 1. Each node is assumed to have exactly four diagonals attached to it. However, perimeter nodes will have fewer diagonals (corner nodes will only have one diagonal and edge nodes two). In such cases, placeholders of 0 are used.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.diagonal_dirs_at_node.shape == (grid.number_of_nodes, 4)
True
>>> grid.diagonal_dirs_at_node
array([[-1,  0,  0,  0], [-1, -1,  0,  0], [ 0, -1,  0,  0],
       [-1,  0,  0,  1], [-1, -1,  1,  1], [ 0, -1,  1,  0],
       [-1,  0,  0,  1], [-1, -1,  1,  1], [ 0, -1,  1,  0],
       [ 0,  0,  0,  1], [ 0,  0,  1,  1], [ 0,  0,  1,  0]],
       dtype=int8)

The lower-right corner node only has one diagonal.

>>> grid.diagonal_dirs_at_node[2]
array([ 0, -1,  0,  0], dtype=int8)

A node on one of the edges has two diagonals.

>>> grid.diagonal_dirs_at_node[3]
array([-1,  0,  0,  1], dtype=int8)
property diagonal_status_at_corner#
property diagonal_status_at_node#
property diagonals_at_corner#

Diagonals attached to corners.

property diagonals_at_node#

Diagonals attached to nodes.

Returns:

Diagonals at each node.

Return type:

ndarray of int, shape (n_nodes, 4)

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.diagonals_at_node.shape == (grid.number_of_nodes, 4)
True
>>> grid.diagonals_at_node
array([[ 0, -1, -1, -1], [ 2,  1, -1, -1], [-1,  3, -1, -1],
       [ 4, -1, -1,  1], [ 6,  5,  0,  3], [-1,  7,  2, -1],
       [ 8, -1, -1,  5], [10,  9,  4,  7], [-1, 11,  6, -1],
       [-1, -1, -1,  9], [-1, -1,  8, 11], [-1, -1, 10, -1]])

Return an (nnodes, X) shape array of link IDs of which links are downwind of each node, according to values (array or field).

X is the maximum downwind links at any node. Nodes with fewer downwind links than this have additional slots filled with bad_index. Links are ordered anticlockwise from east.

Parameters:
  • values (str or array) – Name of variable field defined at links, or array of values at links.

  • bad_index (int) – Index to place in array indicating no link.

Returns:

Array of upwind link IDs

Return type:

ndarray

Examples

>>> import numpy as np
>>> 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,
... ]
>>> rmg.downwind_links_at_node("grad", bad_index=rmg.BAD_INDEX)
array([[ 0,  3],
       [ 1,  4],
       [ 2,  5],
       [ 6, -1],
       [ 7, 10],
       [ 8, 11],
       [ 9, 12],
       [13, -1],
       [14, -1],
       [15, -1],
       [16, -1],
       [-1, -1]])
property ds#
property dual#
property dx#
property dy#
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"))
[]
property event_layers#

EventLayers for each cell.

property extent#

Extent of the grid in the y and x-dimensions.

Return the y and x-dimension of the grid. Because boundary nodes don’t have cells, the dimension of the grid is ((num_rows - 1) * dy, (num_columns - 1) * dx), not (num_rows * dy, num_cols * dx).

Returns:

(y_extent, x_extent) – Length of the grid in the y and x-dimensions.

Return type:

tuple of float

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.extent
(3.0, 4.0)
>>> grid = RasterModelGrid((4, 5), xy_spacing=2.0)
>>> grid.extent
(6.0, 8.0)
>>> grid = RasterModelGrid((4, 5), xy_spacing=(3, 2))
>>> grid.extent
(6.0, 12.0)
property face_dirs_at_corner#
property face_status_at_corner#
property faces_at_cell#

Get the faces that define a cell.

See also

links_at_patch

property faces_at_corner#
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
fields(include='*', exclude=None)#

List of fields held by the grid.

The returned field names are returned as their canonical names. That is, as a string of the for "at_<location>:<name>". This allows for fields with the same name to be defined at different grid locations. You could have, for example, a variable “elevation” defined at both nodes and links.

Both the include and exclude patterns are glob-style expressions, not regular expressions. If either include or exclude are lists, then the patterns are matched using an “or”.

The include filters are applied before the exclude filters.

Parameters:
  • include (str, or iterable of str, optional) – Glob-style pattern for field names to include.

  • exclude (str, or iterable of str, optional) – Glob-style pattern for field names to exclude.

Returns:

Filtered set of canonical field names held by the grid

Return type:

set

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))

Add some fields to the grid. Notice that we have defined a field named “elevation” at both nodes and links.

>>> _ = grid.add_full("elevation", 3.0, at="node")
>>> _ = grid.add_full("elevation", 4.0, at="link")
>>> _ = grid.add_full("temperature", 5.0, at="node")
>>> sorted(grid.fields())
['at_link:elevation', 'at_node:elevation', 'at_node:temperature']
>>> sorted(grid.fields(include="at_node*"))
['at_node:elevation', 'at_node:temperature']
>>> sorted(grid.fields(include="at_node*", exclude="*temp*"))
['at_node:elevation']

Fields can also be defined at layers. In the following example we’ve filtered the results to just return the layer fields.

>>> grid.event_layers.add(1.0, rho=0.5)
>>> sorted(grid.fields(include="at_layer*"))
['at_layer:rho']

If a list, the fields are matched with an “or”.

>>> sorted(grid.fields(include=["at_node*", "*elevation*"]))
['at_link:elevation', 'at_node:elevation', 'at_node:temperature']
find_nearest_node(coords, mode='raise')[source]#

Node nearest a point.

Find the index to the node nearest the given x, y coordinates. Coordinates are provided as numpy arrays in the coords tuple.

Use the mode keyword to specify what to do if the given coordinates are out-of-bounds. See np.ravel_multi_index for a description of possible values for mode. Note that a coordinate is out-of-bounds if it is beyond one half the node spacing from the exterior nodes.

Returns the indices of the nodes nearest the given coordinates.

Parameters:
  • coords (tuple of array-like) – Coordinates of points.

  • mode ({'raise', 'wrap', 'clip'}, optional) – What to do if a point is off the grid.

Returns:

IDs of the nearest nodes.

Return type:

array-like

Notes

For coordinates that are equidistant to two or more nodes, see the rounding rules for numpy.around.

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5))
>>> rmg.find_nearest_node([0.2, 0.2])
0
>>> rmg.find_nearest_node((np.array([1.6, 3.6]), np.array([2.3, 0.7])))
array([12,  9])
>>> rmg.find_nearest_node((-0.4999, 1.0))
5
property fixed_faces#

Get array of fixed faces.

See also

fixed_links

property fixed_gradient_boundary_corner_anchor_corner#

Returns the corner at the other end of the fixed face for a fixed

property fixed_gradient_boundary_corner_fixed_face#

An array of the fixed_faces connected to fixed gradient boundary

property fixed_gradient_boundary_corners#

Get array of fixed gradient boundary corners.

property fixed_gradient_boundary_node_anchor_node#

Returns the node at the other end of the fixed link for a fixed gradient boundary node.

Degenerate NodeStatus.FIXED_GRADIENT nodes (e.g., corners) are handled as in fixed_gradient_boundary_node_fixed_link, by pointing to a neighboring NodeStatus.FIXED_GRADIENT node.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> leftedge = grid.nodes_at_left_edge
>>> grid.status_at_node[leftedge] = grid.BC_NODE_IS_FIXED_GRADIENT
>>> grid.fixed_gradient_boundary_nodes
array([0, 4, 8])
>>> grid.fixed_gradient_boundary_node_fixed_link
array([ 3,  7, 10])
>>> grid.fixed_gradient_boundary_node_anchor_node
array([4, 5, 4])

An array of the fixed_links connected to fixed gradient boundary nodes.

Note that on a raster, some nodes (notably the corners) can be NodeStatus.FIXED_GRADIENT, but not have a true LinkStatus.FIXED neighboring link. In such cases, the link returned will be a closed link joining the corner node to a neighboring NodeStatus.FIXED_GRADIENT node (see example).

An AssertionError will be raised if for some reason a NodeStatus.FIXED_GRADIENT node exists which has neither a NodeStatus.FIXED_GRADIENT neighbor, or a LinkStatus.FIXED.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> leftedge = grid.nodes_at_left_edge
>>> grid.status_at_node[leftedge] = grid.BC_NODE_IS_FIXED_GRADIENT
>>> grid.fixed_gradient_boundary_nodes
array([0, 4, 8])
>>> grid.fixed_gradient_boundary_node_fixed_link
array([ 3,  7, 10])
property fixed_gradient_boundary_nodes#

Get array of fixed gradient boundary nodes.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_FIXED_GRADIENT
>>> mg.fixed_gradient_boundary_nodes
array([15, 16, 17, 18, 19])

Get array of fixed links.

Examples

>>> from landlab import NodeStatus, RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.status_at_node.reshape(grid.shape)
array([[1, 1, 1, 1],
       [1, 0, 0, 1],
       [1, 1, 1, 1]], dtype=uint8)
>>> grid.fixed_links.size
0
>>> grid.status_at_node[:4] = NodeStatus.FIXED_GRADIENT
>>> grid.status_at_node.reshape(grid.shape)
array([[2, 2, 2, 2],
       [1, 0, 0, 1],
       [1, 1, 1, 1]], dtype=uint8)
>>> grid.fixed_links
array([4, 5])
property fixed_value_boundary_corners#

Get array of fixed value boundary corners.

property fixed_value_boundary_nodes#

Get array of fixed value boundary nodes.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))

Initially all the perimeter nodes are fixed value boundary.

>>> grid.fixed_value_boundary_nodes
array([ 0,  1,  2,  3,  4, 5,  9, 10, 14, 15, 16, 17, 18, 19])

Set left, right, and bottom edges to closed.

>>> for edge in (
...     grid.nodes_at_left_edge,
...     grid.nodes_at_right_edge,
...     grid.nodes_at_bottom_edge,
... ):
...     grid.status_at_node[edge] = grid.BC_NODE_IS_CLOSED

Now nodes on just the top edge are fixed.

>>> grid.fixed_value_boundary_nodes
array([16, 17, 18])
freeze()#

Freeze the graph by making arrays read-only.

classmethod from_dataset(dataset)[source]#
classmethod from_dict(params)[source]#

Create a RasterModelGrid from a dictionary.

Parameters:

params (dict_like) – Initialization parameters for a RasterModelGrid.

Returns:

A newly-created grid.

Return type:

RasterModelGrid

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid.from_dict({"shape": (3, 4), "bc": {"top": "closed"}})
>>> grid.number_of_nodes
12
classmethod from_file(file_like)#

Create grid from a file-like object.

File to load either as a file-like object, path to an existing file, or the contents of a file as a string.

Parameters:

file_like – File-like object, filepath, or string.

Examples

>>> from io import StringIO
>>> from landlab import RasterModelGrid
>>> filelike = StringIO(
...     '''
... shape:
...     - 3
...     - 4
... xy_spacing: 2
... '''
... )
>>> grid = RasterModelGrid.from_file(filelike)
>>> grid.x_of_node
array([ 0.,  2.,  4.,  6.,  0.,  2.,  4.,  6.,  0.,  2.,  4.,  6.])
>>> grid.y_of_node
array([ 0.,  0.,  0.,  0.,  2.,  2.,  2.,  2.,  4.,  4.,  4.,  4.])
classmethod from_netcdf(fname)#
property frozen#
grid_coords_to_node_id(row, col, **kwds)[source]#

Convert node indices to node ID.

Returns the ID of the node at the specified row and col of the raster grid. Since this is a wrapper for the numpy ravel_multi_index function, the keyword arguments are the same as that function. In addition, row and col can both be either scalars or arrays (of the same length) to get multiple ids.

As with ravel_multi_index use the mode keyword to change the behavior of the method when passed an out-of-range row or col. The default is to raise ValueError (not IndexError, as you might expect).

Note

The syntax assumes that first row and column are 0, so max entry for a mg with 4 rows and 5 cols is row=3, col=4

Parameters:
  • row (array-like) – Row of node.

  • col (array-like) – Column of node.

Returns:

Node IDs.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.grid_coords_to_node_id(2, 3)
13
>>> mg.grid_coords_to_node_id([2, 0], [3, 4])
array([13,  4])
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
property horizontal_faces#
imshow(*args, **kwds)#

Plot a data field.

This is a wrapper for plot.imshow_grid, and can take the same keywords. See that function for full documentation.

Parameters:

values (str, or array-like) – Name of a field or an array of values to plot.

See also

landlab.plot.imshow_grid

LLCATS

GINF

is_point_on_grid(xcoord, ycoord)[source]#

Check if a point is on the grid.

This method takes x, y coordinates and tests whether they lie within the grid. The limits of the grid are taken to be links connecting the boundary nodes. We perform a special test to detect looped boundaries.

Coordinates can be ints or arrays of ints. If arrays, will return an array of the same length of boolean truth values.

Parameters:
  • xcoord (float or array_like) – The point’s x-coordinate.

  • ycoord (float or array_like) – The point’s y-coordinate.

Returns:

True if the point is on the grid. Otherwise, False.

Return type:

bool

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5), xy_spacing=(1, 2))
>>> grid.is_point_on_grid(1, 1)
True
>>> grid.is_point_on_grid(
...     (1, 1, 1),
...     (1, 3.1, 6.1),
... )
array([ True,  True, False], dtype=bool)
>>> grid.is_point_on_grid((-0.1, 0.1, 3.9, 4.1), (1, 1, 1, 1))
array([False,  True,  True, False], dtype=bool)
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']
property length_of_d8#

Length of links and diagonals.

Return the lengths links and diagonals in the grid. Links are listed first and then diagonals.

Returns:

Link lengths.

Return type:

ndarray of float

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3), xy_spacing=(4, 3))
>>> grid.length_of_link
array([ 4.,  4.,  3.,  3.,  3.,  4.,  4.,  3.,  3.,  3.,  4.,  4.])
>>> grid.length_of_d8
array([ 4.,  4.,  3.,  3.,  3.,
        4.,  4.,  3.,  3.,  3.,
        4.,  4.,  5.,  5.,  5.,
        5.,  5.,  5.,  5.,  5.])
property length_of_diagonal#
property length_of_face#

Get the length of faces.

See also

length_of_link

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 a boolean the same shape as links_at_node which flags links which are downwind of the node as True.

link_at_node_is_downwind 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. It then return a boolean array the same shape as links_at_node flagging these links. e.g., for a raster, the returned array will be shape (nnodes, 4).

Parameters:
  • values (str or array) – Name of variable field defined at links, or array of values at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array. Must be correct shape and boolean dtype.

Returns:

Boolean of which links are downwind at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> 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,
... ]
>>> rmg.link_at_node_is_downwind("grad")
array([[ True,  True, False, False],
       [ True,  True, False, False],
       [ True,  True, False, False],
       [False,  True, False, False],
       [ True,  True, False, False],
       [ True,  True, False, False],
       [ True,  True, False, False],
       [False,  True, False, False],
       [ True, False, False, False],
       [ True, False, False, False],
       [ True, False, False, False],
       [False, False, False, False]], dtype=bool)

Return a boolean the same shape as links_at_node which flags links which are upwind of the node as True.

link_at_node_is_upwind 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. It then return a boolean array the same shape as links_at_node flagging these links. e.g., for a raster, the returned array will be shape (nnodes, 4).

Parameters:
  • values (str or array) – Name of variable field defined at links, or array of values at links.

  • out (ndarray, optional) – Buffer to place mapped values into or None to create a new array. Must be correct shape and boolean dtype.

Returns:

Boolean of which links are upwind at nodes.

Return type:

ndarray

Examples

>>> import numpy as np
>>> 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,
... ]
>>> rmg.link_at_node_is_upwind("grad")
array([[False, False, False, False],
       [False, False,  True, False],
       [False, False,  True, False],
       [False, False,  True, False],
       [False, False, False,  True],
       [False, False,  True,  True],
       [False, False,  True,  True],
       [False, False,  True,  True],
       [False, False, False,  True],
       [False, False,  True,  True],
       [False, False,  True,  True],
       [False, False,  True,  True]], dtype=bool)

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)

Return array of IDs of links with given angle.

Examples

>>> from landlab import HexModelGrid
>>> grid = HexModelGrid((3, 3))
>>> grid.link_with_angle(0.0)
array([  0,  1,  8,  9, 10, 17, 18])
>>> grid.link_with_angle(60.0, in_degrees=True)
array([  3,  5,  7, 11, 13, 15])
>>> grid.link_with_angle(2.0944)  # 120 degrees
array([  2,  4,  6, 12, 14, 16])
>>> len(grid.link_with_angle(0.5236))  # no links at 30 deg
0
>>> grid = HexModelGrid((3, 3), orientation="vertical")
>>> grid.link_with_angle(30.0, in_degrees=True)
array([  1,  3,  8, 10, 15, 17])
>>> grid.link_with_angle(1.5708)  # 90 degrees
array([ 2,  5,  6,  9, 12, 13, 16])
>>> grid.link_with_angle(330.0, in_degrees=True)
array([ 0,  4,  7, 11, 14, 18])
>>> len(grid.link_with_angle(60.0, in_degrees=True))  # none at 60 deg
0

Links with a given node status.

Parameters:
  • status_at_tail (NodeStatus, optional) – Status of the link tail node.

  • status_at_head (NodeStatus, optional) – Status of the link head node.

Returns:

Links with the given tail and head node statuses.

Return type:

array of int

Examples

>>> from landlab import RasterModelGrid, NodeStatus
>>> grid = RasterModelGrid((4, 5))
>>> grid.status_at_node[13] = NodeStatus.FIXED_VALUE
>>> grid.status_at_node[2] = NodeStatus.CLOSED
>>> grid.link_with_node_status(
...     status_at_tail=NodeStatus.CORE, status_at_head=NodeStatus.CORE
... )
array([10, 11, 14, 15, 19])
>>> grid.link_with_node_status(
...     status_at_tail=NodeStatus.CORE, status_at_head=NodeStatus.FIXED_VALUE
... )
array([12, 16, 20, 23, 24])
>>> grid.link_with_node_status(
...     status_at_tail=NodeStatus.FIXED_VALUE, status_at_head=NodeStatus.CORE
... )
array([ 5,  7,  9, 18])
>>> grid.link_with_node_status(status_at_head=NodeStatus.CORE)
array([ 5,  6,  7,  9, 10, 11, 14, 15, 18, 19])
>>> grid.link_with_node_status(status_at_tail=NodeStatus.CORE)
array([10, 11, 12, 14, 15, 16, 19, 20, 23, 24])
>>> grid.link_with_node_status()
array([ 0,  1,  2,  3,  4,  5,  6,  7,  8,  9, 10, 11, 12, 13, 14, 15, 16,
       17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30])

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)#
property looped_neighbors_at_cell#

For each cell in a raster, return the D8 neighboring cells, looping across grid boundaries as necessary.

Returns lists of looped neighbor cell IDs of given cell ids. If cell ids are not given, it returns a 2D array of size (self.number_of_cells, 8). Order or neighbors is [ E, NE, N, NW, W, SW, S, SE ]

Output is looped, regardless of boundary conditions! (see examples)

Returns:

The eight neighbors of each cell.

Return type:

ndarray (num_cells, 8)

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> neighbors = grid.looped_neighbors_at_cell
>>> neighbors[1, :]
array([2, 5, 4, 3, 0, 3, 4, 5])
>>> neighbors[5, :]
array([3, 0, 2, 1, 4, 1, 2, 0])
>>> grid.looped_neighbors_at_cell[np.array([1, 5]), :]
array([[2, 5, 4, 3, 0, 3, 4, 5],
       [3, 0, 2, 1, 4, 1, 2, 0]])
property looped_neighbors_at_patch#

For each patch in a raster, return the D8 neighboring patches, looping

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 (x,y) vector components of data_at_link onto nodes.

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> link_data = np.arange(rmg.number_of_links)
>>> x, y = map_link_vector_components_to_node_raster(rmg, link_data)
>>> x[5:7]
array([ 7.5,  8.5])
>>> y[5:7]
array([ 7.5,  8.5])

Map the vector sum of links around a patch to the patch.

The resulting vector is returned as a length-2 list, with the two items being arrays of the x component and the y component of the resolved vectors at the patches, respectively.

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

  • ignore_inactive_links (bool) – If True, do not incorporate inactive links into calc. If all links are inactive at a patch, record zero if out is None or leave the existing value if out.

  • out (len-2 list of npatches-long arrays, optional) – Buffer to place mapped values into or None to create a new array.

Returns:

[x_component_of_link_vals_at_patch, y_component_of_link_vals_at_patch].

Return type:

len-2 list of arrays

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_link_vector_sum_to_patch
>>> from landlab import HexModelGrid
>>> mg = HexModelGrid((4, 3))
>>> interior_nodes = mg.status_at_node == mg.BC_NODE_IS_CORE
>>> exterior_nodes = mg.status_at_node != mg.BC_NODE_IS_CORE

Add a ring of closed nodes at the edge:

>>> mg.status_at_node[exterior_nodes] = mg.BC_NODE_IS_CLOSED

This gives us 5 core nodes, 7 active links, and 3 present patches

>>> (mg.number_of_core_nodes == 5 and mg.number_of_active_links == 7)
True
>>> A = mg.add_ones("vals", at="link")
>>> A.fill(9.0)  # any old values on the inactive links
>>> A[mg.active_links] = np.array([1.0, -1.0, 1.0, -1.0, -1.0, -1.0, -1.0])

This setup should give present patch 0 pure east, patch 1 zero (vorticity), and patch 2 westwards and downwards components.

>>> xcomp, ycomp = map_link_vector_sum_to_patch(mg, "vals")
>>> xcomp, ycomp = np.round(xcomp, decimals=5), np.round(ycomp, decimals=5)
>>> np.allclose(xcomp[(6, 9, 10),], [2.0, 0.0, -1.0])
True
>>> np.allclose(ycomp[(6, 9, 10),] / np.sqrt(3.0), [0.0, 0.0, -1.0])
True

These are the patches with LinksStatus.INACTIVE on all three sides:

>>> absent_patches = np.array([0, 1, 2, 4, 8, 11, 12, 15, 16, 17, 18])
>>> np.allclose(xcomp[absent_patches], 0.0)
True
>>> np.allclose(ycomp[absent_patches], 0.0)
True

Now demonstrate the remaining functionality:

>>> A = mg.at_link["vals"].copy()
>>> A.fill(1.0)
>>> _ = map_link_vector_sum_to_patch(
...     mg, A, ignore_inactive_links=False, out=[xcomp, ycomp]
... )
>>> np.allclose(xcomp[absent_patches], 0.0)
False
>>> np.allclose(ycomp[absent_patches], 0.0)
False

Map the maximum of links entering a node to the node.

map_max_of_inlinks_to_node takes an array at the links and finds the inlink values for each node in the grid. it finds the maximum value at the the inlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_max_of_inlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_max_of_inlinks_to_node(rmg, "z")
array([  0.,   0.,   1.,   2.,
         3.,   7.,   8.,   9.,
        10.,  14.,  15.,  16.])

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 max of links leaving a node to the node.

map_max_of_outlinks_to_node takes an array at the links and finds the outlink values for each node in the grid. it finds the maximum value at the the outlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_max_of_outlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_max_of_outlinks_to_node(rmg, "z")
array([  3.,   4.,   5.,   6.,  10.,  11.,  12.,  13.,  14.,  15.,  16.,
         0.])
map_max_of_patch_nodes_to_patch(var_name, ignore_closed_nodes=True, out=None)#

Map the maximum value of nodes around a patch to the patch.

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

  • ignore_closed_nodes (bool) – If True, do not incorporate closed nodes into calc. If all nodes are masked at a patch, record zero if out is None or leave the existing value if out.

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

Returns:

Mapped values at patches.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_max_of_patch_nodes_to_patch
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_node["vals"] = [
...     [5.0, 4.0, 3.0, 2.0],
...     [3.0, 4.0, 3.0, 2.0],
...     [3.0, 2.0, 1.0, 0.0],
... ]
>>> map_max_of_patch_nodes_to_patch(rmg, "vals")
array([ 5., 4., 3.,
        4., 4., 3.])
>>> rmg.at_node["vals"] = [
...     [5.0, 4.0, 3.0, 2.0],
...     [3.0, 4.0, 3.0, 2.0],
...     [3.0, 2.0, 1.0, 0.0],
... ]
>>> rmg.status_at_node[rmg.node_x > 1.5] = rmg.BC_NODE_IS_CLOSED
>>> ans = np.zeros(6, dtype=float)
>>> _ = map_max_of_patch_nodes_to_patch(rmg, "vals", out=ans)
>>> ans
array([ 5., 4., 0.,
        4., 4., 0.])

Map the mean of active links in the x direction touching node to the node.

map_mean_of_horizontal_active_links_to_node takes an array at the links and finds the average of all horizontal (x-direction) link neighbor values for each node in the grid. It returns an array at the nodes of the mean of these values. If a link is absent, it is ignored. If a node has no active links, it receives 0. Note that here a positive returned value means flux to the east, and a negative to the west.

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.raster_mappers import (
...     map_mean_of_horizontal_active_links_to_node,
... )
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", -np.arange(17, dtype=float), at="link")
>>> rmg.status_at_node[rmg.nodes_at_left_edge] = rmg.BC_NODE_IS_CLOSED
>>> map_mean_of_horizontal_active_links_to_node(rmg, "z")
array([ 0. ,  0. ,  0. ,  0. ,  0. , -8. , -8.5, -9. ,  0. ,  0. ,  0. ,
        0. ])

Map the mean of links in the x direction touching a node to the node.

map_mean_of_horizontal_links_to_node takes an array at the links and finds the average of all horizontal (x-direction) link neighbor values for each node in the grid. It returns an array at the nodes of the mean of these values. If a link is absent, it is ignored. Note that here a positive returned value means flux to the east, and a negative to the west.

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.raster_mappers import map_mean_of_horizontal_links_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_mean_of_horizontal_links_to_node(rmg, "z")
array([  0. ,   0.5,   1.5,   2. ,   7. ,   7.5,   8.5,   9. ,  14. ,
        14.5,  15.5,  16. ])

Map the mean of links entering a node to the node.

map_mean_of_inlinks_to_node takes an array at the links and finds the inlink values for each node in the grid. It finds the average of the inlinks and returns values at the nodes.

This considers all inactive links to have a value of 0.

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.raster_mappers import map_mean_of_inlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_mean_of_inlinks_to_node(rmg, "z")
array([  0. ,   0. ,   0.5,   1. ,   1.5,   5.5,   6.5,   7.5,   5. ,
        12.5,  13.5,  14.5])

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 mean of links touching a node to the node.

map_mean_all_links_to_node takes an array at the links and finds the average of all ~existing~ link neighbor values for each node in the grid. it returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_mean_of_links_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_mean_of_links_to_node(rmg, "z")
array([  1.5       ,   1.66666667,   2.66666667,   4.        ,
         6.66666667,   7.5       ,   8.5       ,   9.33333333,
        12.        ,  13.33333333,  14.33333333,  14.5       ])

Map the mean of links leaving a node to the node.

map_mean_of_outlinks_to_node takes an array at the links and finds the outlink values for each node in the grid. it finds the average of the outlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_mean_of_outlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_mean_of_outlinks_to_node(rmg, "z")
array([  1.5,   2.5,   3.5,   3. ,   8.5,   9.5,  10.5,   6.5,   7. ,
         7.5,   8. ,   0. ])
map_mean_of_patch_nodes_to_patch(var_name, ignore_closed_nodes=True, out=None)#

Map the mean value of nodes around a patch to the patch.

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

  • ignore_closed_nodes (bool) – If True, do not incorporate closed nodes into calc. If all nodes are masked at a patch, record zero if out is None or leave the existing value if out.

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

Returns:

Mapped values at patches.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_mean_of_patch_nodes_to_patch
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_node["vals"] = [
...     [5.0, 4.0, 3.0, 2.0],
...     [5.0, 4.0, 3.0, 2.0],
...     [3.0, 2.0, 1.0, 0.0],
... ]
>>> map_mean_of_patch_nodes_to_patch(rmg, "vals")
array([ 4.5, 3.5, 2.5,
        3.5, 2.5, 1.5])
>>> rmg.at_node["vals"] = [
...     [5.0, 4.0, 3.0, 2.0],
...     [5.0, 4.0, 3.0, 2.0],
...     [3.0, 2.0, 1.0, 0.0],
... ]
>>> rmg.status_at_node[rmg.node_x > 1.5] = rmg.BC_NODE_IS_CLOSED
>>> ans = np.zeros(6, dtype=float)
>>> _ = map_mean_of_patch_nodes_to_patch(rmg, "vals", out=ans)
>>> ans
array([ 4.5, 4. , 0. ,
        3.5, 3. , 0. ])

Map the mean of active links in the y direction touching node to the node.

map_mean_of_vertical_active_links_to_node takes an array at the links and finds the average of all vertical (y-direction) link neighbor values for each node in the grid. It returns an array at the nodes of the mean of these values. If a link is absent, it is ignored. If a node has no active links, it receives 0. Note that here a positive returned value means flux to the north, and a negative to the south.

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.raster_mappers import (
...     map_mean_of_vertical_active_links_to_node,
... )
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", -np.arange(17, dtype=float), at="link")
>>> rmg.status_at_node[rmg.nodes_at_bottom_edge] = rmg.BC_NODE_IS_CLOSED
>>> map_mean_of_vertical_active_links_to_node(rmg, "z")
array([  0.,   0.,   0.,   0.,   0., -11., -12.,   0.,   0., -11., -12.,
         0.])

Map the mean of links in the y direction touching a node to the node.

map_mean_of_vertical_links_to_node takes an array at the links and finds the average of all vertical (y-direction) link neighbor values for each node in the grid. It returns an array at the nodes of the mean of these values. If a link is absent, it is ignored. Note that here a positive returned value means flux to the north, and a negative to the south.

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.raster_mappers import map_mean_of_vertical_links_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_mean_of_vertical_links_to_node(rmg, "z")
array([  3. ,   4. ,   5. ,   6. ,   6.5,   7.5,   8.5,   9.5,  10. ,
        11. ,  12. ,  13. ])

Map the minimum of links entering a node to the node.

map_min_of_inlinks_to_node takes an array at the links and finds the inlink values for each node in the grid. it finds the minimum value at the the inlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_min_of_inlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_min_of_inlinks_to_node(rmg, "z")
array([  0.,   0.,   0.,   0.,   0.,   4.,   5.,   6.,   0.,  11.,  12.,
        13.])

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

Map the min of links leaving a node to the node.

map_min_of_outlinks_to_node takes an array at the links and finds the outlink values for each node in the grid. It finds the minimum value at the the outlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_min_of_outlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_min_of_outlinks_to_node(rmg, "z")
array([ 0.,  1.,  2.,  0.,  7.,  8.,  9.,  0.,  0.,  0.,  0.,  0.])
map_min_of_patch_nodes_to_patch(var_name, ignore_closed_nodes=True, out=None)#

Map the minimum value of nodes around a patch to the patch.

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

  • ignore_closed_nodes (bool) – If True, do not incorporate closed nodes into calc. If all nodes are masked at a patch, record zero if out is None or leave the existing value if out.

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

Returns:

Mapped values at patches.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_min_of_patch_nodes_to_patch
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> rmg.at_node["vals"] = [
...     [5.0, 4.0, 3.0, 2.0],
...     [5.0, 4.0, 3.0, 2.0],
...     [3.0, 2.0, 1.0, 0.0],
... ]
>>> map_min_of_patch_nodes_to_patch(rmg, "vals")
array([ 4., 3., 2.,
        2., 1., 0.])
>>> rmg.at_node["vals"] = [
...     [5.0, 4.0, 3.0, 2.0],
...     [5.0, 4.0, 3.0, 2.0],
...     [3.0, 2.0, 1.0, 0.0],
... ]
>>> rmg.status_at_node[rmg.node_x > 1.5] = rmg.BC_NODE_IS_CLOSED
>>> ans = np.zeros(6, dtype=float)
>>> _ = map_min_of_patch_nodes_to_patch(rmg, "vals", out=ans)
>>> ans
array([ 4., 4., 0.,
        2., 2., 0.])
map_node_to_cell(var_name, out=None)#

Map values for nodes to cells.

map_node_to_cell iterates across the grid and identifies the all node values of ‘var_name’.

This function takes node values of ‘var_name’ and mapes that value to the corresponding cell area for each node.

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 cells.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_node_to_cell
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(12.0), at="node")
>>> map_node_to_cell(rmg, "z")
array([ 5.,  6.])
>>> values_at_cells = rmg.empty(at="cell")
>>> rtn = map_node_to_cell(rmg, "z", out=values_at_cells)
>>> values_at_cells
array([ 5.,  6.])
>>> rtn is values_at_cells
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 sum of links entering a node to the node.

map_sum_of_inlinks_to_node takes an array at the links and finds the inlink values for each node in the grid. it sums the inlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_sum_of_inlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_sum_of_inlinks_to_node(rmg, "z")
array([  0.,   0.,   1.,   2.,   3.,  11.,  13.,  15.,  10.,  25.,  27.,
        29.])

Map the sum of links leaving a node to the node.

map_sum_of_outlinks_to_node takes an array at the links and finds the outlink values for each node in the grid. it sums the outlinks and returns values at the nodes.

Note

This considers all inactive links to have a value of 0.

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.raster_mappers import map_sum_of_outlinks_to_node
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(17.0), at="link")
>>> map_sum_of_outlinks_to_node(rmg, "z")
array([  3.,  5.,  7.,   6.,  17.,  19.,  21.,  13.,  14.,  15.,  16.,
         0.])

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])
property material_layers#

MaterialLayers for each cell.

merge(dual, node_at_cell=None, nodes_at_face=None)#
property midpoint_of_face#

Get the middle of faces.

See also

midpoint_of_link

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#

Number of spatial dimensions of the grid.

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
property node_at_cell#
property node_at_core_cell#

Get array of nodes associated with core cells.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))

Initially each cell’s node is core.

>>> grid.node_at_core_cell
array([ 6,  7,  8,
       11, 12, 13])

Setting a node to closed causes means its cell is also “closed”.

>>> grid.status_at_node[8] = grid.BC_NODE_IS_CLOSED
>>> grid.node_at_core_cell
array([ 6,  7, 11, 12, 13])

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])
node_axis_coordinates(axis=0)#

Get the coordinates of nodes along a particular axis.

Return node coordinates from a given axis (defaulting to 0). Axis numbering is the same as that for numpy arrays. That is, the zeroth axis is along the rows, and the first along the columns.

Parameters:

axis (int, optional) – Coordinate axis.

Returns:

Coordinates of nodes for a given axis.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.node_axis_coordinates(0).reshape(grid.shape)
array([[ 0.,  0.,  0.,  0.,  0.],
       [ 1.,  1.,  1.,  1.,  1.],
       [ 2.,  2.,  2.,  2.,  2.],
       [ 3.,  3.,  3.,  3.,  3.]])
>>> grid.node_axis_coordinates(1).reshape(grid.shape)
array([[ 0.,  1.,  2.,  3.,  4.],
       [ 0.,  1.,  2.,  3.,  4.],
       [ 0.,  1.,  2.,  3.,  4.],
       [ 0.,  1.,  2.,  3.,  4.]])
property node_corners#
node_has_boundary_neighbor(ids, method='d8')[source]#

Check if nodes have neighbors that are boundary nodes.

Checks to see if one of the eight neighbor nodes of node(s) with id has a boundary node. Returns True if a node has a boundary node, False if all neighbors are interior.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((5, 5))
>>> mg.node_has_boundary_neighbor(6)
True
>>> mg.node_has_boundary_neighbor(12)
False
>>> mg.node_has_boundary_neighbor([12, -1])
array([False,  True], dtype=bool)
>>> mg.node_has_boundary_neighbor(25)
Traceback (most recent call last):
    ...
IndexError: index 25 is out of bounds for axis 0 with size 25
node_is_boundary(ids, boundary_flag=None)#

Check if nodes are boundary nodes.

Check if nodes at given ids are boundary nodes. Use the boundary_flag to specify a particular boundary type status flag.

Parameters:
  • ids (ndarray) – Node IDs to check.

  • boundary_flag (int, optional) – A boundary type to check for.

Returns:

Array of booleans indicating if nodes are boundary nodes.

Return type:

ndarray

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.node_is_boundary([0, 6])
array([ True, False], dtype=bool)
>>> mg.node_is_boundary([0, 6], boundary_flag=mg.BC_NODE_IS_CLOSED)
array([False, False], dtype=bool)
node_vector_to_raster(u, flip_vertically=False)[source]#

Unravel an array of node values.

Converts node vector u to a 2D array and returns it, so that it can be plotted, output, etc.

If the flip_vertically keyword is True, this function returns an array that has the rows in reverse order. This is useful for use in plot commands (such as the image display functions) that puts the first row at the top of the image. In the landlab coordinate system, the first row is thought to be at the bottom. Thus, a flipped matrix will plot in the landlab style with the first row at the bottom.

The returned array is a view of u, not a copy.

See also

RasterModelGrid.nodes

An equivalent property, but without the option to flip the grid.

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5))
>>> u = rmg.zeros(centering="node")
>>> u = u + range(0, len(u))
>>> u
array([ 0.,   1.,   2.,   3.,   4.,   5.,   6.,   7.,   8.,   9.,  10.,
       11.,  12.,  13.,  14.,  15.,  16.,  17.,  18.,  19.])
>>> ur = rmg.node_vector_to_raster(u)
>>> ur
array([[  0.,   1.,   2.,   3.,   4.],
       [  5.,   6.,   7.,   8.,   9.],
       [ 10.,  11.,  12.,  13.,  14.],
       [ 15.,  16.,  17.,  18.,  19.]])
>>> ur = rmg.node_vector_to_raster(u, flip_vertically=True)
>>> ur
array([[ 15.,  16.,  17.,  18.,  19.],
       [ 10.,  11.,  12.,  13.,  14.],
       [  5.,   6.,   7.,   8.,   9.],
       [  0.,   1.,   2.,   3.,   4.]])
property node_x#
property node_y#
property nodes#

A shaped array of node ids.

Returns:

Node IDs in an array shaped as number_of_node_rows by number_of_node_columns.

Return type:

ndarray

nodes_around_point(xcoord, ycoord)[source]#

Get the nodes surrounding a point.

Return IDs of the four nodes of the area around a point with coordinates xcoord, ycoord. Node IDs are returned counter-clockwise order starting from the southwest node.

If either xcoord or ycoord are arrays the usual numpy broadcasting rules apply.

Parameters:
  • xcoord (float, array-like) – x-coordinate of point

  • ycoord (float, array-like) – y-coordinate of point

Returns:

IDs of nodes around the point.

Return type:

(4, N) ndarray

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.nodes_around_point(0.4, 1.2)
array([4, 8, 9, 5])
>>> grid.nodes_around_point([0.9, 1.1], 1.2)
array([[ 4,  5],
       [ 8,  9],
       [ 9, 10],
       [ 5,  6]])
>>> grid = RasterModelGrid((3, 4), xy_spacing=(1, 2))
>>> grid.nodes_around_point(0.5, 1.5)
array([0, 4, 5, 1])
>>> grid = RasterModelGrid((3, 4))
>>> grid.nodes_around_point(0.5, 1.5)
array([4, 8, 9, 5])
property nodes_at_bottom_edge#
property nodes_at_corners_of_grid#

Nodes at corners of grid.

The nodes at at the corners of the grid. The nodes are returned counterclockwise starting with the upper-right.

Returns:

Nodes at the four corners.

Return type:

tuple of int

Examples

>>> from landlab.graph import UniformRectilinearGraph
>>> graph = UniformRectilinearGraph((4, 5))
>>> graph.nodes_at_corners_of_grid
(19, 15, 0, 4)
property nodes_at_d8#
property nodes_at_diagonal#

Nodes at diagonal tail and head.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.nodes_at_diagonal.shape == (grid.number_of_diagonals, 2)
True
>>> grid.nodes_at_diagonal
array([[ 0,  4], [ 1,  3], [ 1,  5], [ 2,  4],
       [ 3,  7], [ 4,  6], [ 4,  8], [ 5,  7],
       [ 6, 10], [ 7,  9], [ 7, 11], [ 8, 10]])

The first column is diagonal tail and the second the head.

>>> grid.diagonals_at_node[3]
array([ 4, -1, -1,  1])
>>> grid.nodes_at_diagonal[(4, 1),]
array([[3, 7],
       [1, 3]])
>>> grid.diagonal_dirs_at_node[3]
array([-1,  0,  0,  1], dtype=int8)
nodes_at_edge(edge)#
property nodes_at_face#
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#
property number_of_active_faces#

Total number of active faces.

Returns:

Total number of active faces in the grid.

Return type:

int

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.number_of_active_faces
7

The number of active faces is updated when a node status changes.

>>> grid.status_at_node[6] = grid.BC_NODE_IS_CLOSED
>>> grid.number_of_active_faces
3

Number of active links.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.number_of_active_links
17
>>> for edge in (
...     mg.nodes_at_left_edge,
...     mg.nodes_at_right_edge,
...     mg.nodes_at_bottom_edge,
... ):
...     mg.status_at_node[edge] = mg.BC_NODE_IS_CLOSED
>>> mg.number_of_active_links
10
property number_of_cell_columns#

Number of cell columns.

Returns the number of columns, including boundaries.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.number_of_cell_columns
3
property number_of_cell_rows#

Number of cell rows.

Returns the number of rows, including boundaries.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.number_of_cell_rows
2
property number_of_cells#

Get the number of cells.

property number_of_cells_present_at_corner#

Return the number of cells at a corner without a closed corner.

property number_of_cells_present_at_face#

Return the number of cells at a face without a closed corner.

property number_of_core_cells#

Number of core cells.

A core cell excludes all boundary cells.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.number_of_core_cells
6
>>> grid.status_at_node[7] = grid.BC_NODE_IS_CLOSED
>>> grid.number_of_core_cells
5
property number_of_core_corners#

Number of core corners.

property number_of_core_nodes#

Number of core nodes.

The number of core nodes on the grid (i.e., excluding all boundary nodes).

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.number_of_core_nodes
6
>>> grid.status_at_node[7] = grid.BC_NODE_IS_CLOSED
>>> grid.number_of_core_nodes
5
property number_of_core_patches#

Number of core patches.

property number_of_corner_columns#
property number_of_corner_rows#
property number_of_corners#

Get total number of corners.

See also

number_of_nodes

property number_of_d8#

Number of links and diagonals.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.number_of_d8
29
property number_of_diagonals#

Number of diagonals in the grid.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.number_of_diagonals
12
number_of_elements(name)#

Number of instances of an element.

Get the number of instances of a grid element in a grid.

Parameters:

name ({'node', 'cell', 'link', 'face', 'core_node', 'core_cell',) – ‘active_link’, ‘active_face’} Name of the grid element.

Returns:

Number of elements in the grid.

Return type:

int

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.number_of_elements("node")
20
>>> mg.number_of_elements("core_cell")
6
>>> mg.number_of_elements("link")
31
>>> mg.number_of_elements("active_link")
17
>>> mg.status_at_node[8] = mg.BC_NODE_IS_CLOSED
>>> mg.number_of_elements("link")
31
>>> mg.number_of_elements("active_link")
13
property number_of_faces#

Get corners at face head.

See also

number_of_links

property number_of_fixed_faces#

Number of fixed faces.

Number of fixed links.

Examples

>>> from landlab import NodeStatus, RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.number_of_fixed_links
0
>>> mg.status_at_node[mg.nodes_at_top_edge] = NodeStatus.FIXED_GRADIENT
>>> mg.number_of_fixed_links
3
property number_of_interior_corners#

Number of interior corners.

property number_of_interior_nodes#

Number of interior nodes.

Returns the number of interior nodes on the grid, i.e., non-perimeter nodes. Compare self.number_of_core_nodes.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 5))
>>> grid.number_of_interior_nodes
6

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_node_columns#
property number_of_node_rows#
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_patch_columns#

Number of patch columns.

property number_of_patch_rows#

Number of patch rows.

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

Return the number of patches at a link without a closed node.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.patches_present_at_link
array([[False,  True],
       [False,  True],
       [ True, False],
       [ True,  True],
       [False,  True],
       [ True, False],
       [ True, False],
       [False, False],
       [False, False],
       [False, False],
       [False, False],
       [False, False]], dtype=bool)
>>> mg.number_of_patches_present_at_link
array([1, 1, 1, 2, 1, 1, 1, 0, 0, 0, 0, 0])
property number_of_patches_present_at_node#

Return the number of patches at a node without a closed node.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.patches_present_at_node
array([[ True, False, False, False],
       [ True,  True, False, False],
       [False,  True, False, False],
       [False, False, False,  True],
       [False, False,  True,  True],
       [False, False,  True, False],
       [False, False, False, False],
       [False, False, False, False],
       [False, False, False, False]], dtype=bool)
>>> mg.number_of_patches_present_at_node
array([1, 2, 1, 1, 2, 1, 0, 0, 0])
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 open_boundary_corners#

Get array of open boundary corners.

property open_boundary_nodes#

Get array of open boundary nodes.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> for edge in (
...     mg.nodes_at_left_edge,
...     mg.nodes_at_right_edge,
...     mg.nodes_at_bottom_edge,
... ):
...     mg.status_at_node[edge] = mg.BC_NODE_IS_CLOSED
>>> mg.open_boundary_nodes
array([16, 17, 18])
property orientation_of_face#

Return array of face orientation codes (one value per face).

Return array of link orientation codes (one value per link).

Orientation codes are defined by LinkOrientation; 1 = E, 2 = ENE, 4 = NNE, 8 = N, 16 = NNW, 32 = ESE (using powers of 2 allows for future applications that might want additive combinations).

property origin#
property parallel_faces_at_face#

Return similarly oriented faces connected to each face.

Return similarly oriented links connected to each link.

Return IDs of links of the same orientation that are connected to each given link’s tail or head node.

The data structure is a numpy array of shape (n_links, 2) containing the IDs of the “tail-wise” (connected to tail node) and “head-wise” (connected to head node) links, or -1 if the link is inactive (e.g., on the perimeter) or it has no attached parallel neighbor in the given direction.

For instance, consider a 3x4 raster, in which link IDs are as shown:

.-14-.-15-.-16-.
|    |    |    |
10  11   12   13
|    |    |    |
.--7-.--8-.--9-.
|    |    |    |
3    4    5    6
|    |    |    |
.--0-.--1-.--2-.

Here’s a mapping of the tail-wise (shown at left or bottom of links) and head-wise (shown at right or top of links) links:

.----.----.----.
|    |    |    |
|    |    |    |
|    4    5    |
.---8.7--9.8---.
|   11   12    |
|    |    |    |
|    |    |    |
.----.----.----.

So the corresponding data structure would be mostly filled with -1, but for the 7 active links, it would look like:

4: [[-1, 11],
5:  [-1, 12],
7:  [-1,  8],
8:  [ 7,  9],
9:  [ 8, -1],
11: [ 4, -1],
12: [ 5, -1]]

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> pll = grid.parallel_links_at_link
>>> pll[4:13, :]
array([[-1, 11],
       [-1, 12],
       [-1, 13],
       [-1,  8],
       [ 7,  9],
       [ 8, -1],
       [ 3, -1],
       [ 4, -1],
       [ 5, -1]])
property patch_area_at_corner#

Cell areas in a ncorners-long array.

property patch_at_corner#
property patch_grid_shape#

Get the shape of the patchular grid (grid with only patches).

See also

cell_grid_shape

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 patches_at_nodes_of_grid#

Get array of patches in patchular grid (grid with only patches) nodes.

A boolean array, False where a patch has a closed node or is missing.

The array is the same shape as patches_at_link, and is designed to mask it.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.patches_at_link
array([[-1,  0],
       [-1,  1],
       [ 0, -1],
       [ 1,  0],
       [-1,  1],
       [ 0,  2],
       [ 1,  3],
       [ 2, -1],
       [ 3,  2],
       [-1,  3],
       [ 2, -1],
       [ 3, -1]])
>>> mg.patches_present_at_link
array([[False,  True],
       [False,  True],
       [ True, False],
       [ True,  True],
       [False,  True],
       [ True, False],
       [ True, False],
       [False, False],
       [False, False],
       [False, False],
       [False, False],
       [False, False]], dtype=bool)
>>> 1 in mg.patches_at_link * mg.patches_present_at_link
True
>>> 2 in mg.patches_at_link * mg.patches_present_at_link
False
property patches_present_at_node#

A boolean array, False where a patch has a closed node or is missing.

The array is the same shape as patches_at_node, and is designed to mask it.

Note that in cases where patches may have more than 3 nodes (e.g., rasters), a patch is considered still present as long as at least 3 open nodes are present.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((3, 3))
>>> mg.status_at_node[mg.nodes_at_top_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.patches_at_node
array([[ 0, -1, -1, -1],
       [ 1,  0, -1, -1],
       [-1,  1, -1, -1],
       [ 2, -1, -1,  0],
       [ 3,  2,  0,  1],
       [-1,  3,  1, -1],
       [-1, -1, -1,  2],
       [-1, -1,  2,  3],
       [-1, -1,  3, -1]])
>>> mg.patches_present_at_node
array([[ True, False, False, False],
       [ True,  True, False, False],
       [False,  True, False, False],
       [False, False, False,  True],
       [False, False,  True,  True],
       [False, False,  True, False],
       [False, False, False, False],
       [False, False, False, False],
       [False, False, False, False]], dtype=bool)
>>> 1 in mg.patches_at_node * mg.patches_present_at_node
True
>>> 2 in mg.patches_at_node * mg.patches_present_at_node
False
property perimeter_corners#
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()#

Resolve the xy-components of links.

Resolves values provided defined on links into the x and y directions. Returns values_along_x, values_along_y

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.
roll_nodes_ud(data_name, shift, interior_only=False)[source]#

Roll (shift) specified data on nodes up or down in a raster grid.

Similar to the Numpy roll() function, in that it shifts node values up by shift rows, wrapping the data in the top row(s) around to the bottom. If the interior_only is set, data along the left and right grid edges are not changed.

Note that the contents of the data_name field are changed.

Parameters:
  • data_name (string) – Name of node-data item attached to grid.

  • shift (int) – Number of rows to shift upward.

  • interior_only (bool, optional) – If True, data along left and right edges are not shifted

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 3))
>>> data = rmg.add_zeros("test_data", at="node")
>>> data[:] = np.arange(12)
>>> rmg.roll_nodes_ud("test_data", 1)
>>> data.reshape(rmg.shape)
array([[  9.,  10.,  11.],
       [  0.,   1.,   2.],
       [  3.,   4.,   5.],
       [  6.,   7.,   8.]])
>>> rmg.roll_nodes_ud("test_data", 2)
>>> data.reshape(rmg.shape)
array([[  3.,   4.,   5.],
       [  6.,   7.,   8.],
       [  9.,  10.,  11.],
       [  0.,   1.,   2.]])
>>> rmg.roll_nodes_ud("test_data", 1, interior_only=True)
>>> data.reshape(rmg.shape)
array([[  3.,   1.,   5.],
       [  6.,   4.,   8.],
       [  9.,   7.,  11.],
       [  0.,  10.,   2.]])
save(path, names=None, format=None, at=None)[source]#

Save a grid and fields.

If more than one field name is specified for names when saving to ARC ascii, multiple files will be produced, suffixed with the field names.

When saving to netCDF (.nc), the fields are incorporated into the single named .nc file.

Parameters:
  • path (str) – Path to output file.

  • names (iterable of strings, optional) – List of field names to save, defaults to all if not specified.

  • format ({'netcdf', 'esri-ascii'}, optional) – Output file format. Guess from file extension if not given.

  • at (str) – Grid element where values are defined.

Examples

>>> from landlab import RasterModelGrid
>>> import os
>>> from tempfile import mkdtemp
>>> grid = RasterModelGrid((4, 5))
>>> fname = os.path.join(mkdtemp(), "mysave.nc")
>>> grid.save(fname)
>>> os.path.isfile(fname)
True
>>> os.remove(fname)
property second_ring_looped_neighbors_at_cell#

Get list of second ring looped neighbor cell IDs (all 16 neighbors).

Returns lists of looped second ring neighbor cell IDs of given cell ids. If cell ids are not given, it returns a 2D array of size ( self.number_of_cells, 16 ).

The cells are the 16 which encircle the nine true neighbor cells. Order of neighbors: Starts with E and goes counter clockwise

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((10, 10))
>>> mg.second_ring_looped_neighbors_at_cell[36, :]
array([38, 46, 54, 53, 52, 51, 50, 42, 34, 26, 18, 19, 20, 21, 22, 30])
>>> mg.second_ring_looped_neighbors_at_cell[8, :]
array([10, 18, 26, 25, 24, 31, 30, 22, 14,  6, 62, 63, 56, 57, 58,  2])

…take a look at the cell grid to understand why:

[56, 57, 58, 59, 60, 61, 62, 63]
[48, 49, 50, 51, 52, 53, 54, 55]
[40, 41, 42, 43, 44, 45, 46, 47]
[32, 33, 34, 35, 36, 37, 38, 39]
[24, 25, 26, 27, 28, 29, 30, 31]
[16, 17, 18, 19, 20, 21, 22, 23]
[ 8,  9, 10, 11, 12, 13, 14, 15]
[ 0,  1,  2,  3,  4,  5,  6,  7]
property second_ring_looped_neighbors_at_patch#

Get list of second ring looped neighbor patch IDs (all 16 neighbors).

set_closed_boundaries_at_grid_edges(right_is_closed, top_is_closed, left_is_closed, bottom_is_closed)[source]#

Set boundary not to be closed.

Sets the status of nodes along the specified side(s) of a raster grid (bottom, right, top, and/or left) to BC_NODE_IS_CLOSED.

Arguments are booleans indicating whether the bottom, left, top, and right are closed (True) or not (False).

For a closed boundary:

  • the nodes are flagged BC_NODE_IS_CLOSED (status type 4)

  • all links that connect to a BC_NODE_IS_CLOSED node are flagged as inactive (so they appear on link-based lists, but not active_link-based lists)

This means that if you call the calc_grad_at_active_link method, links connecting to closed boundaries will be ignored: there can be no gradients or fluxes calculated, because the links that connect to that edge of the grid are not included in the calculation. So, setting a grid edge to BC_NODE_IS_CLOSED is a convenient way to impose a no-flux boundary condition. Note, however, that this applies to the grid as a whole, rather than a particular variable that you might use in your application. In other words, if you want a no-flux boundary in one variable but a different boundary condition for another, then use another method.

Parameters:
  • right_is_closed (boolean) – If True right-edge nodes are closed boundaries.

  • top_is_closed (boolean) – If True top-edge nodes are closed boundaries.

  • left_is_closed (boolean) – If True left-edge nodes are closed boundaries.

  • bottom_is_closed (boolean) – If True bottom-edge nodes are closed boundaries.

Notes

Note that the four corners are treated as follows:

  • bottom left = BOTTOM

  • bottom right = BOTTOM

  • top right = TOP

  • top left = TOP

This scheme is necessary for internal consistency with looped boundaries.

Examples

The following example sets the top and left boundaries as closed in a four-row by five-column grid that initially has all boundaries open and all boundary nodes coded as BC_NODE_IS_FIXED_VALUE (=1):

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5))  # rows, columns, spacing
>>> rmg.number_of_active_links
17
>>> rmg.status_at_node.reshape(rmg.shape)
array([[1, 1, 1, 1, 1],
       [1, 0, 0, 0, 1],
       [1, 0, 0, 0, 1],
       [1, 1, 1, 1, 1]], dtype=uint8)
>>> rmg.set_closed_boundaries_at_grid_edges(True, True, False, False)
>>> rmg.number_of_active_links
12
>>> rmg.status_at_node.reshape(rmg.shape)
array([[1, 1, 1, 1, 1],
       [1, 0, 0, 0, 4],
       [1, 0, 0, 0, 4],
       [4, 4, 4, 4, 4]], dtype=uint8)
set_fixed_value_boundaries_at_grid_edges(right_is_fixed_val, top_is_fixed_val, left_is_fixed_val, bottom_is_fixed_val, value=None, value_of='topographic__elevation')[source]#

Create fixed values boundaries.

Sets the status of nodes along the specified side(s) of a raster grid—bottom, right, top, and/or left—to NODE_IS_FIXED_VALUE

Arguments are booleans indicating whether the bottom, right, top, and left sides are to be set (True) or not (False).

value controls what values are held constant at these nodes. It can be either a float, an array of length number_of_fixed_nodes or number_of_nodes (total), or left blank. If left blank, the values will be set from the those already in the grid fields, according to ‘value_of’.

value_of controls the name of the model field that contains the values. Remember, if you don’t set value, the fixed values will be set from the field values *at the time you call this method*. If no values are present in the field, the module will complain but accept this, warning that it will be unable to automatically update boundary conditions.

The status of links (active or inactive) is automatically updated to reflect the changes.

The following example sets the bottom and right boundaries as fixed-value in a four-row by five-column grid that initially has all boundaries closed (i.e., flagged as node_status=4):

Parameters:
  • bottom_is_fixed_val (boolean) – Set bottom edge as fixed boundary.

  • left_is_fixed_val (boolean) – Set left edge as fixed boundary.

  • top_is_fixed_val (boolean) – Set top edge as fixed boundary.

  • right_is_fixed_val (boolean) – Set right edge as fixed boundary.

  • value (float, array or None (default).) – Override value to be kept constant at nodes.

  • value_of (string.) – The name of the grid field containing the values of interest.

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5), xy_spacing=(1, 1))
>>> rmg.number_of_active_links
17

Put some arbitrary values in the grid fields:

>>> import numpy as np
>>> rmg.at_node["topographic__elevation"] = np.random.rand(20)
>>> rmg.set_closed_boundaries_at_grid_edges(True, True, True, True)
>>> rmg.status_at_node.reshape(rmg.shape)
array([[4, 4, 4, 4, 4],
       [4, 0, 0, 0, 4],
       [4, 0, 0, 0, 4],
       [4, 4, 4, 4, 4]], dtype=uint8)
>>> rmg.set_fixed_value_boundaries_at_grid_edges(True, True, False, False)
>>> rmg.number_of_active_links
12
>>> rmg.status_at_node.reshape(rmg.shape)
array([[4, 4, 4, 4, 4],
       [4, 0, 0, 0, 1],
       [4, 0, 0, 0, 1],
       [1, 1, 1, 1, 1]], dtype=uint8)

Note that the four corners are treated as follows:

  • bottom left = BOTTOM

  • bottom right = BOTTOM

  • top right = TOP

  • top left = TOP

This scheme is necessary for internal consistency with looped boundaries.

set_looped_boundaries(top_bottom_are_looped, sides_are_looped)[source]#

Create wrap-around boundaries.

Handles boundary conditions by setting corresponding parallel grid edges as looped “tracks_cell” (==3) status, linked to each other. If top_bottom_are_looped is True, the top and bottom edges will link to each other. If sides_are_looped is True, the left and right edges will link to each other.

Looped boundaries are experimental, and not as yet well integrated into the Landlab framework. Many functions may not recognise them, or silently create unforeseen errors. Use at your own risk!

Note that because of the symmetries this BC implies, the corner nodes are all paired with the bottom/top edges, not the sides.

Parameters:
  • top_bottom_are_looped (bool) – Top and bottom are wrap-around.

  • sides_are_looped (bool) – Left and right sides are wrap-around.

Examples

>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 5))  # rows, columns, spacing
>>> rmg.number_of_active_links
17
>>> rmg.status_at_node.reshape(rmg.shape)
array([[1, 1, 1, 1, 1],
       [1, 0, 0, 0, 1],
       [1, 0, 0, 0, 1],
       [1, 1, 1, 1, 1]], dtype=uint8)
>>> rmg.add_zeros("topographic__elevation", at="node")
array([ 0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,  0.,
        0.,  0.,  0.,  0.,  0.,  0.,  0.])
>>> rmg.set_looped_boundaries(True, True)
>>> rmg.looped_node_properties["boundary_node_IDs"]
array([ 0,  1,  2,  3,  4,  5,  9, 10, 14, 15, 16, 17, 18, 19])
>>> rmg.looped_node_properties["linked_node_IDs"]
array([10, 11, 12, 13, 14,  8,  6, 13, 11,  5,  6,  7,  8,  9])
set_nodata_nodes_to_closed(node_data, nodata_value)#

Make no-data nodes closed boundaries.

Sets node status to BC_NODE_IS_CLOSED for all nodes whose value of node_data is equal to the nodata_value.

Any links connected to BC_NODE_IS_CLOSED nodes are automatically set to LinkStatus.INACTIVE boundary.

Parameters:
  • node_data (ndarray) – Data values.

  • nodata_value (float) – Value that indicates an invalid value.

Examples

The following example uses the following grid:

*--I--->o------>o------>o
^       ^       ^       ^
I       I       |       |
|       |       |       |
*--I--->*--I--->o------>o
^       ^       ^       ^
I       I       I       I
|       |       |       |
*--I--->*--I--->*--I--->*

Note

Links set to LinkStatus.ACTIVE are not shown in this diagram.

* indicates the nodes that are set to NodeStatus.CLOSED

o indicates the nodes that are set to NodeStatus.CORE

I indicates the links that are set to LinkStatus.INACTIVE

>>> import numpy as np
>>> import landlab as ll
>>> mg = ll.RasterModelGrid((3, 4))
>>> mg.status_at_node.reshape(mg.shape)
array([[1, 1, 1, 1],
       [1, 0, 0, 1],
       [1, 1, 1, 1]], dtype=uint8)
>>> h = np.array(
...     [
...         [-9999, -9999, -9999, -9999],
...         [-9999, -9999, 12345.0, 0.0],
...         [-9999, 0.0, 0.0, 0.0],
...     ]
... ).flatten()
>>> mg.set_nodata_nodes_to_closed(h, -9999)
>>> mg.status_at_node.reshape(mg.shape)
array([[4, 4, 4, 4],
       [4, 4, 0, 1],
       [4, 1, 1, 1]], dtype=uint8)
set_nodata_nodes_to_fixed_gradient(node_data, nodata_value)#

Make no-data nodes fixed gradient boundaries.

Set node status to BC_NODE_IS_FIXED_VALUE for all nodes whose value of node_data is equal to nodata_value.

Any links between BC_NODE_IS_FIXED_GRADIENT nodes and BC_NODE_IS_CORE are automatically set to LinkStatus.FIXED boundary status.

Parameters:
  • node_data (ndarray) – Data values.

  • nodata_value (float) – Value that indicates an invalid value.

Examples

The following examples use this grid:

*--I--->*--I--->*--I--->*--I--->*--I--->*--I--->*--I--->*--I--->*
^       ^       ^       ^       ^       ^       ^       ^       ^
I       I       I       X       X       X       X       X       I
|       |       |       |       |       |       |       |       |
*--I--->*--I--->*--X--->o       o       o       o       o--X--->*
^       ^       ^       ^       ^       ^       ^       ^       ^
I       I       I       |       |       |       |       |       I
|       |       |       |       |       |       |       |       |
*--I--->*--I--->*--X--->o       o       o       o       o--X--->*
^       ^       ^       ^       ^       ^       ^       ^       ^
I       I       I       X       X       X       X       X       I
|       |       |       |       |       |       |       |       |
*--I--->*--I--->*--I--->*--I--->*--I--->*--I--->*--I--->*--I--->*

Note

Links set to LinkStatus.ACTIVE are not shown in this diagram.

X indicates the links that are set to LinkStatus.FIXED

I indicates the links that are set to LinkStatus.INACTIVE

o indicates the nodes that are set to NodeStatus.CORE

* indicates the nodes that are set to

BC_NODE_IS_FIXED_GRADIENT

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 9))
>>> rmg.status_at_node.reshape(rmg.shape)
array([[1, 1, 1, 1, 1, 1, 1, 1, 1],
       [1, 0, 0, 0, 0, 0, 0, 0, 1],
       [1, 0, 0, 0, 0, 0, 0, 0, 1],
       [1, 1, 1, 1, 1, 1, 1, 1, 1]], dtype=uint8)
>>> z = np.array(
...     [
...         [-99.0, -99.0, -99.0, -99.0, -99.0, -99.0, -99.0, -99.0, -99.0],
...         [-99.0, -99.0, -99.0, 0.0, 0.0, 0.0, 0.0, 0.0, -99.0],
...         [-99.0, -99.0, -99.0, 0.0, 0.0, 0.0, 0.0, 0.0, -99.0],
...         [-99.0, -99.0, -99.0, -99.0, -99.0, -99.0, -99.0, -99.0, -99.0],
...     ]
... ).flatten()
>>> rmg.set_nodata_nodes_to_fixed_gradient(z, -99)
>>> rmg.status_at_node.reshape(rmg.shape)
array([[2, 2, 2, 2, 2, 2, 2, 2, 2],
       [2, 2, 2, 0, 0, 0, 0, 0, 2],
       [2, 2, 2, 0, 0, 0, 0, 0, 2],
       [2, 2, 2, 2, 2, 2, 2, 2, 2]], dtype=uint8)
>>> rmg.status_at_link[rmg.horizontal_links].reshape((4, 8))
array([[4, 4, 4, 4, 4, 4, 4, 4],
       [4, 4, 2, 0, 0, 0, 0, 2],
       [4, 4, 2, 0, 0, 0, 0, 2],
       [4, 4, 4, 4, 4, 4, 4, 4]], dtype=uint8)
>>> rmg.status_at_link[rmg.vertical_links].reshape((3, 9))
array([[4, 4, 4, 2, 2, 2, 2, 2, 4],
       [4, 4, 4, 0, 0, 0, 0, 0, 4],
       [4, 4, 4, 2, 2, 2, 2, 2, 4]], dtype=uint8)
set_open_nodes_disconnected_from_watershed_to_closed(node_data, outlet_id=None, nodata_value=-9999.0, adjacency_method='D8')[source]#

Identifys all non-closed nodes that are disconnected from the node given in.

outlet_id and sets them as closed.

If outlet_id is not given, the outlet will be identified as the node for which the status at the node is BC_NODE_IS_FIXED_VALUE. If more than one node has this value, the algorithm will fail.

If outlet_id is given, the algorithm will check that it is not a node with status of BC_NODE_IS_CLOSED.

The method supports both D4 and D8 (default) neighborhood evaluation in determining if a node is connected. This can be modified with the flag adjacency_method.

This function can be run directly, or by setting the flag remove_disconnected to True in set_watershed_boundary_condition

Parameters:
  • node_data (field name or ndarray) – At-node field name or at-node data values to use for identifying watershed location.

  • outlet_id (one element numpy array, optional.) – The node ID of the outlet that all open nodes must be connected to. If a node ID is provided, it does not need have the status BC_NODE_IS_FIXED_VALUE. However, it must not have the status of BC_NODE_IS_CLOSED.

  • nodata_value (float, optional, default is -9999.) – Value that indicates an invalid value.

  • adjacency_method (string, optional. Default is 'D8'.) – Sets the connection method.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> mg1 = RasterModelGrid((4, 6))
>>> z1 = np.array(
...     [
...         [-9999.0, -9999.0, -9999.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, 67.0, 67.0, -9999.0, 50.0, -9999.0],
...         [-9999.0, 67.0, 0.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, -9999.0, -9999.0, -9999.0, -9999.0, -9999.0],
...     ]
... ).flatten()
>>> mg2 = RasterModelGrid((4, 6))
>>> z2 = np.array(
...     [
...         [-9999.0, -9999.0, -9999.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, 67.0, 67.0, -9999.0, 50.0, -9999.0],
...         [-9999.0, 67.0, 0.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, -9999.0, -9999.0, -9999.0, -9999.0, -9999.0],
...     ]
... ).flatten()
>>> mg1.set_watershed_boundary_condition(z1, remove_disconnected=True)
>>> mg2.set_watershed_boundary_condition(z2)
>>> mg2.status_at_node.reshape(mg2.shape)
array([[4, 4, 4, 4, 4, 4],
       [4, 0, 0, 4, 0, 4],
       [4, 0, 1, 4, 4, 4],
       [4, 4, 4, 4, 4, 4]], dtype=uint8)
>>> mg2.set_open_nodes_disconnected_from_watershed_to_closed(z2)
>>> np.allclose(mg1.status_at_node, mg2.status_at_node)
True
>>> np.allclose(z1, z2)
True
>>> mg2.status_at_node.reshape(mg2.shape)
array([[4, 4, 4, 4, 4, 4],
       [4, 0, 0, 4, 4, 4],
       [4, 0, 1, 4, 4, 4],
       [4, 4, 4, 4, 4, 4]], dtype=uint8)
>>> z1.reshape(mg1.shape)
array([[-9999., -9999., -9999., -9999., -9999., -9999.],
       [-9999.,    67.,    67., -9999., -9999., -9999.],
       [-9999.,    67.,     0., -9999., -9999., -9999.],
       [-9999., -9999., -9999., -9999., -9999., -9999.]])
set_status_at_node_on_edges(right=None, top=None, left=None, bottom=None)#

Set node status on grid edges.

Parameters:
  • right (int, optional) – Node status along right edge.

  • top (int, optional) – Node status along top edge.

  • left (int, optional) – Node status along left edge.

  • bottom (int, optional) – Node status along bottom edge.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.status_at_node.reshape(grid.shape)
array([[1, 1, 1, 1],
       [1, 0, 0, 1],
       [1, 1, 1, 1]], dtype=uint8)
>>> grid.set_status_at_node_on_edges(right=grid.BC_NODE_IS_CLOSED)
>>> grid.status_at_node.reshape(grid.shape)
array([[1, 1, 1, 4],
       [1, 0, 0, 4],
       [1, 1, 1, 4]], dtype=uint8)
>>> grid = RasterModelGrid((3, 4))

The status of a corner is set along with its clockwise edge. That is, if setting the status for the top and right edges, the upper-right corner has the status of the right edge.

>>> grid.set_status_at_node_on_edges(
...     top=grid.BC_NODE_IS_CLOSED, right=grid.BC_NODE_IS_FIXED_GRADIENT
... )
>>> grid.status_at_node.reshape(grid.shape)
array([[1, 1, 1, 2],
       [1, 0, 0, 2],
       [4, 4, 4, 2]], dtype=uint8)

In the above example, if you wanted the corner to have the status of the top edge, you need to make two calls to set_status_at_node_on_edges,

>>> grid = RasterModelGrid((3, 4))
>>> grid.set_status_at_node_on_edges(right=grid.BC_NODE_IS_FIXED_GRADIENT)
>>> grid.set_status_at_node_on_edges(top=grid.BC_NODE_IS_CLOSED)
>>> grid.status_at_node.reshape(grid.shape)
array([[1, 1, 1, 2],
       [1, 0, 0, 2],
       [4, 4, 4, 4]], dtype=uint8)

An example that sets all of the edges shows how corners are set.

>>> grid.set_status_at_node_on_edges(right=1, top=2, left=3, bottom=4)
>>> grid.status_at_node.reshape(grid.shape)
array([[3, 4, 4, 4],
       [3, 0, 0, 1],
       [2, 2, 2, 1]], dtype=uint8)
set_watershed_boundary_condition(node_data, nodata_value=-9999.0, return_outlet_id=False, remove_disconnected=False, adjacency_method='D8')[source]#

Finds the node adjacent to a boundary node with the smallest value. This node is set as the outlet. The outlet node must have a data value. Can return the outlet id as a one element numpy array if return_outlet_id is set to True.

All nodes with nodata_value are set to BC_NODE_IS_CLOSED (grid.status_at_node == 4). All nodes with data values are set to BC_NODE_IS_CORE (grid.status_at_node == 0), with the exception that the outlet node is set to a BC_NODE_IS_FIXED_GRADIENT (grid.status_at_node == 1).

Note that the outer ring (perimeter) of the raster is set to BC_NODE_IS_CLOSED, even if there are nodes that have values. The only exception to this would be if the outlet node is on the perimeter, which is acceptable.

This routine assumes that all of the nodata_values are on the outside of the data values. In other words, there are no islands of nodata_values surrounded by nodes with data.

This also assumes that the grid has a single watershed (that is a single outlet node).

If you are considering running flow routing algorithms on this model grid, it is recommended that you verify the absence of non-closed nodes surrounded only by closed nodes. The presence of these isolated non- closed nodes may result from clipping a raster to a polygon and will prevent the flow router from functioning as intended.

This can be acomplished by setting the parameter: remove_disconnected to True (default is False).

This will run the function: set_open_nodes_disconnected_from_watershed_to_closed which will find any isolated open nodes that have no neigbors in the main watershed and set them to closed. The adjacency method used to assess connectivity can be set to either ‘D8’(default) or ‘D4’ using the flag adjacency_method.

Finally, the developer has seen cases in which DEM data that has been filled results in a different outlet from DEM data which has not been filled. Be aware that if you identify an outlet on a filled DEM, make sure that the filled DEM is what is being used for your modeling. Otherwise, this may find a different outlet. To force the outlet location, use either set_watershed_boundary_condition_outlet_coords or set_watershed_boundary_condition_outlet_id.

Parameters:
  • node_data (field name or ndarray) – At-node field name or at-node data values to use for identifying watershed location.

  • nodata_value (float, optional) – Value that indicates an invalid value.

  • return_outlet_id (boolean, optional) – Indicates whether or not to return the id of the found outlet

  • remove_disconnected (boolean, optional) – Indicates whether to search for and remove disconnected nodes.

  • adjacency_method (string, optional. Default is 'D8'.) – Sets the connection method for use if remove_disconnected==True

Returns:

outlet_loc – Array of size 1 containing id of outlet location

Return type:

array

Examples

The first example will use a 4,4 grid with node data values as illustrated:

-9999. -9999. -9999. -9999.
-9999.    67.     0. -9999.
-9999.    67.    67. -9999.
-9999. -9999. -9999. -9999.

The second example will use a 4,4 grid with node data values as illustrated:

-9999. -9999. -9999. -9999.
-9999.    67.     0. -9999.
-9999.    67.     67.   -2.
-9999. -9999. -9999. -9999.
>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 4))
>>> node_data = np.array(
...     [
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, 67.0, 67.0, -9999.0],
...         [-9999.0, 67.0, 0.0, -9999.0],
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...     ]
... ).flatten()
>>> out_id = rmg.set_watershed_boundary_condition(node_data, -9999.0, True)
>>> out_id
array([10])
>>> rmg.status_at_node
array([4, 4, 4, 4, 4, 0, 0, 4, 4, 0, 1, 4, 4, 4, 4, 4], dtype=uint8)
>>> rmg2 = RasterModelGrid((4, 4))
>>> node_data2 = np.array(
...     [
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, 67.0, 67.0, -2.0],
...         [-9999.0, 67.0, 0.0, -9999.0],
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...     ]
... ).flatten()
>>> rmg2.set_watershed_boundary_condition(node_data2, -9999.0)
>>> rmg2.status_at_node
array([4, 4, 4, 4, 4, 0, 0, 1, 4, 0, 0, 4, 4, 4, 4, 4], dtype=uint8)

The node data can also be provided as a model grid field.

>>> rmg = RasterModelGrid((4, 4))
>>> node_data = [
...     [-9999.0, -9999.0, -9999.0, -9999.0],
...     [-9999.0, 67.0, 67.0, -9999.0],
...     [-9999.0, 67.0, 0.0, -9999.0],
...     [-9999.0, -9999.0, -9999.0, -9999.0],
... ]
>>> _ = rmg.add_field("topographic__elevation", node_data, at="node")
>>> out_id = rmg.set_watershed_boundary_condition(
...     "topographic__elevation", -9999.0, True
... )
>>> out_id
array([10])
>>> rmg.status_at_node
array([4, 4, 4, 4, 4, 0, 0, 4, 4, 0, 1, 4, 4, 4, 4, 4], dtype=uint8)
set_watershed_boundary_condition_outlet_coords(outlet_coords, node_data, nodata_value=-9999.0)[source]#

Set the boundary conditions for a watershed. All nodes with nodata_value are set to BC_NODE_IS_CLOSED (grid.status_at_node == 4). All nodes with data values are set to CORE_NODES (grid.status_at_node == 0), with the exception that the outlet node is set to a BC_NODE_IS_FIXED_VALUE (grid.status_at_node == 1).

Note that the outer ring of the raster is set to BC_NODE_IS_CLOSED, even if there are nodes that have values. The only exception to this would be if the outlet node is on the boundary, which is acceptable.

Assumes that outlet is already known.

This assumes that the grid has a single watershed. If this is not the case this will not work.

This must be passed the values of the outlet_row and outlet_column. Also takes node_data and optionally, nodata_value.

Parameters:
  • outlet_coords (list - two integer values) – row, column of outlet, NOT THE ABSOLUTE X AND Y LOCATIONS

  • node_data (field name or ndarray) – At-node field name or at-node data values to use for identifying watershed location.

  • nodata_value (float, optional) – Value that indicates an invalid value.

Examples

The example will use a 4,4 grid with node data values as illustrated:

-9999. -9999. -9999. -9999.
-9999.    67.     0. -9999.
-9999.    67.    67. -9999.
-9999. -9999. -9999. -9999.
>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 4))
>>> rmg.status_at_node
array([1, 1, 1, 1, 1, 0, 0, 1, 1, 0, 0, 1, 1, 1, 1, 1], dtype=uint8)
>>> node_data = np.array(
...     [
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, 67.0, 67.0, -9999.0],
...         [-9999.0, 67.0, 0.0, -9999.0],
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...     ]
... ).flatten()
>>> rmg.set_watershed_boundary_condition_outlet_coords(
...     (2, 2), node_data, -9999.0
... )
>>> rmg.status_at_node.reshape(rmg.shape)
array([[4, 4, 4, 4],
       [4, 0, 0, 4],
       [4, 0, 1, 4],
       [4, 4, 4, 4]], dtype=uint8)
set_watershed_boundary_condition_outlet_id(outlet_id, node_data, nodata_value=-9999.0)[source]#

Set the boundary conditions for a watershed. All nodes with nodata_value are set to BC_NODE_IS_CLOSED (4). All nodes with data values are set to BC_NODE_IS_CORE (0), with the exception that the outlet node is set to a BC_NODE_IS_FIXED_VALUE (1).

Note that the outer ring of the raster is set to BC_NODE_IS_CLOSED, even if there are nodes that have values. The only exception to this would be if the outlet node is on the boundary, which is acceptable.

Assumes that the id of the outlet is already known.

This assumes that the grid has a single watershed. If this is not the case this will not work.

Parameters:
  • outlet_id (integer) – id of the outlet node

  • node_data (field name or ndarray) – At-node field name or at-node data values to use for identifying watershed location.

  • nodata_value (float, optional) – Value that indicates an invalid value.

Examples

The example will use a 4,4 grid with node data values as illustrated:

-9999. -9999. -9999. -9999. -9999. 67. 0. -9999. -9999. 67. 67. -9999. -9999. -9999. -9999. -9999.

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((4, 4))
>>> rmg.status_at_node.reshape(rmg.shape)
array([[1, 1, 1, 1],
       [1, 0, 0, 1],
       [1, 0, 0, 1],
       [1, 1, 1, 1]], dtype=uint8)
>>> node_data = np.array(
...     [
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...         [-9999.0, 67.0, 67.0, -9999.0],
...         [-9999.0, 67.0, 0.0, -9999.0],
...         [-9999.0, -9999.0, -9999.0, -9999.0],
...     ]
... ).flatten()
>>> outlet = rmg.set_watershed_boundary_condition_outlet_id(
...     10, node_data, -9999.0
... )
>>> rmg.status_at_node.reshape(rmg.shape)
array([[4, 4, 4, 4],
       [4, 0, 0, 4],
       [4, 0, 1, 4],
       [4, 4, 4, 4]], dtype=uint8)
property shape#
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.

property spacing#
property status_at_corner#

Get array of the boundary status for each corner.

See also

status_at_node

property status_at_d8#
property status_at_diagonal#

Status at diagonals.

Examples

>>> from landlab import LinkStatus, NodeStatus, RasterModelGrid
>>> grid = RasterModelGrid((4, 3))

An inactive link (or diagonal) is one that joins two boundary nodes or has one end that is a closed boundary.

>>> grid.status_at_node
array([1, 1, 1,
       1, 0, 1,
       1, 0, 1,
       1, 1, 1], dtype=uint8)
>>> NodeStatus.CORE, NodeStatus.FIXED_VALUE
(<NodeStatus.CORE: 0>, <NodeStatus.FIXED_VALUE: 1>)

Diagonals that connect two boundary nodes are inactive.

>>> grid.status_at_diagonal
array([0, 4, 4, 0,
       0, 0, 0, 0,
       4, 0, 0, 4], dtype=uint8)
>>> LinkStatus.ACTIVE, LinkStatus.INACTIVE
(<LinkStatus.ACTIVE: 0>, <LinkStatus.INACTIVE: 4>)

By setting a node to closed that wasn’t before, a new link becomes inactive.

>>> grid.status_at_node[5] = NodeStatus.CLOSED
>>> grid.status_at_diagonal
array([0, 4, 4, 0,
       0, 0, 0, 4,
       4, 0, 0, 4], dtype=uint8)
property status_at_face#

Get array of the status of all faces.

See also

status_at_link

Get array of the status of all links.

Examples

>>> from landlab import RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.status_at_node[mg.nodes_at_left_edge] = mg.BC_NODE_IS_CLOSED
>>> mg.status_at_node[mg.nodes_at_right_edge] = mg.BC_NODE_IS_FIXED_GRADIENT
>>> mg.status_at_link
array([4, 4, 4, 4, 4, 0, 0, 0, 4, 4, 0, 0, 2, 4, 0, 0, 0, 4, 4, 0, 0,
       2, 4, 0, 0, 0, 4, 4, 4, 4, 4], dtype=uint8)
property status_at_node#

Get array of the boundary status for each node.

Examples

>>> import numpy as np
>>> from landlab import LinkStatus, NodeStatus, RasterModelGrid
>>> mg = RasterModelGrid((4, 5))
>>> mg.status_at_node.reshape((4, 5))
array([[1, 1, 1, 1, 1],
       [1, 0, 0, 0, 1],
       [1, 0, 0, 0, 1],
       [1, 1, 1, 1, 1]], dtype=uint8)
>>> np.any(mg.status_at_link == LinkStatus.FIXED)
False
>>> mg.status_at_node[mg.nodes_at_left_edge] = NodeStatus.FIXED_GRADIENT
>>> mg.status_at_node.reshape((4, 5))
array([[2, 1, 1, 1, 1],
       [2, 0, 0, 0, 1],
       [2, 0, 0, 0, 1],
       [2, 1, 1, 1, 1]], dtype=uint8)
>>> np.any(mg.status_at_link == LinkStatus.FIXED)  # links auto-update
True
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).

property unit_vector_at_corner#

Get a unit vector for each corner.

property unit_vector_at_face#

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

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 unit_vector_sum_xcomponent_at_corner#

Get array of x-component of unit vector sums at each corner.

property unit_vector_sum_xcomponent_at_node#

Get array of x-component of unit vector sums at each node.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3))
>>> len(grid.unit_vector_sum_xcomponent_at_node) == grid.number_of_nodes
True
>>> grid.unit_vector_sum_xcomponent_at_node
array([ 1.,  2.,  1.,  1.,  2.,  1.,  1.,  2.,  1.])
property unit_vector_sum_ycomponent_at_corner#

Get array of y-component of unit vector sums at each corner.

property unit_vector_sum_ycomponent_at_node#

Get array of y-component of unit vector sums at each node.

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 3))
>>> len(grid.unit_vector_sum_ycomponent_at_node) == grid.number_of_nodes
True
>>> grid.unit_vector_sum_ycomponent_at_node
array([ 1.,  1.,  1.,  2.,  2.,  2.,  1.,  1.,  1.])

Return an (nnodes, X) shape array of link IDs of which links are upwind of each node, according to values (field or array).

X is the maximum upwind links at any node. Nodes with fewer upwind links than this have additional slots filled with bad_index. Links are ordered anticlockwise from east.

Parameters:
  • values (str or array) – Name of variable field defined at links, or array of values at links.

  • bad_index (int) – Index to place in array indicating no link.

Returns:

Array of upwind link IDs

Return type:

ndarray

Examples

>>> import numpy as np
>>> 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,
... ]
>>> rmg.upwind_links_at_node("grad", bad_index=-1)
array([[-1, -1],
       [ 0, -1],
       [ 1, -1],
       [ 2, -1],
       [ 3, -1],
       [ 7,  4],
       [ 8,  5],
       [ 9,  6],
       [10, -1],
       [14, 11],
       [15, 12],
       [16, 13]])
property vertical_faces#
property x_of_corner#

Get x-coordinate of corner.

See also

x_of_node

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_cell#

Get the centroid of each cell.

See also

xy_of_patch

property xy_of_corner#

Get x and y-coordinates of corner.

See also

xy_of_node

property xy_of_face#
property xy_of_lower_left#

Return (x, y) of the reference point.

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_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 RasterModelGrid
>>> rmg = RasterModelGrid((4, 5), xy_of_reference=(12345, 678910))
>>> rmg.xy_of_reference
(12345, 678910)
>>> rmg.xy_of_reference = (98765, 43210)
>>> rmg.xy_of_reference
(98765, 43210)
property y_of_corner#

Get y-coordinate of corner.

See also

y_of_node

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"))
[]
grid_edge_is_closed_from_dict(boundary_conditions)[source]#

Get a list of closed-boundary status at grid edges.

Get a list that indicates grid edges that are closed boundaries. The returned list provides a boolean that gives the boundary condition status for edges order as [bottom, left, top, right].

boundary_conditions is a dict whose keys indicate edge location (as “bottom”, “left”, “top”, “right”) and values must be one of “open”, or “closed”. If an edge location key is missing, that edge is assumed to be open.

Parameters:

boundary_conditions (dict) – Boundary condition for grid edges.

Returns:

List of booleans indicating if an edge is a closed boundary.

Return type:

list

Examples

>>> from landlab.grid.raster import grid_edge_is_closed_from_dict
>>> grid_edge_is_closed_from_dict(dict(bottom="closed", top="open"))
[False, False, False, True]
>>> grid_edge_is_closed_from_dict({})
[False, False, False, False]