High Level API#

These high-level functions attempt to “do what the user means.” They group other, lower-level functions into common sets of operations. The “default workflow” is nearly entirely built from calls to this namespace. If the user needs to deviate from the default workflow, they can then call the lower level functions.

This top level module provides functionality for getting shapes and rasters representing watershed boundaries, river networks, digital elevation models, and other GIS datasets and then processing those data sets for use in simulations.

Functions here that get shapes and rasters should be used instead of directly using the source managers, because they additionally do things like convert coordinate systems, get the data into common data structures, etc.

watershed_workflow.elevate(m2: Mesh2D, dem: DataArray, **kwargs) None[source]#

Elevate a mesh onto the provided dem, in place.

Parameters:
  • mesh_crs (crs-type) – Mesh coordinate system.

  • dem (xr.DataArray) – 2D array forming an elevation raster.

  • algorithm (str, optional) – Algorithm used for interpolation. One of: * “nearest” for nearest-neighbor pixels * “piecewise bilinear” for interpolation (default)

Returns:

out – Array of triangle vertices, including a z-dimension.

Return type:

np.array((n_points, 3), ‘d’)

watershed_workflow.findHUC(source: Any, shape: Polygon, in_crs: CRS, hint: str, shrink_factor: float = 1e-05) str[source]#

Finds the smallest HUC containing shape.

Parameters:
  • source (source-type) – Source object providing a getShapes() method that gets HUCs by ID.

  • shape (shapely.geometry.Polygon) – Find this shape in a HUC.

  • in_crs (watershed_workflow.crs.CRS) – CRS of shape.

  • hint (str) – HUC in which to start searching. This should be at least as long as the indexing file size – HUC 2 or longer for WBD, 4 or longer for NHD Plus, or 8 or longer for NHD.

  • shrink_factor (float, optional) – A fraction of the radius of shape to shrink prior for checking containment within HUCs. This fixes cases where shape is on a HUC boundary with potentially some numerical error.

Returns:

out – The smallest containing HUC.

Return type:

str

watershed_workflow.getDatasetOnMesh(m2: Mesh2D, data: DataArray, **kwargs) ndarray[source]#

Interpolate xarray data onto cell centroids of a mesh.

watershed_workflow.getShapePropertiesOnMesh(m2: Mesh2D, df: GeoDataFrame, column: str, resolution: float, nodata: int | float | None = None, **kwargs) ndarray[source]#

Intepolate shape data onto cell centroids of a mesh.

watershed_workflow.reduceRivers(rivers: List[River], ignore_small_rivers: int = 0, keep_n_rivers: int = -1, prune_by_area: float = 0.0, area_property: str = 'drainage_area_sqkm', remove_diversions: bool = False, remove_braided_divergences: bool = False, tol: float | None = 0.1) List[River][source]#

Reduce the extent of the river network through a variety of methods.

Parameters:
  • rivers (list(river_tree.River)) – A list of rivers to reduce.

  • ignore_small_rivers (int, optional) – If provided and positive, removes rivers whose number of reaches is less than this value. If negative, keeps the N biggest (in number of reaches) rivers, where N is the negative of the provided value (e.g. -2 keeps the biggest 2 rivers).

  • prune_by_area (float, optional) – If provided, remove reaches whose total contributing area is less than this tol. NOTE: only valid for reaches that include a contributing area property (e.g. NHDPlus).

  • area_property (str, optional='DivergenceRoutedDrainAreaSqKm') – Name of the area property to use for determining reach CA. Note that this defines the units of prune_by_area value.

  • remove_diversions (bool, optional=False) – If true, remove diversions (see documentation of modify_rivers_remove_divergences()).

  • remove_braided_divergences (bool, optional=False) – If true, remove braided divergences (see documentation of modify_rivers_remove_divergences()).

  • tol (float, optional=0.1) – Defines what close is in the case of method == ‘geometry’

Returns:

out – A list of rivers, as River objects.

Return type:

list(river_tree.River)

