Meshing

Mesh

Tools for turning data into meshes, then writing them to file.

Works with and assumes all polyhedra cells (and polygon faces).

Requires building a reasonably recent version of Exodus to get the associated exodus.py wrappers in python3.

Note that this is typically done in your standard ATS installation, assuming you have built your Amanzi TPLs with shared libraries (the default through bootstrap).

In that case, simply ensure that ${AMANZI_TPLS_DIR}/SEACAS/lib is in your PYTHONPATH.

class watershed_workflow.mesh.SideSet(name: str, setid: int, cell_list: list[int], side_list: list[int])[source]

A collection of faces in cells.

class watershed_workflow.mesh.LabeledSet(name: str, setid: int, entity: str, ent_ids: list[int], to_extrude: bool = False)[source]

A generic collection of entities.

class watershed_workflow.mesh.Mesh2D(coords, conn, labeled_sets=_Nothing.NOTHING, crs=None, eps=0.001, check_handedness=True, validate=False)[source]

A 2D mesh class.

Parameters:
  • coords (np.ndarray(NCOORDS,NDIMS)) – Array of coordinates of the 2D mesh.

  • conn (list(lists)) – List of lists of indices into coords that form the cells.

  • labeled_sets (list(LabeledSet), optional) – List of labeled sets to add to the mesh.

  • crs (CRS) – Keep this as a property for future reference.

  • eps (float, optional=0.01) – A small measure of length between coords.

  • check_handedness (bool, optional=True) – If true, ensure all cells are oriented so that right-hand rule ordering of the nodes points up.

  • validate (bool, optional=False) – If true, validate coordinates and connections post-construction.

  • Note ((coords, conn) may be output provided by a) –

  • call. (watershed_workflow.triangulation.triangulate()) –

property conn

Note that conn is immutable because changing this breaks all the other properties, which may be cached. To change topology, one must construct a new mesh!

static edge_hash(i, j)[source]

Hashes edges in a direction-independent way.

property boundary_edges

Return edges in the boundary of the mesh, ordered around the boundary.

check_handedness()[source]

Ensures all cells are oriented via the right-hand-rule, i.e. in the +z direction.

next_available_labeled_setid()[source]

Returns next available LS id.

compute_centroid(c)[source]

Computes, based on coords, the centroid of a cell with ID c.

Note this ALWAYS recomputes the value, not using the cache.

property centroids

Calculate surface mesh centroids.

clear_geometry_cache()[source]

If coordinates are changed, any computed, cached geometry must be recomputed. It is the USER’s responsibility to call this function if any coords are changed!

plot(color=None, ax=None)[source]

Plot the flattened 2D mesh.

write_vtk(filename)[source]

Writes to VTK.

transform(mat=None, shift=None)[source]

Transform a 2D mesh

classmethod read_VTK(filename)[source]

Constructor from a VTK file.

classmethod read_VTK_Unstructured(filename)[source]

Constructor from an unstructured VTK file.

classmethod read_VTK_Simplices(filename)[source]

Constructor from an structured VTK file.

classmethod from_Transect(x, z, width=1, **kwargs)[source]

Creates a 2D surface strip mesh from transect data

to_dual()[source]

Creates a 2D surface mesh from a primal 2D triangular surface mesh.

Returns:

  • dual_nodes (np.ndarray)

  • dual_conn (list(lists)) – The nodes and cell_to_node_conn of a 2D, polygonal, nearly Voronoi mesh that is the truncated dual of self. Here we say nearly because the boundary triangles (see note below) are not Voronoi.

  • dual_from_primal_mapping (np.array( (len(dual_conn),), ‘i’)) – Mapping from a dual cell to the primal node it is based on. Without the boundary, this would be simply arange(0,len(dual_conn)), but beacuse we added triangles on the boundary, this mapping is useful.

Note

  • Nodes of the dual are on circumcenters of the primal triangles.

  • Centers of the dual are numbered by nodes of the primal mesh, modulo the boundary.

  • At the boundary, the truncated dual polygon may be non-convex. To avoid this issue, we triangulate boundary polygons.

