landlab

API for landlab.grid.base

Python implementation of ModelGrid, a base class used to create and manage grids for 2D numerical models.

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.

BAD_INDEX_VALUE = -1

Indicates an index is, in some way, bad.

class ModelGrid(**kwds)[source]

Bases: landlab.field.graph_field.GraphFields, landlab.layers.eventlayers.EventLayersMixIn, landlab.layers.materiallayers.MaterialLayersMixIn

Base class for 2D structured or unstructured grids for numerical models.

The idea is to have at least two inherited classes, RasterModelGrid and DelaunayModelGrid, that can create and manage grids. To this might be added a GenericModelGrid, which would be an unstructured polygonal grid that doesn’t necessarily obey or understand the Delaunay triangulation, but rather simply accepts an input grid from the user. Also a HexModelGrid for hexagonal.

at_node

Values at nodes.

Type:dict-like
at_cell

Values at cells.

Type:dict-like

Values at links.

Type:dict-like
at_face

Values at faces.

Type:dict-like
at_grid

Global values

Type:dict-like
BAD_INDEX

Indicates a grid element is undefined.

Type:int
BC_NODE_IS_CORE

Indicates a node is core.

Type:int
BC_NODE_IS_FIXED_VALUE

Indicates a boundary node has a fixed value.

Type:int
BC_NODE_IS_FIXED_GRADIENT

Indicates a boundary node has a fixed gradient.

Type:int
BC_NODE_IS_LOOPED

Indicates a boundary node is wrap-around.

Type:int
BC_NODE_IS_CLOSED

Indicates a boundary node is closed

Type:int

Indicates a link is active, and can carry flux.

Type:int

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

Type:int

Indicates a link is inactive, and cannot carry flux.

Type:int
Other Parameters:
 
  • axis_name (tuple, optional) – Name of axes
  • axis_units (tuple, optional) – Units of coordinates
BAD_INDEX = -1

Indicates a node is bad index.

BC_LINK_IS_ACTIVE = 0

Indicates a link is active, and can carry flux

BC_LINK_IS_FIXED = 2

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

BC_LINK_IS_INACTIVE = 4

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.

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

LLCATS: NINF CONN BC

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

LLCATS: FINF BC

1=incoming flux, -1=outgoing flux, 0=no flux. Note that inactive links receive zero, 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 rows by max number of links per node. A zero indicates no link at this position.
Return type:(NODES, LINKS) ndarray of int

Examples

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((4, 3))
>>> grid.status_at_node[grid.nodes_at_left_edge] = grid.BC_NODE_IS_CLOSED
>>> grid.active_link_dirs_at_node # doctest: +NORMALIZE_WHITESPACE
array([[ 0,  0,  0,  0], [ 0, -1,  0,  0], [ 0,  0,  0,  0],
       [ 0,  0,  0,  0], [-1, -1,  0,  1], [ 0,  0,  1,  0],
       [ 0,  0,  0,  0], [-1, -1,  0,  1], [ 0,  0,  1,  0],
       [ 0,  0,  0,  0], [ 0,  0,  0,  1], [ 0,  0,  0,  0]],
       dtype=int8)

LLCATS: NINF LINF CONN

Type:Link flux directions at each node

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

LLCATS: LINF BC

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. / np.pi
>>> angles[0, :4]
array([  0.,  90.,  90.,  90.])
>>> angles[0, ::4]
array([ 0.,  0.,  0.])
>>> angles[0, ::5]
array([  0.,  45.,  45.])

LLCATS: NINF MEAS

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

LLCATS: NINF MEAS

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

LLCATS: LINF MEAS

as_dataset(include='*', exclude=None)[source]

Create an xarray Dataset representation of a grid.

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

at_cell = {}
at_corner = {}
at_face = {}
at_grid = {}
at_link = {}
at_node = {}
at_patch = {}
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')

LLCATS: GINF

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., 2.))
>>> 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')

LLCATS: GINF

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

LLCATS: NINF BC

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., 2.))
>>> 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.

LLCATS: NINF SURF

