Mapping data between different grid elements#

Map values from one grid element to another.

Grid mapping functions#

map_link_head_node_to_link(grid, var_name[, out])

Map values from a link head nodes to links.

map_link_tail_node_to_link(grid, var_name[, out])

Map values from a link tail nodes to links.

map_min_of_link_nodes_to_link(grid, var_name)

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

map_max_of_link_nodes_to_link(grid, var_name)

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

map_mean_of_link_nodes_to_link(grid, var_name)

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

map_value_at_min_node_to_link(grid, ...[, out])

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_max_node_to_link(grid, ...[, out])

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_node_to_cell(grid, var_name[, out])

Map values for nodes to cells.

map_min_of_node_links_to_node(grid, var_name)

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

map_max_of_node_links_to_node(grid, var_name)

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

map_upwind_node_link_max_to_node(grid, var_name)

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

map_downwind_node_link_max_to_node(grid, ...)

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

map_upwind_node_link_mean_to_node(grid, var_name)

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

map_downwind_node_link_mean_to_node(grid, ...)

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

map_value_at_upwind_node_link_max_to_node(...)

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_value_at_downwind_node_link_max_to_node(...)

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_link_vector_components_to_node(grid, ...)

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

map_node_to_link_linear_upwind(grid, v, u[, out])

Assign values to links from upstream nodes.

map_node_to_link_lax_wendroff(grid, v, c[, out])

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

Each link has a tail and head node. The tail nodes are located at the start of a link, while the head nodes are located at end of a link.

Below, the numbering scheme for links in RasterModelGrid is illustrated with an example of a four-row by five column grid (4x5). In this example, each * (or X) is a node, the lines represent links, and the ^ and > symbols indicate the direction and head of each link. Link heads in the RasterModelGrid always point in the cardinal directions North (N) or East (E).:

 *--27-->*--28-->*--29-->*--30-->*
 ^       ^       ^       ^       ^
22      23      24      25      26
 |       |       |       |       |
 *--18-->*--19-->*--20-->*--21-->*
 ^       ^       ^       ^       ^
 13      14      15      16     17
 |       |       |       |       |
 *---9-->*--10-->X--11-->*--12-->*
 ^       ^       ^       ^       ^
 4       5       6       7       8
 |       |       |       |       |
 *--0--->*---1-->*--2--->*---3-->*

For example, node X has four link-neighbors. From south and going clockwise, these neighbors are [6, 10, 15, 11]. Both link 6 and link 10 have node X as their head node, while links 15 and 11 have node X as their tail node.

cartesian_to_polar(x, y)[source]#

Return 2d polar coordinates (r, theta) equivalent to given cartesian coordinates (x, y).

Examples

>>> r, theta = cartesian_to_polar(1.0, 1.0)
>>> int(r * 1000)
1414
>>> int(theta * 1000)
785

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

Map values from a link head nodes to links.

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

Map values from a link tail nodes to links.

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

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

Examples

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

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

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

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

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

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

Returns:

[x_component_of_link_vals_at_patch, y_component_of_link_vals_at_patch].

Return type:

len-2 list of arrays

Examples

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

Add a ring of closed nodes at the edge:

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

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

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

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

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

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

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

Now demonstrate the remaining functionality:

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

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at patches.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at patches.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at patches.

Return type:

ndarray

Examples

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

Map values for nodes to cells.

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

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

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

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

Returns:

Mapped values at cells.

Return type:

ndarray

Examples

>>> import numpy as np
>>> from landlab.grid.mappers import map_node_to_cell
>>> from landlab import RasterModelGrid
>>> rmg = RasterModelGrid((3, 4))
>>> _ = rmg.add_field("z", np.arange(12.0), at="node")
>>> map_node_to_cell(rmg, "z")
array([ 5.,  6.])
>>> values_at_cells = rmg.empty(at="cell")
>>> rtn = map_node_to_cell(rmg, "z", out=values_at_cells)
>>> values_at_cells
array([ 5.,  6.])
>>> rtn is values_at_cells
True

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

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

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

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

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

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

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

Values at links when c = -0.2:

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

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

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

Examples

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

Set values for the middle row of horizontal links.

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

Assign values to links from upstream nodes.

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

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

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

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

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

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

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

If u < 0, the link values should be:

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

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

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

Examples

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

Set values for the middle row of horizontal links.

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

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

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

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

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

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

Returns:

Mapped values at links.

Return type:

ndarray

Examples

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

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

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

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

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

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

Returns:

Mapped values at nodes.

Return type:

ndarray

Examples

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

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

Examples

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