We do not simply construct a mesh because likely the user needs to elevate, generate material IDs, labeled sets, etc, prior to mesh construction.

class watershed_workflow.mesh.Mesh3D(coords, face_to_node_conn, cell_to_face_conn, labeled_sets=_Nothing.NOTHING, side_sets=_Nothing.NOTHING, material_ids=None, crs=None, eps=0.001, validate=False)[source]

A 3D mesh class.

Parameters:
  • coords (np.ndarray(NCOORDS,3)) – Array of coordinates of the 3D mesh.

  • face_to_node_conn (list(lists)) – List of lists of indices into coords that form the faces.

  • cell_to_face_conn (list(lists)) – List of lists of indices into face_to_node_conn that form the cells.

  • side_sets (list(SideSet), optional) – List of side sets to add to the mesh.

  • labeled_sets (list(LabeledSet), optional) – List of labeled sets to add to the mesh.

  • material_ids (np.array((len(cell_to_face_conn),),'i'), optional) – Array of length num_cells that specifies material IDs

  • crs (CRS) – Keep the coordinate system for reference.

  • eps (float, optional=0.01) – A small measure of length between coords.

  • (coords (Note that) –

  • a (conn) may be output provided by) –

  • call. (watershed_workflow.triangulation.triangulate()) –

validate()[source]

Checks the validity of the mesh, or throws an AssertionError.

next_available_labeled_setid()[source]

Returns next available LS id.

write_vtk(filename)[source]

Writes to VTK.

Note, this just writes the topology/geometry information, for WEDGE type meshes (extruded triangles). No labeled sets are written. Prefer to use write_exodus() for a fully featured mesh.

write_exodus(filename, face_block_mode='one block')[source]

Write the 3D mesh to ExodusII using arbitrary polyhedra spec

static summarize_extrusion(layer_types, layer_data, ncells_per_layer, mat_ids, surface_cell_id=0)[source]

Summarizes extruded data by printing info to log file.

This is useful in rapidly debugging and understanding the layering before you do the extrusion process.

classmethod extruded_Mesh2D(mesh2D, layer_types, layer_data, ncells_per_layer, mat_ids)[source]

Uniformly extrude a 2D mesh to make a 3D mesh.

Layers of potentially multiple sets of cells are extruded downward in the vertical. The cell dz is uniform horizontally and vertically through the layer in all but the ‘snapped’ case, where it is uniform vertically but not horizontally.

Parameters:
  • mesh2D (Mesh2D) – The 2D mesh to extrude

  • layer_types (str, list[str]) – One of [‘snapped’, ‘function’, ‘node’, ‘cell’,] or a list of these strings, describing the extrusion method for each layer. If only a single string is supplied, all layers are of this type. See layer_data below.

  • layer_data (list[data]) –

    Data required, one for each layer. The data is specific to the layer type:

    • ’constant’ : float, layer thickness

    • ’snapped’ : float, the bottom z-coordinate of the layer. Cells are extruded uniformly to whatever thickness required to match this bottom coordinate.

    • ’function’ : function, layer_thickness = func(x,y)

    • ’node’ : np.array((mesh2D.num_nodes,), float) layer thickness for each node

    • ’cell’ : np.array((mesh2D.num_cells,), float) interpolates nodal layer thickness from neighboring cell thicknesses

  • ncells_per_layer (int, list[int]) – Either a single integer (same number of cells in all layers) or a list of number of cells in each layer

  • mat_ids (int, list[int], np.ndarray((num_layers, mesh2D.num_cells), dtype=int)) – Material ID for each cell in each layer. If an int, all cells in all layers share thesame material ID. If a list or 1D array, one material ID per layer. If a 2D array, provides all material IDs explicitly.

Returns:

The extruded, 3D mesh.

Return type:

Mesh3D

watershed_workflow.mesh.telescope_factor(ncells, dz, layer_dz)[source]

Calculates a telescoping factor to fill a given layer.

Calculates a constant geometric factor, such that a layer of thickness layer_dz is perfectly filled by ncells in the vertical, where the top cell is dz in thickness and each successive cell grows by a factor of that factor.