Calculate differences of node values over links.

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
>>> rmg = RasterModelGrid((3, 3))
>>> z = np.zeros(9)
>>> z[4] = 1.
>>> rmg.calc_diff_at_link(z)
array([ 0.,  0.,  0.,  1.,  0.,  1., -1.,  0., -1.,  0.,  0.,  0.])

LLCATS: LINF GRAD

calc_distances_of_nodes_to_point(coord, get_az=None, node_subset=None, out_distance=None, out_azimuth=None)[source]

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

LLCATS: NINF MEAS

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.

LLCATS: NINF GRAD

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.

LLCATS: NINF GRAD

Calculate gradients of node values at links.

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

Parameters:
  • node_values (ndarray or field name (x number of nodes)) – Values at grid nodes.
  • out (ndarray, optional (x number of links)) – Buffer to hold the result.
Returns:

Gradients across active links.

Return type:

ndarray

Examples

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

LLCATS: LINF GRAD

calc_grad_at_patch(elevs='topographic__elevation', ignore_closed_nodes=True, unit_normal=None, slope_magnitude=None)

Calculate the components of the gradient at each patch.

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

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.
  • unit_normal (array with shape (num_patches, 3) (optional)) – The unit normal vector to each patch, if already known.
  • slope_magnitude (array with size num_patches (optional)) – The slope of each patch, if already known.
Returns:

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

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.)
True
>>> np.allclose(x_grad, 0.)
True

LLCATS: PINF GRAD

calc_hillshade_at_node(alt=45.0, az=315.0, slp=None, asp=None, unit='degrees', elevs='topographic__elevation')[source]

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
  • slp and asp are both not specified, 'elevs' must be provided as (If) –
  • grid field name (defaults to 'topographic__elevation') or an (a) –
  • array of elevation values. In this case, the method will (nnodes-long) –
  • local slopes and aspects internally as part of the hillshade (calculate) –
  • 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.)
>>> z = mg.x_of_node * np.tan(60. * np.pi / 180.)
>>> mg.calc_hillshade_at_node(elevs=z, alt=30., az=210.)
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])

LLCATS: NINF SURF

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.

LLCATS: NINF GRAD

calc_slope_at_node(elevs='topographic__elevation', method='patch_mean', ignore_closed_nodes=True, return_components=False, **kwds)

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. It averages the gradient on each of the patches surrounding the node, creating a value for node slope that better incorporates nonlocal elevation information. Directional information can still be returned through use of the return_components keyword.

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  # not always true

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

Parameters:
  • elevs (str or ndarray, optional) – Field name or array of node values.
  • method ({'patch_mean', 'Horn'}) – By equivalence to the raster version, ‘patch_mean’ returns a scalar mean on the patches; ‘Horn’ returns a vector mean on the patches.
  • 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((4, 5))
>>> z = mg.node_x
>>> slopes = mg.calc_slope_at_node(elevs=z)
>>> np.allclose(slopes, 45. / 180. * np.pi)
True
>>> mg = RasterModelGrid((4, 5))
>>> z = - mg.node_y
>>> slope_mag, cmp = mg.calc_slope_at_node(elevs=z,
...                                        return_components=True)
>>> np.allclose(slope_mag, np.pi / 4.)
True
>>> np.allclose(cmp[0], 0.)
True
>>> np.allclose(cmp[1], - np.pi / 4.)
True
>>> mg = RadialModelGrid(n_rings=3)
>>> z = mg.radius_at_node
>>> slope_at_node = np.round(mg.calc_slope_at_node(elevs=z), decimals=5)
>>> nodes_at_ring = [
...     np.where(np.isclose(mg.radius_at_node, radius)) for radius in range(3)
... ]
>>> slope_at_node[nodes_at_ring[0]]
array([ 0.85707])
>>> slope_at_node[nodes_at_ring[1]]
array([ 0.79417,  0.79417,  0.79417,  0.79417,  0.79417,  0.79417])
>>> slope_at_node[nodes_at_ring[2]]
array([ 0.77542,  0.78453,  0.78453,  0.77542,  0.77542,  0.78453,
        0.78453,  0.77542,  0.77542,  0.78453,  0.78453,  0.77542])

