rivers

river_directionality

Created on Tue Nov 6 14:31:01 2018

@author: Jon

rivgraph.rivers.river_directionality.dir_centerline(links, nodes, meshpolys, meshlines, Imask, gt, pixlen)[source]

Guess flow directions of links in a braided river channel.

Guesses the flow direction of links in a braided river channel network by exploiting a “valleyline” centerline. Two metrics are computed to help guess the correct direction. The first is the number of centerline transects (meshlines) that the link crosses. The second is the local angle of the centerline compared to the link’s angle. These metrics are appended to the links dictionary as links[‘cldist’] and links[‘clangs’].

Parameters:
  • links (dict) – Network links and associated properties.

  • nodes (dict) – Network nodes and associated properties.

  • meshpolys (list) – List of shapely.geometry.Polygons that define the valleyline mesh.

  • meshlines (list) – List of shapely.geometry.LineStrings that define the valleyline mesh.

  • Imask (np.array) – Binary mask of the network.

  • gt (tuple) – gdal-type GeoTransform of the original binary mask.

  • pixlen (float) – Length resolution of each pixel.

Returns:

links – Network links and associated properties with ‘cldists’ and ‘clangs’ attributes appended.

Return type:

dict

Guess link directions based on link widths.

Guesses each link’s direction based on link widths at the endpoints. The node with the larger width is guessed to be the upstream node. The ratio of guessed upstream to guessed downstream node widths is appended to links. If the width at a link’s endpoint is 0, the percent is set to 1000. This ratio is appended to the links dictionary as links[‘wid_pctdiff’].

Parameters:

links (dict) – Network links and associated properties.

Returns:

links – Network links and associated properties with ‘wid_pctdiff’ property appended.

Return type:

dict

rivgraph.rivers.river_directionality.directional_info(links, nodes, Imask, pixlen, exit_sides, gt, meshlines, meshpolys, Idt)[source]

Compute information for link direction setting.

Computes all the information required for link directions to be set for a river channel network.

Parameters:
  • links (dict) – Network links and associated properties.

  • nodes (dict) – Network nodes and associated properties.

  • Imask (np.array) – Binary mask of the network.

  • pixlen (float) – Length resolution of each pixel.

  • exit_sides (str) – Two-character string of cardinal directions denoting the upstream and downsteram sides of the image that the network intersects (e.g. ‘SW’).

  • gt (tuple) – gdal-type GeoTransform of the original binary mask.

  • meshlines (list) – List of shapely.geometry.LineStrings that define the valleyline mesh.

  • meshpolys (list) – List of shapely.geometry.Polygons that define the valleyline mesh.

  • Idt (np.array()) – Distance transform of Imask.

Returns:

  • links (dict) – Network links and associated properties including directional info.

  • nodes (dict) – Network nodes and associated properties including directional info.

rivgraph.rivers.river_directionality.fix_river_cycle(links, nodes, cyclelinks, cyclenodes, imshape)[source]

Attempt to fix a single cycle.

Attempts to resolve a single cycle within a river network. The general logic is that all link directions of the cycle are un-set except for those set by highly-reliable algorithms, and a modified direction-setting recipe is implemented to re-set these algorithms. This was developed according to the most commonly-encountered cases for real braided rivers, but could certainly be improved.

Parameters:
  • links (dict) – Network links and associated properties.

  • nodes (dict) – Network nodes and associated properties.

  • cyclelinks (list) – List of link ids that comprise a cycle.

  • cyclenodes (list) – List of node ids taht comprise a cycle.

  • imshape (tuple) – Shape of binary mask as (nrows, ncols).

Returns:

  • links (dict) – Network links and associated properties with the cycle fixed if possible.

  • nodes (dict) – Network nodes and associated properties with the cycle fixed if possible.

  • fixed (int) – 1 if the cycle was resolved, else 0.

rivgraph.rivers.river_directionality.fix_river_cycles(links, nodes, imshape)[source]