Parameters:
  • ncells (int) – Number of cells (in the vertical) needed.

  • dz (float) – Top cell’s thickness in the vertical.

  • layer_dz (float) – Thickness of the total layer.

Returns:

The telescoping factor.

Return type:

float

watershed_workflow.mesh.optimize_dzs(dz_begin, dz_end, thickness, num_cells, p_thickness=1000, p_dz=10000, p_increasing=1000, p_smooth=10, tol=1)[source]

Tries to optimize dzs

watershed_workflow.mesh.create_submesh(m2, shp)[source]

Given a shape that contains some cells of m2, create the submesh.

watershed_workflow.mesh.merge_meshes(meshes)[source]

Combines multiple 2D meshes into a single mesh.

It is assumed that the meshes to be combined have common vertices on the shared edge (no steiner points). labeledsets should be added after merging the mesh. Option of merging pre-existing labeledsets in the to-be-merged meshes will be added soon

Parameters:

meshes (list(mesh.Mesh2D)) – The list of meshes to be merged

Returns:

combined mesh

Return type:

mesh.Mesh2D

watershed_workflow.mesh.merge_mesh(mesh1, mesh2)[source]

merge two meshes (mesh.Mesh2D objects)

River-aligned Mesh

creates river mesh using quad, pentagon and hexagon elements

exception watershed_workflow.river_mesh.IntersectionError[source]
watershed_workflow.river_mesh.sort_children_by_angle(tree, reverse=False)[source]

Sorts the children of a given segment by their angle with respect to that segment.

watershed_workflow.river_mesh.create_rivers_meshes(rivers, widths=8, enforce_convexity=True, ax=None, label=True)[source]

Returns list of elems and river corridor polygons for a given list of river trees

Parameters:

rivers: list(watershed_workflow.river_tree.River object)

List of river tree along which river meshes are to be created

widths: 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”

enforce_convexity: boolean

If true, enforce convexity of the pentagons/hexagons at the junctions.

axmatplotlib Axes object, optional

For debugging – plots troublesome reaches as quad elements are generated to find tricky areas.

labelbool, optional = True

If true and ax is provided, animates the debugging plot with reach ID labels as the user hovers over the plot. Requires a widget backend for matplotlib.

returns:
  • elems (list(list)) – List of river elements, each element a list of indices into corr.coords.

  • corrs (list(shapely.geometry.Polygon)) – List of river corridor polygons, one per river, storing the coordinates used in elems.

watershed_workflow.river_mesh.create_river_mesh(river, widths=8, enforce_convexity=True, gid_shift=0, dilation_width=4, ax=None)[source]

Returns list of elems and river corridor polygons for a given river tree

Parameters:

river: watershed_workflow.river_tree.River object)

River tree along which mesh is to be created

widths: 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”

enforce_convexity: boolean

If true, enforce convexity of the pentagons/hexagons at the junctions.

gid_shift: int

All the node-ids used in the element defination are shifted by this number to make it consistant with the global ids in the m2 mesh, necessary in the case of multiple rivers

dilation_width: float

This is used for initial buffering of the river tree into the river corridor polygon. For a typical watershed the 4 m default should work well; for smaller domains, setting smaller initial dilation_width might be desirable (much smaller than expected quad element length).

axmatplotlib Axes object, optional

For debugging – plots troublesome reaches as quad elements are generated to find tricky areas.

returns:
  • elems (List(List)) – List of river elements, each element a list of indices into corr.coords.

  • corr (shapely.geometry.Polygon) – The polygon of all elements, stores the coordinates.

watershed_workflow.river_mesh.create_river_corridor(river, width)[source]

Returns a polygon representing the river corridor with a fixed dilation width.

Parameters:
  • river (watershed_workflow.river_tree.River object) – River tree along which corridor polygon is to be created.

  • width (float) – Width to dilate the river.

Returns:

Resulting polygon.

Return type:

shapely.geometry.Polygon

watershed_workflow.river_mesh.to_quads(river, corr, width, gid_shift=0, ax=None)[source]