LLCATS: NINF GRAD SURF

calc_slope_at_patch(elevs='topographic__elevation', ignore_closed_nodes=True, unit_normal=None)

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

If ignore_closed_nodes is True, closed nodes do not affect slope calculations. If a 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.
  • unit_normal (array with shape (num_patches, 3) (optional)) – The unit normal vector to each patch, if already known.
Returns:

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

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

LLCATS: PINF GRAD

calc_unit_normal_at_patch(elevs='topographic__elevation')

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

Parameters:elevs (str or ndarray, optional) – Field name or array of node values.
Returns:nhat – The unit normal vector <a, b, c> to each patch.
Return type:num-patches x length-3 array

Examples

>>> from landlab import HexModelGrid
>>> mg = HexModelGrid((3, 3))
>>> z = mg.node_x * 3. / 4.
>>> mg.calc_unit_normal_at_patch(z)
array([[-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8],
       [-0.6,  0. ,  0.8]])

LLCATS: PINF GRAD

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

LLCATS: CINF NINF CONN

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

LLCATS: NINF BC

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

LLCATS: CINF BC

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

LLCATS: NINF BC

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'] = np.array([-1., -2., -1.,
...                                 -2., -3., -4., -5.,
...                                 -1., -2., -1.,
...                                 -1., -2., -3., -4.,
...                                 -1., -2., -1.])
>>> 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]])

LLCATS: LINF NINF CONN

fields(include='*', exclude=None)[source]

List of fields held by the grid.

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

LLCATS: NINF BC

Get array of fixed links.

Examples

>>> from landlab import NodeStatus, RasterModelGrid
>>> grid = RasterModelGrid((3, 4))
>>> grid.status_at_node # doctest: +NORMALIZE_WHITESPACE
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 # doctest: +NORMALIZE_WHITESPACE
array([2, 2, 2, 2,
       1, 0, 0, 1,
       1, 1, 1, 1], dtype=uint8)
>>> grid.fixed_links
array([4, 5])

LLCATS: LINF BC

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

LLCATS: NINF BC

classmethod from_dict(params)[source]

Create grid from dictionary.

Parameters:params (dictionary) – Dictionary of required parameters to create a model grid.

Examples

>>> from landlab import RasterModelGrid
>>> params = {"shape": (3,4), "xy_spacing": 2}
>>> grid = RasterModelGrid.from_dict(params)
>>> 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_file(file_like)[source]

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

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'] = np.array([-1., -2., -1.,
...                                 -2., -3., -4., -5.,
...                                 -1., -2., -1.,
...                                 -1., -2., -3., -4.,
...                                 -1., -2., -1.])
>>> 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)

LLCATS: LINF NINF CONN

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'] = np.array([-1., -2., -1.,
...                                 -2., -3., -4., -5.,
...                                 -1., -2., -1.,
...                                 -1., -2., -3., -4.,
...                                 -1., -2., -1.])
>>> 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)

LLCATS: LINF NINF CONN

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

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'] = np.array([-1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.])
>>> 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

LLCATS: NINF LINF MAP

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'] = np.array([-1., -2., -1.,
...                                 -2., -3., -4., -5.,
...                                 -1., -2., -1.,
...                                 -1., -2., -3., -4.,
...                                 -1., -2., -1.])
>>> 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

LLCATS: NINF LINF MAP

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

LLCATS: NINF LINF MAP

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

LLCATS: NINF LINF MAP

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

Examples

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

Map the 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., -1.,  1., -1., -1., -1., -1.])

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

LLCATS: PINF LINF MAP

Map data defined on links to nodes.

Given a variable defined on links, breaks it into x and y components and assigns values to nodes by averaging each node’s attached links.

Parameters:q (ndarray of floats (1D, length = number of links in grid)) – Variable defined on links
Returns:x and y components of variable mapped to nodes (1D, length = number of nodes)
Return type:ndarray, ndarray

See also

_create_link_unit_vectors
sets up unit vectors at links and unit-vector sums at nodes

Notes

THIS ALGORITHM IS NOT CORRECT AND NEEDS TO BE CHANGED!