watershed_workflow.simplify(hucs: SplitHUCs, rivers: List[River], reach_segment_target_length: float, huc_segment_target_length: float | None = None, river_close_distance: float = 100.0, river_far_distance: float = 500.0, resample_by_reach_property: bool = False, min_angle: float = 20, junction_min_angle: float = 20, snap_triple_junctions_tol: float | None = None, plot_diagnostics: bool = False, keep_points: bool = False) None[source]#

Simplifies, in place, the HUC and river shapes to create constrained, discrete segments.

Parameters:
  • hucs (SplitHUCs) – A split-form HUC object containing all reaches.

  • rivers (list(River)) – A list of river objects.

  • huc_segment_target_length (float) – Target length of a typical triangle edge away from the river.

  • reach_segment_target_length (float) – Target length of a typical triangle edge at the river.

  • waterbodies (list(shply), optional) – A list of waterbodies.

  • simplify_hucs (float, optional) – If provided, simply the hucs by moving points at most this many units (see also shapely.simplify). Units are that of the CRS of shapes.

  • simplify_rivers (float, optional) – If provided, simply the rivers by moving points at most this many units (see also shapely.simplify). Units are that of the CRS of shapes. If not provided, use the simplify_hucs value. Provide 0 to make no changes to the rivers.

  • simplify_waterbodies (float, optional) – Simplify the waterbodies. If not provided, uses the simplify_hucs value. Provide 0 to make no changes to the waterbodies.

  • prune_tol (float, optional = simplify_rivers) – Prune leaf reaches that are smaller than this tolerance. If not provided, uses simplify_rivers value. Provide 0 to not do this step.

  • merge_tol (float, optional = simplify_rivers) – Merges reaches that are smaller than this tolerance with their downstream parent reach. Note that if there is a branchpoint at the downstream node of the small reach, it will get moved to the upstream node. If not provided, uses simplify_rivers value. Provide 0 to not do this step.

  • snap_tol (float, optional = 0.75 * simplify_rivers) – Tolerance used for snapping rivers to nearby huc boundaries. Provide 0 to not snap.

  • snap_triple_junctions_tol (float, optional = 3 * snap_tol) – Tolerance used for snapping river triple junctions to huc triple junctions.

  • snap_reach_endpoints_tol (float, optional = 2 * snap_tol) – Tolerance used for snapping reach junctions to huc boundaries.

  • snap_waterbodies_tol (float, optional = snap_tol) – If not 0, snaps waterbody and HUC nodes that are nearly coincident to be discretely coincident. Note this is not recommended; prefer to include major waterbodies in the HUC network.

  • cut_intersections (bool, optional = False) – If true, force intersections of the river network and the HUC boundary to occur at a coincident node by adding nodes as needed.

watershed_workflow.tessalateRiverAligned(hucs: SplitHUCs, rivers: List[River], river_width: Any, internal_boundaries: List[River | BaseGeometry] | None = None, as_mesh: bool = True, debug: bool = False, **kwargs) Tuple[ndarray, List[List[int]]] | Mesh2D | Tuple[ndarray, List[List[int]], GeoDataFrame] | Tuple[ndarray, List[List[int]], ndarray, ndarray] | Tuple[Mesh2D, ndarray, ndarray][source]#

Tessalate HUCs using river-aligned quads along the corridor and triangles away from it.

Parameters:
  • hucs (SplitHUCs) – The huc geometry to tessalate. Note this will be adjusted if required by the river corridor.

  • rivers (list[River]) – The rivers to mesh with quads

  • river_width (float or dict or callable or boolean) – Width of the quads, either a float or a dictionary providing a {StreamOrder : width} mapping. Or a function (callable) that computer width using node properties Or boolean, where True means, width for each reach is explicitely provided properties as “width”

  • river_n_quads (int, optional) – Number of quads across the river. Currently only 1 is supported (the default).

  • hole_points (list(shapely.Point), optional) – List of points inside the polygons to be left as holes/voids (excluded from mesh)

  • internal_boundaries (list[shapely.Polygon], optional) – List of internal boundaries to embed in the domain, e.g. waterbodies.

  • diagnostics (bool, optional) – If true, prints extra diagnostic info.

  • ax (matplotlib Axes object, optional) – For debugging – plots troublesome reaches as quad elements are generated to find tricky areas.

  • kwargs – All other arguments are passed to the triangulation function for refinement.