Iterate over the rivers, creating quads and pentagons forming the corridor.

The global_id_adjustment is to keep track of node_id in elements w.r.t to global id in mesh mainly relevant for multiple river

Parameters:
  • rivers (watershed_workflow.river_tree.River) – river tree

  • corr (shapely.geometry.Polygon) – a river corridor polygon for the river

  • width (float) – fixed width used for dilating the river

  • gid_shift (int) – All the node-ids used in the element defination are shifted by this number to make it consistant with the global ids in the m2 mesh, necessary in the case of multiple rivers.

Returns:

elems – List of river elements

Return type:

list(list)

watershed_workflow.river_mesh.set_width(river, corr, widths=8, dilation_width=8, gid_shift=0)[source]

Adjust the river-corridor polygon width by order.

Parameters:
  • river (watershed_workflow.river_tree.River) – river tree along which mesh is to be created

  • corr (shapely.geometry.Polygon) – a river corridor polygon for the river

  • widths (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”

Returns:

river corridor polygon with adjusted width

Return type:

shapely.geometry.Polygon

watershed_workflow.river_mesh.move_to_target_separation(p1, p2, target, dilation_width=8)[source]

Returns the points after moving them to a target separation from each other

watershed_workflow.river_mesh.width_by_order(width_dict, order)[source]

Returns the reach width based using the {order:width dictionary}

watershed_workflow.river_mesh.convexity_enforcement(river, corr, gid_shift)[source]

Ensure convexity of each river-corridor element.

Moves nodes onto the convex hull of the element if needed.

Parameters:
Returns:

The new polygon with all convex elements.

Return type:

shapely.geometry.Polygon

watershed_workflow.river_mesh.make_convex_using_hull(points)[source]

Snaps non-convex points to the corresonding edge on the convex hull for the non-convex element

watershed_workflow.river_mesh.make_convex_by_nudge(points)[source]

Nudges the neck-points of the junction until the element is convex.

Used if efficient convexity does not work

watershed_workflow.river_mesh.hucsegs_at_intersection(point, hucs)[source]

For a given intersection point, return a list of indices for huc.segments touching this point

watershed_workflow.river_mesh.angle_rivers_segs(ref_seg, seg)[source]

Returns the angle of incoming-river-segment or huc-segment w.r.t outgoing river.

The angle is measured clockwise; this is useful to sort orientation wise and add river corridor points at junction

watershed_workflow.river_mesh.rc_points_for_rt_point(rt_point, node, river_corr)[source]

Returns the points (list of indices of coords) on the river-corridor-polygon for a given junction point on river tree.

watershed_workflow.river_mesh.adjust_seg_for_rc(seg, river_corr, new_seg_point, integrate_rc=False)[source]

Return the modified segment accomodating river-corridor-polygon (exclude river corridor or integrate with it).

watershed_workflow.river_mesh.adjust_hucs_for_river_corridors(hucs, rivers, river_corrs, integrate_rc=True, ax=None)[source]

Adjusts hucs to accomodate river corridor polygons.

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

  • rivers (list(watershed_workflow.river_tree.River)) – A list of river tree object

  • river_corrs (list(shapely.geometry.Polygons)) – A list of river corridor polygons for each river

  • integrate_rc (bool, optional) – if false, will leave gap in the huc whereever rc crosses huc except at the overall outlet; hence hucs.polygons() will break, this mode is to be used during triangulation to creates NodesEdges object with a rc as a whole if true, will extend the hucs-segments alogn the edge of quads

watershed_workflow.river_mesh.adjust_hucs_for_river_corridor(hucs, river, river_corr, gid_shift=0, integrate_rc=True, ax=None)[source]

Adjusts hucs to accomodate river corridor polygon.

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

  • river (watershed_workflow.river_tree.River object) – river tree

  • river_corr (shapely.geometry.Polygons) – A river corridor polygon for given river

  • integrate_rc (bool, optional) – If false, this will leave gap in the huc whereever rc crosses huc except at the overall outlet; hence hucs.polygons() will break, this mode is to be used during triangulation to creates NodesEdges object with a rc as a whole. If true, will extend the hucs-segments alogn the edge of quads.

Triangulation

Triangulates polygons

class watershed_workflow.triangulation.Nodes(decimals=3)[source]

A collection of nodes that are indexed in the order they are added.

Note this uses round() for an efficient solution, which is potentially fragile if numbers were not generated consistently. In this use case, however, it should be safe – numbers were originally rounded (in watershed_workflow.config), and then any function which inserted new points or moved points were always careful to ensure that all duplicates were always assigned the identical values – i.e. all math was done BEFORE assignment. So duplicates have identical representations in floating point. This suggests we shouldn’t have to round here at all, but we do anyway to keep things cleaner.

watershed_workflow.triangulation.oneway_trip_connect(inds)[source]

Connect indices in edges in a oneway fashion

watershed_workflow.triangulation.round_trip_connect(inds)[source]

Connect indices in edges in a round-trip fashion

class watershed_workflow.triangulation.NodesEdges(objlist=None)[source]

A collection of nodes and edges.

add(obj)[source]

Adds nodes and edges from obj into collection.

check(tol)[source]

Checks consistency of the interal representation.

watershed_workflow.triangulation.triangulate(hucs, river_corrs=None, internal_boundaries=None, hole_points=None, tol=1, **kwargs)[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()

  • river_corrs (list(shapely.geometry.Polygons), optional) – A list of river corridor polygons for each river

  • refinement_polygon (list, optional) – List of RiverTrees or other iterable collections of coordinates used to refine the mesh given the distance function refinement.

  • internal_boundaries (list, optional) – List of shapely objects or RiverTrees or other iterable collections of coordinates used as internal boundaries that must be included in the mesh.

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

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

  • meshpy.triangle.build() (Additional keyword arguments include all options for) –

watershed_workflow.triangulation.refine_from_max_area(max_area)[source]

Returns a refinement function based on max area, for use with Triangle.

watershed_workflow.triangulation.refine_from_river_distance(near_distance, near_area, away_distance, away_area, rivers)[source]

Returns a graded refinement function based upon a distance function from rivers, for use with Triangle.

Triangle area must be smaller than near_area when the triangle centroid is within near_distance from the river network. Area must be smaller than away_area when the triangle centroid is at least away_distance from the river network. Area must be smaller than a linear interpolant between near_area and away_area when between near_distance and away_distance from the river network.

watershed_workflow.triangulation.refine_from_max_edge_length(edge_length)[source]

Returns a refinement function based on max edge length, for use with Triangle.

watershed_workflow.triangulation.pick_hole_point(poly)[source]

A function to pick a point inside a polygon

Condition

watershed_workflow.condition.fill_pits(m2, outlet=None, algorithm=3)[source]

Conditions a mesh, IN PLACE, by removing pits.

Starts at outlet and works through all coordinates in the mesh, ensuring that there is a path to all nodes of the mesh from the outlet that monotonically increases in elevation.

Available algorithms (likely all should be equivalent):

1: original, 2-pass algorithm 2: refactored single-pass algorithm based on sorted lists 3: boundary marching method. Should be the fastest.

Parameters:
  • m2 (mesh.Mesh2D object) – The mesh to condition.

  • outlet (int, optional) – If provided, the ID of the point to start conditioning from. If not provided, will use the boundary node with minimal elevation.

  • algorithm (int) – See above, defaults to 3.

watershed_workflow.condition.fill_pits1(points, outletID=None)[source]

This is the origional, 2-pass algorithm, and is likely inefficient.

watershed_workflow.condition.fill_pits2(points, outletID)[source]

This is a refactored, single pass algorithm that leverages a sorted list.

watershed_workflow.condition.fill_pits3(points, outletID)[source]

This algorithm is based on a boundary marching method

watershed_workflow.condition.fill_pits_dual(m2, is_waterbody=None, outlet_edge=None, eps=0.001)[source]

Conditions the dual of the mesh, IN PLACE, by filling pits.

This ensures the property that, starting with an outlet cell, there is a path to every cell by way of faces that is monotonically increasing in elevation (except in cells which are a part of waterbodies).

Parameters:
  • m2 (mesh.Mesh2D object) – The mesh to condition.

  • is_waterbody (np.ndarray(int, shape=ncells), optional) – A boolean/integer mask of length m2.num_cells, with True values indicating that the given cell is a part of a waterbody that may be a depression. This is important for reservoirs/lakes/etc where bathymetry is known and pits are physical.

  • outlet_edge ((int,int), optional) – If provided, the point to start conditioning from. If not provided, will use the edge on the boundary of m2 with the lowest elevation.

  • eps (float, optional=1e-3) – A small vertical displacement for soft tolerances.

watershed_workflow.condition.identify_local_minima(m2)[source]

For all cells, identify if their centroid elevation is lower than the elevation of all neighbors.

Parameters:

m2 (mesh.Mesh2D object) – The mesh to check.

Returns:

Array of 0s and 1s, where 1 indicates a local minima.

Return type:

np.array

watershed_workflow.condition.smooth(img_in, algorithm='gaussian', **kwargs)[source]

Smooths an image according to an algorithm, passing kwargs on to that algorithm.

watershed_workflow.condition.condition_river_meshes(m2, rivers, smooth=False, use_parent=False, lower=False, use_nhd_elev=False, treat_banks=False, depress_upstream_by=None, network_burn_in_depth=None, ignore_in_sweep=None)[source]

For multile rivers, condition, IN PLACE, the elevations of stream-corridor elements to ensure connectivity throgh culverts, skips ponds, maintain monotonicity, or otherwise enforce depths of constructed channels.

watershed_workflow.condition.condition_river_mesh(m2, river, smooth=False, use_parent=False, lower=False, use_nhd_elev=False, treat_banks=False, depress_upstream_by=None, network_burn_in_depth=None, ignore_in_sweep=None)[source]

Condition, IN PLACE, the elevations of stream-corridor elements to ensure connectivity throgh culverts, skips ponds, maintain monotonicity, or otherwise enforce depths of constructed channels.

m2: watershed_workflow.mesh.Mesh2D object

2D mesh with 3D coordinates.

river: watershed_workflow.river_tree.River object

River tree with node.elements added for quads

smooth: boolean, optional

If true, smooth the profile of each reach using a gaussian filter (mainly to pass through railroads and avoid reservoirs).

use_parent: boolean, optional

If true, use the segment from the original parent reach while smoothing (seems to be not making a huge difference).

lower: boolean, optional

If true, lower the smoothed bed profile to match the lower points on the raw bed profile. This is useful particularly for narrow ag. ditches where NHDPLus flowlines often do not coincide with the DEM depressions and so stream-elements intermitently fall into them.

use_nhd_elev: boolean, optional

If true, enforce maximum and minimum elevation for each reach provided in NHDPlus.

treat_banks: boolean, optional

Where the river is passing right next to the reservoir or NHDline is misplaced into the reservoir, banks may fall into the reservoir. If true, this will enforce that the bank node is at a higher elevation than the stream bed elevation.

depress_upstream_by: float, optional

If the depression is not captured well in the DEM, the river-mesh elements (streambed) headwater reaches may be lowered by this number. The effect of propogated downstream only upto where it is needed to maintain topographic gradients on the network scale in the network sweep step.

network_burn_in_depth: float, dict, or function

Like depress_upstream_by, this also lowers river-mesh elements by this value, but this variant lowers all reaches. The depth may be provided as a float (uniform lowering), dictionary {stream order : depth to depress by}, or as a function of drainage area.

ignore_in_sweep: list, optional

If provided, a list of IDs to not be burned in via the network sweep.

watershed_workflow.condition.elevate_rivers(rivers, crs, dem, dem_profile)[source]

Elevate the river using dem and store reach-bed-profile as node properties.

Regions

This namespace hosts functions to add labelsets to define regions

watershed_workflow.regions.add_nlcd_labeled_sets(m2, nlcd_colors, nlcd_names=None)[source]

Add labeled sets to a mesh – one per unique color.

Parameters:
  • m2 (mesh.Mesh2D object) – The mesh to label.

  • nlcd_colors (np.array) – 1D array of len m2.num_cells of cell colors.

  • nlcd_names (dict, optional) – Dictionary mapping colors to color name.

watershed_workflow.regions.add_polygonal_regions(m2, polygons, labels=None, kind='watershed', volume=False)[source]

Label m2 with region(s) for each polygon.

Always adds a surface region; if volume also adds a volume region for extrusion in 3D.

Parameters:
  • m2 (mesh.Mesh2D) – The mesh to label.

  • polygons (iterable[shapely.Polygon]) – The polygons covering each watershed.

  • labels (iterable[str], optional) – Names of the polygons (often HUC strings)

  • kind (str, optional) – The kind of polygon to add – this is only used to generate labels. Default is ‘watershed’.

  • volume (bool, optional) – If true, also add the volumetric region below the polygon that will be extruded in a 3D mesh eventually.

Returns:

partitions – A list of length polygons, each entry of which is the list of cell indices in that polygon.

Return type:

list[list[int]]

watershed_workflow.regions.add_watershed_regions(m2, polygons, labels=None)[source]

Deprecated – kept for backward compatibility.

watershed_workflow.regions.add_watershed_regions_and_outlets(m2, hucs, outlets=None, outlet_width=300, labels=None, exterior_outlet=True)[source]

Add four labeled sets to m2 for each polygon:

  • cells in the polygon, to be extruded

  • cells in the polygon, to be kept as faces upon extrusion

  • boundary of the polygon (edges, kept as faces)

  • outlet of the polygon (edges, kept as faces, within outlet width of the outlet)

Parameters:
  • m2 (mesh.Mesh2D) – The mesh to label.

  • hucs (iterable[shapely.Polygon] or SplitHUCs object) – Watershed polygons.

  • outlets (iterable[shapely.Point], optional) – If provided, the outlet points. If SplitHUCs are provided, outlets are in that object and this must be None.

  • outlet_width (float, optional) – How wide should the outlet region be? Note this should include not just the river width, but more realistically the best resolved floodplain, or 1-2 face-widths, whichever is bigger.

  • labels (iterable[str], optional) – Name of the polygons. If SplitHUCs are provided, HUC names are used as the default.

  • exterior_outlet (bool, optional) – If true, find the outlet point that intersects the boundary and include regions around that outlet as well.

watershed_workflow.regions.add_river_corridor_regions(m2, rivers, labels=None)[source]

Add labeled sets to m2 for each river corridor.

Parameters:

m2: watershed_workflow.mesh.Mesh2D object

2D mesh elevated on DEMs

rivers: list(watershed_workflow.river_tree.RiverTree)

List of rivers used to create the river corridors.

labels: list(str), optional

List of names, one per river.

watershed_workflow.regions.add_regions_by_stream_order_rivers(m2, rivers, labels=None)[source]

Add labeled sets to m2 for reaches of each stream order for each river.

Parameters:

m2: watershed_workflow.mesh.Mesh2D object

2D mesh elevated on DEMs

rivers: list(watershed_workflow.river_tree.RiverTree)

List of rivers used to create the river corridors.

labels: list(str), optional

List of names, one per river.

watershed_workflow.regions.add_regions_by_stream_order(m2, river, river_id=0)[source]

Add labeled sets to m2 for reaches of each stream order .

Parameters:

m2: watershed_workflow.mesh.Mesh2D object

2D mesh elevated on DEMs

rivers: list(watershed_workflow.river_tree.RiverTree)

List of rivers used to create the river corridors.

river_id: str/int, optional

river identifier/name.

watershed_workflow.regions.add_region_by_reach_id(m2, river, reach_ids=None, labels=None)[source]

Add labeled sets to m2 for reaches of each stream order .

Parameters:

m2: watershed_workflow.mesh.Mesh2D object

2D mesh elevated on DEMs

river: watershed_workflow.river_tree.RiverTree

List of rivers used to create the river corridors.

reaches: list(str)

list of NHDID IDs to be labeled.

watershed_workflow.regions.getNode(nhd_id, rivers)[source]

return node given NHDID