The concept here is that q contains a vector variable that is defined at each link. The magnitude is given by the value of q, and the direction is given by the orientation of the link, as described by its unit vector.

To map the link-vector values to the nodes, we break the values into x- and y-components according to each link’s unit vector. The x-component of q at a node is a weighted sum of the x-components of the links that are attached to that node. A good way to appreciate this is by example. Consider a 3x4 raster grid:

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

Imagine that for each node, we were to add up the unit vector components for each connected link; in other words, add up all the x components of the unit vectors associated with each link, and add up all the y components. Here’s what that would look like for the above grid (“vsx” and “vsy” stand for “vector sum x” and “vector sum y”):

  • Corner nodes (0, 3, 8, 11): vsx = 1, vsy = 1
  • Bottom and top nodes (1-2, 9-10): vsx = 2, vsy = 1
  • Left and right nodes (4, 7): vsx = 1, vsy = 2
  • All others: vsx = 2, vsy = 2

The process of creating unit-vector sums at nodes is handled by ModelGrid._create_link_unit_vectors() (and, for raster grids, by the overriding method RasterModelGrid._create_link_unit_vectors()). The node unit-vector sums are then stored in self.node_unit_vector_sum_x and self.node_unit_vector_sum_y.

How would you use this? Suppose you have a vector variable q defined at links. What’s the average at the nodes? We’ll define the average as follows. The terminology here is: \(q = (u,v)\) represents the vector quantity defined at links, \(Q = (U,V)\) represents its definition at nodes, \((m,n)\) represents the unit vector components at a link, and \((S_x,S_y)\) represents the unit-vector sum at a given node.

\[U_i = \sum_{j=1}^{L_i} q_j m_j / S_{xi} V_i = \sum_{j=1}^{L_i} q_j n_j / S_{yi}\]

Suppose that the vector q is uniform and equal to one. Then, at node 0 in the above grid, this works out to:

U_0 = (q_0 m_0) / 1 + (q_8 m_8) / 1 = (1 0)/ 1 + (1 1)/1 = 1
V_0 = (q_0 n_0) / 1 + (q_8 n_8) / 1 = (1 1) / 1 + (1 0) / 1 = 1

At node 1, in the bottom row but not a corner, we add up the values of q associated with THREE links. The x-vector sum of these links is 2 because there are two horizontal links, each with an x- unit vector value of unity. The y-vector sum is 1 because only one of the three (link #1) has a non-zero y component (equal to one). Here is how the numbers work out:

U_1 = (q_1 m_1) / 2 + (q_8 m_8) / 2 + (q_9 m_9) / 2
    = (1 0) / 2 + (1 1) / 2 + (1 1) / 2 = 1
V_1 = (q_1 n_1) / 1 + (q_8 n_8) / 1 + (q_9 n_9) / 1
    = (1 1) / 1 + (1 0) / 1 + (1 0) / 1 = 1

At node 5, in the interior, there are four connected links (two in-links and two out-links; two horizontal and two vertical). So, we add up the q values associated with all four:

U_5 = (q_1 m_1) / 2 + (q_5 m_5) / 2 + (q_11 m_11) / 2 + (q_12 m_12) / 2
    = (1 0) / 2 + (1 0) / 2 + (1 1) / 2 + (1 1) / 2 = 1

V_5 = (q_1 n_1) / 2 + (q_5 n_5) / 2 + (q_11 n_11) / 2 + (q_12 n_12) / 2
    = (1 1) / 2 + (1 1) / 2 + (1 0) / 2 + (1 0) / 2 = 1

To do this calculation efficiently, we use the following algorithm:

FOR each row in _node_inlink_matrix (representing one inlink @ each
node)
    Multiply the link's q value by its unit x component ...
    ... divide by node's unit vector sum in x ...
    ... and add it to the node's total q_x
    Multiply the link's q value by its unit y component ...
    ... divide by node's unit vector sum in y ...
    ... and add it to the node's total q_y

Examples

Example 1

q[:] = 1. Vector magnitude is \(\sqrt{2}\), direction is \((1,1)\).

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((3, 4), xy_spacing=(2., 2.))
>>> grid.unit_vector_at_node
array([[ 1.,  1.],
       [ 2.,  1.],
       [ 2.,  1.],
       [ 1.,  1.],
       [ 1.,  2.],
       [ 2.,  2.],
       [ 2.,  2.],
       [ 1.,  2.],
       [ 1.,  1.],
       [ 2.,  1.],
       [ 2.,  1.],
       [ 1.,  1.]])
>>> q = grid.ones(at='link')
>>> grid.map_link_vector_to_nodes(q)
array([[ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.],
       [ 1.,  1.]])

Example 2

Vector magnitude is 5, angle is 30 degrees from horizontal, forming a 3-4-5 triangle.

>>> import numpy as np
>>> q = np.array([4., 4., 4., 3., 3., 3., 3.,
...               4., 4., 4., 3., 3., 3., 3.,
...               4., 4., 4])
>>> grid.map_link_vector_to_nodes(q)
array([[ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.],
       [ 4.,  3.]])

..todo:

Fix and finish example 3 below.

Example 3: Hexagonal grid with vector as above. Here, q is pre-calculated to have the right values to represent a uniform vector with magnitude 5 and orientation 30 degrees counter-clockwise from horizontal.

LLCATS: NINF LINF CONN MAP

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

LLCATS: NINF LINF MAP

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

LLCATS: NINF LINF MAP

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'] = np.array([5., 4., 3., 2.,
...                                 3., 4., 3., 2.,
...                                 3., 2., 1., 0.])
>>> map_max_of_patch_nodes_to_patch(rmg, 'vals')
array([ 5., 4., 3.,
        4., 4., 3.])
