watershed module

Functions that implement raster watershed analyses

Function

Description

condition

Conditions a DEM by filling pits and depressions, and resolving flats

flow

Computes D8 flow directions from a conditioned DEM

slopes

Computes D8 flow slopes

relief

Computes vertical relief to the nearest ridge cell

accumulation

Computes basic, weighted, or masked flow accumulation

catchment

Returns the catchment mask for a watershed pixel

network

Returns the stream segments as a list of shapely.LineString objects

Georeferencing

The slopes function requires the input DEM to have both a CRS and affine Transform. Most workflows will require flow slopes, so we recommend using a properly georeferenced DEM whenever possible.

Tip

Most users will not need the catchment or network functions, as these are implemented internally by the Segments class.

Warning

The functions in this module rely on the pysheds library, which assigns a default NoData of 0 to any raster lacking a NoData value. This can cause unexpected results when a raster has valid 0 values and lacks NoData. Consider using a placeholder NoData when this is the case.

Workflow

Typical workflow is to first use the condition function to condition a DEM (i.e. filling pits and resolving flats). Then, use the flow function to compute D8 flow directions from a DEM. These flow directions are an essential input to all other watershed functions, as well as the Segments class. With the flow directions, users can compute flow accumulation (also referred to as upslope area), D8 flow slopes, and the vertical relief of watershed pixels.

Flow direction style

This module produces flow directions in the TauDEM style:

\[\begin{split}\begin{matrix} 4 & 3 & 2\\ 5 & \mathrm{X} & 1\\ 6 & 7 & 8\\ \end{matrix}\end{split}\]

where X is the current pixel, and integers indicate flow in a particular direction. So for example, if pixel X flows into the next pixel to the left, then X will be marked with a flow direction of 5. But if X flows into the pixel to the right, then its flow direction will be 1.

Note

All pfdf routines that use flow directions expect values in the TauDEM style.


pfdf.watershed.condition(dem, *, fill_pits=True, fill_depressions=True, resolve_flats=True)

Conditions a DEM to resolve pits, depressions, and/or flats

Conditioning a DEM
condition(dem)

Conditions a DEM by filling pits, filling depressions, and then resolving flats. A pit is defined as a single cell lower than all its surrounding neighbors. When a pit is filled, its elevation is raised to match that of its lowest neighbor. A depression consists of multiple cells surrounded by higher terrain. When a depression is filled, the elevations of all depressed cells are raised to match the elevation of the lowest pixel on the border of the depression. Flats are sets of adjacent cells with the same elevation, and often result from filling pits and depressions (although they may also occur naturally). When a flat is resolved the elevations of the associated cells are minutely adjusted so that their elevations no longer match.

Skipping steps
condition(dem, *, fill_pits=False)
condition(dem, *, fill_depressions=False)
condition(dem, *, resolve_flats=False)

Allows you to skip specific steps of the conditioning algorithm. Setting an option to False will disable the associated conditioning step. Raises a ValueError if you attempt to skip all three steps.

Inputs:
  • dem (Raster-like) – A digital elevation model raster

  • fill_pits (bool) – True (default) to fill pits. False to disable this step

  • fill_depressions (bool) – True (default) to fill depressions. False to disable this step

  • resolve_flats (bool) – True (default) to resolve flats. False to disable this step

Outputs:

Raster – A conditioned DEM raster

pfdf.watershed.flow(dem)

Compute D8 flow directions from a conditioned DEM

flow(dem)

Computes D8 flow directions for a conditioned DEM. Flow direction numbers follow the TauDEM style. Values of 0 indicate NoData - these may result from NoData values in the original DEM, as well as any unresolved pits, depressions, or flats.

Inputs:
  • dem (Raster-like) – A conditioned digital elevation model raster

Outputs:

Raster – The D8 flow directions for the DEM

pfdf.watershed.slopes(dem, flow, dem_per_m=1, check_flow=True)

Computes D8 flow slopes for a watershed

Computing slopes
slopes(dem, flow)
slopes(dem, flow, dem_per_m)