Attempt to resolve cycles in the network.

Attempts to resolve all cycles within the river network. This function is essentially a wrapper for fix_river_cycle(), which is where the heavy lifting is actually done. This function finds cycles, calls fix_river_cycle() on each one, then aggregates the results.

Parameters:
  • links (dict) – Network links and associated properties.

  • nodes (dict) – Network nodes and associated properties.

  • imshape (tuple) – Shape of binary mask as (nrows, ncols).

Returns:

  • links (dict) – Network links and associated properties with all possible cycles fixed.

  • nodes (dict) – Network nodes and associated properties with all possible cycles fixed.

  • cantfix_links (list of lists) – Contains link ids of unresolvable cycles. Length is equal to number of unresolvable cycles.

  • cantfix_nodes (TYPE) – Contains node ids of unresolvable cycles. Length is equal to number of unresolvable cycles.

rivgraph.rivers.river_directionality.re_set_linkdirs(links, nodes, imshape)[source]

Reset link directions.

Resets the link directions for a braided river channel network. This function is called to reset directions of links that belong to a cycle.

Parameters:
  • links (dict) – Network links and associated properties.

  • nodes (dict) – Network nodes and associated properties.

  • imshape (tuple) – Shape of binary mask as (nrows, ncols).

Returns:

  • links (dict) – Network links and associated properties with the directions re-set.

  • nodes (dict) – Network nodes and associated properties with the directions re-set.

rivgraph.rivers.river_directionality.set_directionality(links, nodes, Imask, exit_sides, gt, meshlines, meshpolys, Idt, pixlen, manual_set_csv)[source]

Set direction of each link.

This function sets the direction of each link within the network. It calls a number of helping functions and uses a somewhat-complicated logic to achieve this. The algorithms and logic is described in this open-access paper: https://esurf.copernicus.org/articles/8/87/2020/esurf-8-87-2020.pdf

Every time this is run, all directionality information is reset and recomputed. This includes checking for manually set links via the provided csv.

Parameters:
  • links (dict) – Network links and associated properties.

  • nodes (dict) – Network nodes and associated properties.

  • Imask (np.array) – Binary mask of the network.

  • exit_sides (str) – Two-character string of cardinal directions denoting the upstream and downsteram sides of the image that the network intersects (e.g. ‘SW’).

  • gt (tuple) – gdal-type GeoTransform of the original binary mask.

  • meshlines (list) – List of shapely.geometry.LineStrings that define the valleyline mesh.

  • meshpolys (list) – List of shapely.geometry.Polygons that define the valleyline mesh.

  • Idt (np.array()) – Distance transform of Imask.

  • pixlen (float) – Length resolution of each pixel.

  • manual_set_csv (str, optional) – Path to a user-provided csv file of known link directions. The default is None.

Returns:

  • links (dict) – Network links and associated properties with all directions set.

  • nodes (dict) – Network nodes and associated properties with all directions set.

rivgraph.rivers.river_directionality.set_unknown_cluster_by_widthpct(links, nodes)[source]

Set unknown links based on width differences at endpoints.

(flow goes wide->narrow)

river_utils

Created on Tue Nov 6 14:29:10 2018

@author: Jon

rivgraph.rivers.river_utils.centerline_mesh(coords, width_chan, meshwidth, grid_spacing, smoothing_param=1)[source]

Generate a centerline mesh.

Generates a centerline mesh. Differs from valleyline_mesh() in that it draws perpendiculars rather than offsetting the valley line to compute mesh polygons. This method is more effective for narrower channels that don’t require an exceptionally wide mesh (i.e. not much change).

Parameters:
  • coords – 2xN list, tuple, np.array (xs, ys) of coordinates defining centerline

  • width_chan – width of the river in same units of coords

  • meshwidth – how wide should the mesh be, in same units of coords

  • grid_spacing – how far apart should mesh cells be, in same units of coords