>>> rmg.at_node['vals'] = np.array([5., 4., 3., 2.,
...                                 3., 4., 3., 2.,
...                                 3., 2., 1., 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 # doctest: +NORMALIZE_WHITESPACE
array([ 5., 4., 0.,
        4., 4., 0.])

LLCATS: PINF NINF MAP

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

LLCATS: NINF LINF MAP

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'] = np.array([5., 4., 3., 2.,
...                                 5., 4., 3., 2.,
...                                 3., 2., 1., 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'] = np.array([5., 4., 3., 2.,
...                                 5., 4., 3., 2.,
...                                 3., 2., 1., 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 # doctest: +NORMALIZE_WHITESPACE
array([ 4.5, 4. , 0. ,
        3.5, 3. , 0. ])

LLCATS: PINF NINF MAP

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

LLCATS: NINF LINF MAP

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

LLCATS: NINF LINF MAP

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'] = np.array([5., 4., 3., 2.,
...                                 5., 4., 3., 2.,
...                                 3., 2., 1., 0.])
>>> map_min_of_patch_nodes_to_patch(rmg, 'vals')
array([ 4., 3., 2.,
        2., 1., 0.])
>>> rmg.at_node['vals'] = np.array([5., 4., 3., 2.,
...                                 5., 4., 3., 2.,
...                                 3., 2., 1., 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 # doctest: +NORMALIZE_WHITESPACE
array([ 4., 4., 0.,
        2., 2., 0.])

LLCATS: PINF NINF MAP

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

LLCATS: CINF NINF MAP

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'] = np.array([-1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.])
>>> map_upwind_node_link_max_to_node(rmg, 'grad')
array([ 0.,  1.,  2.,  1.,
        0.,  1.,  2.,  1.,
        0.,  1.,  2.,  1.])
>>> 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
array([ 0.,  1.,  2.,  1.,
        0.,  1.,  2.,  1.,
        0.,  1.,  2.,  1.])
>>> rtn is values_at_nodes
True

LLCATS: NINF LINF MAP

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'] = np.array([-1., -2., -1.,
...                                 -2., -3., -4., -5.,
...                                 -1., -2., -1.,
...                                 -1., -2., -3., -4.,
...                                 -1., -2., -1.])
>>> 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

LLCATS: NINF LINF MAP

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'] = np.array([-1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.])
>>> 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

LLCATS: NINF LINF MAP

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

LLCATS: NINF LINF MAP

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