Returns:

  • vertices (np.array((n_vertices, 2), ‘d’)) – Array of triangle vertices.

  • elems (list[list[int]]) – For each element, an ordered list of indices into the vertices array that make up that cell.

  • areas (_only if diagnostics=True_, np.array((n_cell_vertices), ‘d’)) – Array of areas.

watershed_workflow.triangulate(hucs: SplitHUCs, rivers: List[River] | None = None, internal_boundaries: List[shapely.geometry.BaseGeometry | River] | None = None, hole_points: List[shapely.geometry.Point] | None = None, diagnostics: bool = False, verbosity: int = 1, tol: float = 1.0, refine_max_area: float | None = None, refine_distance: Tuple[float, float, float, float] | None = None, refine_polygons: Tuple[List[shapely.geometry.Polygon], List[float]] | None = None, refine_max_edge_length: float | None = None, refine_min_angle: float | None = None, enforce_delaunay: bool = False, as_mesh: bool = True) Tuple[np.ndarray, np.ndarray] | watershed_workflow.mesh.Mesh2D | Tuple[np.ndarray, np.ndarray, np.ndarray, np.ndarray] | Tuple[watershed_workflow.mesh.Mesh2D, np.ndarray, np.ndarray][source]#

Triangulates HUCs and rivers.

Note, refinement of a given triangle is done if any of the provided criteria is met.

Parameters:
  • hucs (SplitHUCs) – A split-form HUC object from, e.g., get_split_form_hucs()

  • rivers (list[watershed_workflow.river_tree.River], optional) – List of rivers, used to refine the triangulation in conjunction with refine_distance.

  • internal_boundaries (list[shapely.geometry.Polygon, watershed_workflow.river_tree.River], optional) – List of objects, whose boundary (in the case of polygons/waterbodies) or reaches (in the case of River) will be present in the edges of the triangulation.

  • hole_points (list(shapely.Point), optional) – List of points inside the polygons to be left as holes/voids (excluded from mesh)

  • diagnostics (bool, optional) – Plot diagnostics graphs of the triangle refinement.

  • tol (float, optional) – Set tolerance for minimum distance between two nodes. The unit is the same as that of the watershed’s CRS. The default is 1.

  • refine_max_area (float, optional) – Refine a triangle if its area is greater than this area.

  • refine_polygons ([list(shapely.geometry.Polygon), list(float)], optional) – Refine a triangle if it falls within the polygons and its area is greater than the area limit for the polygon

  • refine_distance (list(float), optional) –

    Refine a triangle if its area is greater than a function of its centroid’s distance from the nearest point on the river network. The argument is given by:

    [near_distance, near_area, far_distance, far_area]

    Defining d as the distance from triangle centroid to the nearest point on the river network and area as the area of the triangle in question, refinement occurs if:

    • d < near_distance and area > near_area

    • d > far_distance and area > far_area

    • otherwise, defining d’ = (d - near_distance) / (far_distance - near_distance), refining occurs if area > near_area + (far_area - near_area) * d’

    Effectively this simply writes a piecewise linear function of triangle distance from centroid and uses that as a max area criteria.

  • refine_max_edge_length (float, optional) – Refine a triangle if its max edge length is greater than this length.

  • refine_min_angle (float, optional) – Try to ensure that all triangles have a minimum edge length greater than this value.

  • enforce_delaunay (bool,optional, experimental) –

    Attempt to ensure all triangles are proper Delaunay triangles.

Returns:

  • vertices (np.array((n_points, 2), ‘d’)) – Array of triangle vertices.

  • triangles (np.array((n_tris, 3), ‘i’)) – For each triangle, a list of 3 indices into the vertex array that make up that triangle.

  • areas (_only if diagnostics=True_, np.array((n_tris), ‘d’)) – Array of triangle areas.

  • distances (_only if diagnostics=True_, np.array((n_tris), ‘d’)) – Array of triangle distances from the river network.