Returns D8 flow slopes for a watershed. Computes these slopes using a DEM raster, and TauDEM-style D8 flow directions. The DEM must have both a CRS and an affine Transform. Note that the DEM may be a raw DEM - it does not need to resolve pits and flats, although you may wish to use a resolved DEM for consistency across your analysis. The returned slopes will be unitless gradients. The rise is determined using the DEM, and the run is determined from the CRS and transform. If the CRS base unit is not meters, converts the run to meters before computing slope gradients.

By default, this routine assumes that the DEM is in units of meters. If this is not the case, use the “dem_per_m” to specify a conversion factor (number of DEM units per meter).

Disable flow validation
slopes(..., check_flow=False)

Disables validation checking of the flow directions raster. Validation is not necessary for flow directions directly output by the “watershed.flow” function, and disabling the validation can improve runtimes for large rasters.

Warning

This option may produce unexpected results if the flow directions raster contains invalid values.

Inputs:
  • dem (Raster-like) – A digital elevation model raster

  • flow (Raster-like) – A raster with TauDEM-style D8 flow directions

  • dem_per_m (scalar) – A conversion factor from DEM units to meters

  • check_flow (bool) – True (default) to validate the flow directions raster. False to disable validation checks.

Outputs:

Raster – The computed D8 flow slopes for the watershed

pfdf.watershed.relief(dem, flow, check_flow=True)

Computes vertical relief to the highest ridge cell

Computing Relief
relief(dem, flow)

Computes the vertical relief for each watershed pixel. Here, vertical relief is defined as the change in elevation between each pixel and its nearest ridge cell. (A ridge cell is an upslope cell with no contributing flow from other pixels). Computes these slopes using a DEM raster, and TauDEM-style D8 flow directions.

Note

The DEM can be a raw DEM (as opposed to a conditioned DEM). It does not need to resolve pits and flats.

Disable flow validation
relief(..., check_flow=False)

Disables validation checking of the flow directions raster. Validation is not necessary for flow directions directly output by the flow function, and disabling the validation can improve runtimes for large rasters.

Warning

This option may produce unexpected results if the flow directions raster contains invalid values.

Inputs:
  • dem (Raster-like) – A digital elevation model raster

  • flow (Raster-like) – A TauDEM-style D8 flow direction raster

  • check_flow (bool) – True (default) to validate the flow directions raster. False to disable validation checks.

Outputs:

Raster – The vertical relief of the nearest ridge cell.

pfdf.watershed.accumulation(flow, weights=None, mask=None, *, times=None, omitnan=False, check_flow=True)

Computes basic, weighted, or masked flow accumulation

Computing Accumulation
accumulation(flow)

Uses D8 flow directions to compute basic flow accumulation. In this setup, each pixel is given a value of 1, so the accumulation for each pixel indicates the number of upslope pixels. Note that each pixel is included in its own accumulation, so the minimum valid accumulation is 1. NoData values are indicated by NaN. Flow directions should follow the TauDEM style.

Weighted Accumulation
accumulation(flow, weights)
accumulation(flow, weights, *, omitnan=True)

Computes weighted accumulations. Here, the value of each pixel is set by the input “weights” raster, so the accumulation for each pixel is a sum over itself and all upslope pixels. The weights raster must have the same shape, transform, and crs as the flow raster.

In the default case, NaN and NoData values in the weights raster are set to propagate through the accumulation. So any pixel that is downstream of a NaN or a NoData weight will have its accumulation set to NaN. Setting omitnan=True will change this behavior to instead ignore NaN and NoData values. Effectively, NaN and NoData pixels will be given weights of 0.

Masking
accumulation(..., mask)

Computes a masked accumulation. In this syntax, only the True elements of the mask are included in accumulations. All False elements are given a weight of 0. NoData elements in the mask are interpreted as False. The accumulation for each pixel is thus the sum over all catchment pixels included in the mask. If weights are not specified, then all included pixels are given a weight of 1. Note that the mask raster must have the same shape, transform, and crs as the flow raster.

Multiplicative Factor
accumulation(..., *, times)

Returns accumulation multiplied by the indicated scalar value. This option is often set to the area of a raster pixel in order to return accumulation in specific units, rather than pixel counts.

Disable flow validation
accumulation(..., *, check_flow=False)

Disables validation checking of the flow directions raster. Validation is not necessary for flow directions directly output by the watershed.flow function, and disabling the validation can improve runtimes for large rasters.