LLCATS: NINF LINF MAP

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'] = np.array([-1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.,
...                                  0.,  0.,  0.,  0.,
...                                 -1., -2., -1.])
>>> 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

LLCATS: NINF LINF MAP

ndim

Number of spatial dimensions of the grid.

LLCATS: GINF

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

LLCATS: NINF CINF BC CONN

node_axis_coordinates(axis=0)[source]

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) # doctest: +NORMALIZE_WHITESPACE
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) # doctest: +NORMALIZE_WHITESPACE
array([ 0., 1., 2., 3., 4.,
        0., 1., 2., 3., 4.,
        0., 1., 2., 3., 4.,
        0., 1., 2., 3., 4.])

LLCATS: GINF NINF MEAS

node_has_boundary_neighbor()[source]

Check if ModelGrid 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.

Parameters:ids (int, or iterable of int) – ID of node to test.
Returns:True if node has a neighbor with a boundary ID, False otherwise.
Return type:boolean

Examples

    0,  1,  2,  3,
  4,  5,  6,  7,  8,
9, 10,  11, 12, 13, 14,
  15, 16, 17, 18, 19,
    20, 21, 22, 23
>>> from landlab import HexModelGrid
>>> grid = HexModelGrid((5, 4))
>>> grid.node_has_boundary_neighbor()
array([ True,  True,  True,  True,  True,  True,  True,  True,  True,
        True,  True, False, False,  True,  True,  True,  True,  True,
        True,  True,  True,  True,  True,  True], dtype=bool)
>>> grid.node_has_boundary_neighbor()[6]
True
>>> grid.node_has_boundary_neighbor()[12]
False
>>> grid.node_has_boundary_neighbor()[((12, 0),)]
array([False,  True], dtype=bool)

LLCATS: NINF CONN BC

node_is_boundary(ids, boundary_flag=None)[source]

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)

LLCATS: NINF BC

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

LLCATS: FINF BC

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

LLCATS: LINF BC

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

LLCATS: CINF BC

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

LLCATS: NINF BC

number_of_elements(name)[source]

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

LLCATS: GINF

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

LLCATS: LINF BC

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

LLCATS: PINF LINF BC

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

LLCATS: PINF NINF BC

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

LLCATS: NINF BC

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

LLCATS: PINF LINF

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

LLCATS: PINF NINF

reset_status_at_node()[source]

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

LLCATS: LINF

set_nodata_nodes_to_closed(node_data, nodata_value)[source]

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
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., -9999, 0., 0., 0.])
>>> mg.set_nodata_nodes_to_closed(h, -9999)
>>> mg.status_at_node
array([4, 4, 4, 4, 4, 4, 0, 1, 4, 1, 1, 1], dtype=uint8)

LLCATS: BC NINF

set_nodata_nodes_to_fixed_gradient(node_data, nodata_value)[source]

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 # doctest: +NORMALIZE_WHITESPACE
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 = rmg.zeros(at='node')
>>> z = np.array([
...     -99., -99., -99., -99., -99., -99., -99., -99., -99.,
...     -99., -99., -99.,   0.,   0.,   0.,   0.,   0., -99.,
...     -99., -99., -99.,   0.,   0.,   0.,   0.,   0., -99.,
...     -99., -99., -99., -99., -99., -99., -99., -99., -99.])
>>> rmg.set_nodata_nodes_to_fixed_gradient(z, -99)
>>> rmg.status_at_node # doctest: +NORMALIZE_WHITESPACE
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 # doctest: +NORMALIZE_WHITESPACE
array([4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 2, 2, 2, 2, 2, 4,
       4, 4, 2, 0, 0, 0, 0, 2, 4, 4, 4, 0, 0, 0, 0, 0, 4,
       4, 4, 2, 0, 0, 0, 0, 2, 4, 4, 4, 2, 2, 2, 2, 2, 4,
       4, 4, 4, 4, 4, 4, 4, 4], dtype=uint8)

LLCATS: BC NINF

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 # doctest: +NORMALIZE_WHITESPACE
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)

LLCATS: BC LINF

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

LLCATS: NINF BC

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