Returns:

  • transects (list of shapely.LineStrings) – the “perpendiculars” to the centerline used to generate the mesh

  • polys (list of shapely.Polygons) – coordinates of the polygons representing the grid cells of the mesh

  • cl_smooth (shapely.LineString) – the smoothed centerline used to compute transects

rivgraph.rivers.river_utils.compute_eBI(path_meshlines, path_links, method='local')[source]

method can be ‘local’ or ‘avg’

rivgraph.rivers.river_utils.find_inlet_outlet_nodes(links, nodes, exit_sides, Iskel)[source]

Append inlet and outlet nodes to the node dictionary.

Appends the inlet and outlet nodes to the nodes dictionary. Only works for rivers; deltas must be treated differently.

rivgraph.rivers.river_utils.mask_to_centerline(Imask, es)[source]

Extract centerline from a river mask.

This function takes an input binary mask of a river and extracts its centerline. If there are multiple channels (and therefore islands) in the river, they will be filled before the centerline is computed.

Note

The input mask should have the following properties:

  1. There should be only one “blob” (connected component)

  2. Where the blob intersects the image edges, there should be only one channel. This avoids ambiguity in identifying inlet/outlet links

Parameters:
  • Imask (ndarray) – the mask image (numpy array)

  • es (str) – two-character string comprinsed of “n”, “e”, “s”, or “w”. Exit sides correspond to the sides of the image that the river intersects. Upstream should be first, followed by downstream.

Returns:

  • dt.tif (geotiff) – geotiff of the distance transform of the binary mask

  • skel.tif (geotiff) – geotiff of the skeletonized binary mask

  • centerline.shp (shp) – shapefile of the centerline, arranged upstream to downstream

  • cl.pkl (pkl) – pickle file containing centerline coords, EPSG, and paths dictionary

rivgraph.rivers.river_utils.max_valley_width(Imask)[source]

Computes the maximum valley width of the input mask. Finds the single largest blob in the mask, fills its holes, then uses the distance transform to find the largest width.

Parameters:

Imask (np.array) – Binary mask from which the centerline was computed.

Returns:

max_valley_width – Maximum width of the channel belt, useful for computing a mesh. Units are pixels, so be careful to re-convert.

Return type:

float

rivgraph.rivers.river_utils.mirror_line_ends(xs, ys, npad)[source]

Reflect both ends of a line.

Reflects both ends of a line defined by x and y coordinates. The mirrored distance is set by npad, which refers to the number of vertices along the line to mirror.

rivgraph.rivers.river_utils.prune_river(links, nodes, exit_sides, Iskel, gdobj)[source]

Prune river network.

rivgraph.rivers.river_utils.valleyline_mesh(coords, avg_chan_width, buf_halfwidth, grid_spacing, smoothing=0.15)[source]

Generate a mesh over an input river centerline.

This function generates a mesh over an input river centerline. The mesh is generated across the valley, not just the channel extents, in order to perform larger-scale spatial analyses. With the correct parameter combinations, it can also be used to generate a mesh for smaller-scale analysis, but it is optimized for larger and strange behavior may occur.

Many plotting commands are commented out throughout this script as it’s still somewhat in beta mode.

Parameters:
  • coords – Nx2 list, tuple, or np.array of x,y coordinates. Coordinates MUST be in projected CRS for viable results.

  • width_chan – estimated width. Units MUST correspond to those of the input coordinates

  • buf_width – distance between centerline and left or right bufferline, in units of coords

  • grid_spacing – fraction of input centerline length that should be used for smoothing to create the valley centerline (between 0 and 1)

  • smoothing – fraction of centerline length to use for smoothing window

Returns:

  • transects (list of shapely.LineStrings) – the “perpendiculars” to the centerline used to generate the mesh

  • polys (list of shapely.Polygons) – coordinates of the polygons representing the grid cells of the mesh

  • cl_smooth (shapely.LineString) – the smoothed centerline used to compute transects