Warning

This option may produce unexpected results if the flow directions raster contains invalid values.

Inputs:
  • flow (Raster-like) – A D8 flow direction raster in the TauDEM style

  • weights (Raster-like) – A raster indicating the value of each pixel

  • mask (Raster-like) – A raster whose True elements indicate pixels that should be included in the accumulation.

  • times (scalar) – A multiplicative constant applied to the computed accumulation.

  • omitnan (bool) – True to ignore NaN and NoData values in the weights raster. False (default) propagates these values as NaN to all downstream pixels.

  • check_flow (bool) – True (default) to validate the flow directions raster. False to disable validation checks.

Outputs:

Raster – The computed flow accumulation

pfdf.watershed.catchment(flow, row, column, check_flow=True)

Returns the catchment mask for a DEM pixel

Locate a catchment
catchment(flow, row, column)

Determines the extent of the catchment upstream of the DEM pixel at the indicated row and column. Returns a mask for this catchment extent. The mask will have the same shape as the input flow directions raster - True values indicate pixels that are in the upstream catchment extent, False values are outside of the catchment. Any NoData values in the flow directions will become False values in the catchment mask.

Disable flow validation
catchment(..., check_flow=False)

Disables validation checking of the flow directions raster. Validation is not necessary for flow directions directly output by the watershed.flow function, and disabling the validation can improve runtimes for large rasters.

Warning

This option may produce unexpected results if the flow directions raster contains invalid values.

Inputs:
  • flow (Raster-like) – D8 flow directions for the DEM (in the TauDEM style)

  • row (int) – The row index of the queried pixel in the DEM

  • column (int) – The column index of the queried pixel in the DEM

  • check_flow (bool) – True (default) to validate the flow directions raster. False to disable validation checks.

Outputs:

Raster – The upstream catchment mask for the queried pixel

pfdf.watershed.network(flow, mask, max_length=None, units='meters', check_flow=True)

Returns a list of stream segment LineStrings

Delineate a network
network(flow, mask)

Calculates a stream segment network and returns the segments as a list of shapely.LineString’’ objects. These stream segments approximate the river beds in a drainage basin - they are not the full catchment basin.

The stream segment network is determined using a TauDEM-style D8 flow direction raster and a raster mask. The mask is used to indicate the pixels under consideration as stream segments. True pixels may possibly be assigned to a stream segment, False pixels will never be assigned to a stream segment. The mask typically screens out pixels with low flow accumulations, and may include other screenings - for example, to remove pixels in large bodies of water, or pixels below developed areas.

Maximum length
network(..., max_length)
network(..., max_length, units)

Also specifies a maximum length for the segments in the network. Any segment longer than this length will be split into multiple pieces. The split pieces will all have the same length, which will be <= max_length. By default, the command interprets the units of max_length as meters. Use the units option to specify max_length in different units instead. Unit options include: “base” (CRS/Transform base unit), “meters” (default), “kilometers”, “feet”, and “miles”. Note that the meters/kilometers/feet/miles options all require the flow raster to have a CRS. The units="base" option relaxes this requirement.

Disable flow validation
network(..., check_flow=False)

Disables validation checking of the flow directions raster. Validation is not necessary for flow directions directly output by the watershed.flow function, and disabling the validation can improve runtimes for large rasters.

Warning

This option may produce unexpected results if the flow directions raster contains invalid values.

Inputs:
  • flow (Raster-like) – A TauDEM-style D8 flow direction raster

  • mask (Raster-like) – A raster whose True values indicate the pixels that may potentially belong to a stream segment.

  • max_length (scalar) – A maximum allowed length for segments in the network. Units should be the same as the base units of the flow raster CRS

  • units (str) – Indicates the units of the max_length input. Options include: “base” (CRS/Transform base unit), “meters” (default), “kilometers”, “feet”, and “meters”.

  • check_flow (bool) – True (default) to validate the flow directions raster. False to disable validation checks.

Outputs:

list[shapely.LineString] – The stream segments in the network, represented by shapely.LineString objects. The coordinates of each LineString proceed from upstream to downstream. Coordinates are relative to the flow raster CRS (rather than raster pixel indices).