LLCATS: NINF MEAS

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

LLCATS: NINF MEAS

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'] = np.array([-1., -2., -1.,
...                                 -2., -3., -4., -5.,
...                                 -1., -2., -1.,
...                                 -1., -2., -3., -4.,
...                                 -1., -2., -1.])
>>> 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]])

LLCATS: LINF NINF CONN

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)

Separate a pair of links with vector values into x and y components.

The concept here is that a pair of adjacent links attached to a node are projections of a ‘true’ but unknown vector. This function finds and returns the x and y components of this true vector. The trivial case is the situation in which the two links are orthogonal and aligned with the grid axes, in which case the vectors of these two links are the x and y components.

Parameters:
  • L2 (L1,) – Values (magnitudes) associated with the two links
  • b1y, b2x, b2y (b1x,) – Unit vectors of the two links
Returns:

ax, ay – x and y components of the ‘true’ vector

Return type:

float

Notes

The function does an inverse vector projection. Suppose we have a given ‘true’ vector \(a\), and we want to project it onto two other lines with unit vectors (b1x,b1y) and (b2x,b2y). In the context of Landlab, the ‘true’ vector is some unknown vector quantity, which might for example represent the local water flow velocity. The lines represent two adjacent links in the grid.

Let \(\mathbf{a}\) be the true vector, \(\mathbf{B}\) be a different vector with unit vector \(\mathbf{b}\), and \(L\) be the scalar projection of a onto B. Then,

..math:

L = \mathbf{a} \dot \mathbf{b} = a_x b_x + a_y b_y,

where \((a_x,a_y)\) are the components of a and \((b_x,b_y)\) are the components of the unit vector b.

In this case, we know b (the link unit vector), and we want to know the x and y components of a. The problem is that we have one equation and two unknowns (\(a_x\) and \(a_y\)). But we can solve this if we have two vectors, both of which are projections of a. Using the subscripts 1 and 2 to denote the two vectors, we can obtain equations for both \(a_x\) and \(a_y\):

..math:

a_x = L_1 / b_{1x} - a_y b_{1y} / b_{1x}

a_y = L_2 / b_{2y} - a_x b_{2x} / b_{2y}

Substituting the second into the first,

..math:

a_x = [L_1/b_{1x}-L_2 b_{1y}/(b_{1x} b_{2y})] / [1-b_{1y} b_{2x}/(b_{1x} b_{2y})]

Hence, we find the original vector \((a_x,a_y)\) from two links with unit vectors \((b_{1x},b_{1y})\) and \((b_{2x},b_{2y})\) and associated values \(L_1\) and \(L_2\).

Note that the above equations require that \(b_{1x}>0\) and \(b_{2y}>0\). If this isn’t the case, we invert the order of the two links, which requires \(b_{2x}>0\) and \(b_{1y}>0\). If none of these conditions is met, then we have a degenerate case.

Examples

The following example represents the active links in a 7-node hexagonal grid, with just one core node. The ‘true’ vector has a magnitude of 5 units and an orientation of 30 degrees, pointing up and to the right (i.e., the postive-x and postive-y quadrant), so that its vector components are 4 (x) and 3 (y) (in other words, it is a 3-4-5 triangle). The values assigned to L below are the projection of that true vector onto the six link vectors. The algorithm should recover the correct vector component values of 4 and 3. The FOR loop examines each pair of links in turn.

>>> import numpy as np
>>> from landlab.grid.base import find_true_vector_from_link_vector_pair
>>> bx = np.array([0.5, -0.5, -1., -0.5, 1., 0.5])
>>> by = np.array([0.866, 0.866, 0., -0.866, 0., -0.866])
>>> L = np.array([4.6, 0.6, -4., -4.6, 4., -0.6])
>>> for i in range(5):
...     ax, ay = find_true_vector_from_link_vector_pair(
...         L[i], L[i+1], bx[i], by[i], bx[i+1], by[i+1])
...     round(ax,1), round(ay,1)
(4.0, 3.0)
(4.0, 3.0)
(4.0, 3.0)
(4.0, 3.0)
(4.0, 3.0)