Function API Reference

accumulate_flow(d8_fdr: Union[xarray.DataArray, str, Path], engine: SupportsAccumulateFlow = 'pysheds', upstream_pour_points: Optional[PourPointValuesDict] = None, weights: Optional[xarray.DataArray] = None, out_path: Optional[Union[str, Path]] = None, **kwargs) xarray.DataArray[source]

Create a Flow Accumulation Cell (FAC) raster from a D8 Flow Direction Raster.

Parameters
  • d8_fdr – A D8 Flow Direction Raster (dtype=Int).

  • engine – A terrain engine class that supports flow accumulation.

  • upstream_pour_points – A list of lists each with with coordinate tuples as the first item [0], and updated cell values as the second [1]. This allows the FAC to be made with boundary conditions such as upstream basin pour points.

  • weights – A grid defining the value to accumulate from each cell. Default is a grid of 1s.

  • out_path – Defines a path to save the output raster.

  • **kwargs – keyword arguments, specific options depend on the engine being used.

Returns

The output Flow Accumulation Cells (FAC) raster.

accumulate_parameter(d8_fdr: Union[xarray.DataArray, str, Path], parameter_raster: Union[xarray.DataArray, str, Path], engine: SupportsAccumulateParameter = 'pysheds', upstream_pour_points: Optional[PourPointValuesDict] = None, out_path: Optional[Union[str, Path]] = None, **kwargs) xarray.DataArray[source]

Create a parameter accumulation raster from a D8 Flow Direction Raster and a parameter raster.

A key aspect of this function is that the output DataArray will have dimensions matching param:parameter_raster.

Parameters
  • d8_fdr – A D8 Flow Direction Raster (dtype=Int).

  • parameter_raster – A parameter raster aligned via tools.align_raster() with the us_fdr. This can be multi-dimensional (i.e. f(x, y, t)), and if so, a multi-dimensional output is returned.

  • engine – A terrain engine class that supports parameter accumulation.

  • upstream_pour_points – A list of lists each with with coordinate tuples as the first item [0], and updated cell values as the second [1]. This allows the FAC to be made with boundary conditions such as upstream basin pour points.

  • out_path – Defines a path to save the output raster.

  • **kwargs – keyword arguments, specific options depend on the engine being used.

Returns

The output parameter accumulation raster.

adjust_parameter_raster(parameter_raster: Union[xarray.DataArray, str, Path], d8_fdr: Union[xarray.DataArray, str, Path], upstream_pour_points: PourPointValuesDict, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Adds values to a parameter raster’s at specific pour points / coordinates for use in accumulate_parameter().

The main utility of this function is to enable cascading accumulation values from one basin or raster to another.

Parameters
  • parameter_raster – Input parameter raster to update.

  • d8_fdr – a D8 Flow Direction Raster.

  • upstream_pour_points – A dictionary with coordinates to update values at, and the values to add.

  • out_path – Defines a path to save the output raster.

Returns

The updated parameter raster.

align_raster(in_raster: Union[xarray.DataArray, str, Path], match_raster: Union[xarray.DataArray, str, Path], resample_method: str = 'nearest', out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Aligns the projection/CRS, spatial extent, and resolution of one raster to another.

Parameters
  • in_raster – Input raster that needs to be aligned.

  • match_raster – Raster to align to.

  • resample_method – A valid resampling method from rasterio.enums.Resample (default=’nearest’). NOTE: Do not use any averaging resample methods when working with a categorical raster!

  • out_path – Defines a path to save the output raster.

Returns

The output aligned raster.

binarize_categorical_raster(cat_raster: Union[xarray.DataArray, str, Path], categories_dict: Optional[Dict[int, str]] = None, ignore_categories: Optional[List] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Converts a single band categorical raster into a binary multi-band binary raster.

Each output band represent a unique category, where 1 encodes cells belonging the the category, and 0 encodes cells belonging to any other category. This function is used to prep a categorical raster for parameter_accumulate().

Parameters
  • cat_raster – A categorical (dtype=int) raster with N unique categories (i.e., land cover classes).

  • categories_dict – A dictionary assigning string names (values) to integer raster values (keys).

  • ignore_categories – Category cell values not include in the output raster.

  • out_path – Defines a path to save the output raster.

Returns

A N-band multi-dimensional raster as a xarray DataArray object.

binarize_nodata(in_raster: Union[xarray.DataArray, str, Path], nodata_value: Optional[Union[int, float]] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Creates an output binary raster based on an input where nodata values -> 1, and valued cells -> 0.

Note: while param:inverse=True this can be used with pyfunc:apply_mask() to match nodata cells between rasters.

Parameters
  • in_raster – Input raster.

  • nodata_value – If the nodata value for param:in_raster is not in the metadata, set this parameter to equal the cell value storing nodata (i.e., np.nan or -999).

  • out_path – Defines a path to save the output raster.

Returns

The output binary mask raster.

check_function_kwargs(function: callable, engine: str) Dict[str, Union[str, int, float]][source]

Provides a dictionary of allowed kwargs keywords + input types for a function.

NOTE: This function will raise a ValueError if a non-terrain_engine

function is provided as the input.

Parameters
  • function – A function belonging to a terrain_engine class.

  • engine – The name of the engine being used.

Returns

A dictionary with allowed kwargs as keys, and allowed input types as values.

clip(in_raster: Union[xarray.DataArray, str, Path], match_raster: Optional[Union[xarray.DataArray, str, Path]] = None, match_shapefile: Optional[Union[geopandas.GeoDataFrame, str, Path]] = None, custom_bbox: Optional[List] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Clips a raster to the rectangular extent (aka bounding box) of another raster (or shapefile).

Parameters
  • in_raster – Input raster.

  • match_raster – If defined, in_raster is clipped to match the extent of match_raster.

  • match_shapefile – A shapefile that is used to define the output extent if match_raster == None.

  • out_path – Defines a path to save the output raster.

  • custom_bbox – A list with bounding box coordinates that define the output extent if match_raster == None. Note: Coordinates must be of the form [minX, minY, maxX, maxY]. Note: Using this parameter assumes that coordinates match the CRS of param:in_raster.

Returns

The input raster clipped by the desired geometry.

convert_fdr_formats(d8_fdr: Union[xarray.DataArray, str, Path], out_format: str, in_format: Optional[str] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Converts the D8 encoding of Flow Direction Rasters (FDR).

Parameters
  • d8_fdr – The input D8 Flow Direction Raster (FDR).

  • out_format – A valid D8 flow direction format name in custom_types.D8ConversionDicts.keys().

  • in_format – A valid D8 flow direction format name in custom_types.D8ConversionDicts.keys() that overrides the auto-recognized format from param:d8_fdr. Note: manually inputting param:in_format will improve performance.

  • out_path – Defines a path to save the output raster.

Returns

The re-encoded D8 Flow Direction Raster (FDR).

d8_to_dinfinity(d8_fdr: Union[xarray.DataArray, str, Path], out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Converts a D8 Flow Direction Raster to D-Infinity.

Parameters
  • d8_fdr – A D8 Flow Direction Raster (dtype=int).

  • out_path – Defines a path to save the output raster.

Returns

The output D-Inf Flow Direction Raster.

decay_accumulation(d8_fdr: Union[xarray.DataArray, str, Path], decay_raster: Union[xarray.DataArray, str, Path], engine: SupportsDecayAccumulation = 'taudem', upstream_pour_points: Optional[PourPointValuesDict] = None, parameter_raster: Optional[Union[xarray.DataArray, str, Path]] = None, out_path: Optional[Union[str, Path]] = None, **kwargs) xarray.DataArray[source]

Creates a D-Infinity based accumulation raster (parameter or cell accumulation) while applying decay via a multiplier_raster.

NOTE: Replaces tools.decayAccum() from V1 FCPGtools.

Parameters
  • dinf_fdr – A flow direction raster in D-Infinity format. This input can be made with tools.d8_to_dinfinity().

  • decay_raster – A decay ‘multiplier’ raster calculated from distance to stream via tools.make_decay_raster().

  • engine – A terrain engine class that supports decayed accumulation.

  • upstream_pour_points – A list of lists each with with coordinate tuples as the first item [0], and updated cell values as the second [1]. This allows the FAC to be made with boundary conditions such as upstream basin pour points.

  • parameter_raster – A parameter raster aligned via tools.align_raster() with the us_fdr. This can be multi-dimensional (i.e. f(x, y, t)), and if so, a multi-dimensional output is returned.

  • out_path – Defines a path to save the output raster.

  • **kwargs – keyword arguments, specific options depend on the engine being used.

Returns

The output decayed accumulation raster.

distance_to_stream(d8_fdr: Union[xarray.DataArray, str, Path], fac_raster: Union[xarray.DataArray, str, Path], accum_threshold: int, engine: SupportsDistanceToStream = 'taudem', out_path: Optional[Union[str, Path]] = None, **kwargs) xarray.DataArray[source]

Calculates distance each cell is from a stream (as defined by a cell accumulation threshold).

NOTE: Replaces tools.dist2stream() from V1 FCPGtools.

Parameters
  • d8_fdr – A D8 Flow Direction Raster (dtype=Int).

  • fac_raster – A Flow Accumulation Cell (FAC) raster output from accumulate_flow().

  • accum_threshold – The # of upstream/accumulated cells to consider a cell a stream.

  • engine – A terrain engine class that supports calculating distance to stream.

  • out_path – Defines a path to save the output raster.

  • **kwargs – keyword arguments, specific options depend on the engine being used.

Returns

A raster with values of D8 flow distance from each cell to the nearest stream.

extreme_upslope_values(d8_fdr: Union[xarray.DataArray, str, Path], parameter_raster: Union[xarray.DataArray, str, Path], engine: SupportsExtremeUpslopeValues = 'taudem', mask_streams: Optional[Union[xarray.DataArray, str, Path]] = None, out_path: Optional[Union[str, Path]] = None, get_min_upslope: bool = False, **kwargs) xarray.DataArray[source]

Finds the max (or min if get_min_upslope=True) value of a parameter grid upstream from each cell in a D8 FDR raster.

NOTE: Replaces tools.ExtremeUpslopeValue() from V1 FCPGtools.

Parameters
  • d8_fdr – A flow direction raster .

  • parameter_raster – A parameter raster to find the max values from.

  • engine – A terrain engine class that supports finding extreme upslope values.

  • mask_streams – A stream mask raster from tools.mask_streams(). If provided, the output will be masked to only stream cells.

  • out_path – Defines a path to save the output raster.

  • get_min_upslope – If True, the minimum upslope value is assigned to each cell.

  • **kwargs – keyword arguments, specific options depend on the engine being used.

Returns

A raster with max (or min) upstream value of the parameter grid as each cell’s value.

find_basin_pour_points(fac_raster: Union[xarray.DataArray, str, Path], basins_shp: str, basin_id_field: str = 'HUC12', use_huc4: bool = True) PourPointLocationsDict[source]

Find pour points (aka outflow cells) in a FAC raster by basin using a shapefile.

Parameters
  • fac_raster – A Flow Accumulation Cell raster (FAC).

  • basins_shp – A .shp shapefile containing basin geometries.

  • basin_id_field – Default behavior is for each GeoDataFrame row to be a unique basin. However, if one wants to use a higher level basin id that is shared across rows, this should be set to the column header storing the higher level basin id.

  • use_huc4 – Either ‘HUC4’ or ‘HUC12’.

Returns

A dictionary with keys (i.e., basin IDs) storing coordinates as a tuple(x, y).

find_fac_pour_point(fac_raster: Union[xarray.DataArray, str, Path], basin_name: Optional[Union[str, int]] = None) PourPointLocationsDict[source]

Find pour points (aka outflow cells) in a FAC raster by basin using a shapefile.

Parameters
  • fac_raster – A Flow Accumulation Cell raster (FAC).

  • basin_name – Allows a name to be given to the FAC.

Returns

A dictionary with keys (i.e., basin IDs) storing coordinates as a tuple(x, y).

get_pour_point_values(pour_points_dict: PourPointLocationsDict, accumulation_raster: Union[xarray.DataArray, str, Path]) PourPointValuesDict[source]

Get the accumulation raster values from downstream pour points.

Note: This function is intended to feed into accumulate_flow() or parameter_accumulate() param:upstream_pour_points.

Parameters
  • pour_points_dict – A dictionary of form fcpgtools.custom_types.PourPointValuesDict.

  • accumulation_raster – A Flow Accumulation Cell raster (FAC) or a parameter accumulation raster.

Returns

A list of tuples (one for each pour point) storing their coordinates [0] and accumulation value [1].

load_raster(in_raster: Union[xarray.DataArray, str, Path]) xarray.DataArray[source]

Loads a raster into a xarray.DataArray object.

load_shapefile(in_shapefile: Union[geopandas.GeoDataFrame, str, Path]) geopandas.GeoDataFrame[source]

Loads a shapefile into a geopandas.GeoDataFrame

make_decay_raster(distance_to_stream_raster: Union[xarray.DataArray, str, Path], decay_factor: Union[int, float], out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Creates a decay raster based on distance to stream for use in decay_accumulate().

Grid cell values are computed as the inverse number of grid cells, np.exp((-1 * distance_to_stream * cell_size) / (cell_size ** k)), where k is a constant applied to the cell size values.

Parameters
  • distance_to_stream_raster – A distance to stream raster output from distance_to_stream().

  • decay_factor – Dimensionless constant applied to decay factor denominator. NOTE: Set k to 2 for ‘moderate’ decay; greater than 2 for slower decay; or less than 2 for faster decay.

  • out_path – Defines a path to save the output raster.

Returns

The output decay raster for use in decay_accumulate().

make_fac_weights(parameter_raster: Union[xarray.DataArray, str, Path], fdr_raster: Union[xarray.DataArray, str, Path], out_of_bounds_value: Union[float, int], out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Preps param:parameter_raster for parameter accumulation by matching the nodata boundary to param:fdr_raster.

NOTE: Only this function AFTER aligning the parameter raster to the FDR raster via tools.align_raster()!

Parameters
  • parameter_raster – A parameter raster.

  • fdr_raster – A Flow Direction Raster (either D8 or D-inf).

  • out_of_bounds_value – The value to give param:parameter_raster cells outside the data boundary of param:fdr_raster. Note that this is automatically set within terrain_engine functions.

  • out_path – defines a path to save the output raster.

Returns

The prepped parameter grid.

make_fcpg(param_accum_raster: Union[xarray.DataArray, str, Path], fac_raster: Union[xarray.DataArray, str, Path], out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Creates a Flow Conditioned Parameter Grid raster by dividing a parameter accumulation raster by a Flow Accumulation Cell (FAC) raster.

Parameters
  • param_accum_raster – (xr.DataArray or str raster path)

  • fac_raster – (xr.DataArray or str raster path) input FAC raster.

  • out_path – (str or pathlib.Path, default=None) defines a path to save the output raster.

Returns

The output FCPG raster as a xarray DataArray object.

mask_streams(fac_raster: Union[xarray.DataArray, str, Path], accumulation_threshold: Union[int, float], out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

A simplified version of tools.value_mask() that outputs np.nan for non-‘stream’ cells below the accumulation threshold.

Parameters
  • fac_raster – A single band flow accumulation (FAC) raster.

  • accumulation_threshold – The flow accumulation threshold.

Returns

The output raster with cells below the threshold encoded as np.nan

reproject_raster(in_raster: Union[xarray.DataArray, str, Path], out_crs: Optional[Union[xarray.DataArray, str, Path, geopandas.GeoDataFrame]] = None, resolution: Optional[Union[float, Tuple[float, float]]] = None, wkt_string: Optional[str] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Reprojects a raster to match another shapefile/raster’s Coordinate Reference System (CRS), or a custom CRS.

Parameters
  • in_raster – Input raster.

  • out_crs – A raster or shapefile with the desired CRS to reproject to.

  • resolution – Allows the output resolution to be overridden.

  • wkt_string – Allows the user to define the output CRS using a custom valid WKT string.

  • out_path – Defines a path to save the output raster.

Returns

The input raster reprojected to match the desired Coordinate Reference System (CRS).

reproject_shapefile(in_shapefile: Union[geopandas.GeoDataFrame, str, Path], out_crs: Optional[Union[xarray.DataArray, str, Path, geopandas.GeoDataFrame]] = None, wkt_string: Optional[str] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Reprojects a shapefile to match another shapefile/raster’s Coordinate Reference System (CRS), or a custom CRS.

Parameters
  • in_shapefile – Input shapefile/GeoDataFrame.

  • out_crs – A raster or shapefile with the desired CRS to reproject to.

  • wkt_string – Allows the user to define the output CRS using a custom valid WKT string.

  • out_path – Defines a path to save the output shapefile.

Returns

The input shapefile reprojected to match the desired Coordinate Reference System (CRS).

resample(in_raster: Union[xarray.DataArray, str, Path], match_raster: Union[xarray.DataArray, str, Path], method: str = 'nearest', out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Resamples a raster to match another raster’s cell size.

Parameters
  • in_raster – Input raster.

  • match_raster – A raster to match the resolution / cell size of.

  • method – A valid resampling method from rasterio.enums.Resample. NOTE: Do not use any averaging resample methods when working with a categorical raster!

  • out_path – Defines a path to save the output raster.

Returns

The input raster resampled to the desired resolution / cell size.

save_raster(out_raster: xarray.DataArray, out_path: Union[str, Path]) None[source]

Saves an xarray.DataArray to a .tif raster file at location param:out_path

save_shapefile(out_shapefile: geopandas.GeoDataFrame, out_path: Union[str, Path]) None[source]

Saves an geopandas.GeoDataFrame to a .shp file at location param:out_path

spatial_mask(in_raster: Union[xarray.DataArray, str, Path], mask_shp: Union[geopandas.GeoDataFrame, str], inverse: bool = False, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

Applies a spatial mask via a shapefile to a raster, converting out of bounds values by default to nodata.

Primarily for masking rasters (i.e., FAC) by basin shapefiles, converting out-of-mask raster values to NoData. A cell value can also be used to create a mask for integer rasters. Note: default behavior (inverse=False) will make it so cells NOT COVERED by mask_shp -> NoData.

Parameters
  • in_raster – Input raster.

  • mask_shp – Shapefile used for masking.

  • inverse – If True, cells that ARE COVERED by mask_shp -> NoData.

  • out_path – Defines a path to save the output raster.

Returns

The output spatially masked raster.

value_mask(in_raster: Union[xarray.DataArray, str, Path], thresh: Optional[Union[int, float]] = None, greater_than: bool = True, equals: Optional[int] = None, inverse_equals: bool = False, in_mask_value: Optional[int] = None, out_mask_value: Optional[int] = None, out_path: Optional[Union[str, Path]] = None) xarray.DataArray[source]

“Mask a raster via a value threshold.

Primary use case is to identify high accumulation zones / stream cells. Cells included in the mask are given a value of 1, all other cells are given a value of 0 (unless out_mask_value!=None). Note: this function generalizes V1:pyfunc:makeStreams() functionality.

Parameters
  • in_raster – Input raster.

  • thresh – The threshold to apply the value mask with.

  • greater_than – If False, only values less than param:thresh are included in the mask.

  • equals – Only cells matching the value of param:equals are included in the mask.

  • inverse_equals – Is True and param:equals==True, values NOT equal to :param:thresh are masked out.

  • in_mask_value – Allows included cells to be given a non-zero integer value.

  • out_mask_value – Allows non-included cells to be given a non-zero integer value.

  • out_path – Defines a path to save the output raster.

Returns

out_mask_value.

Return type

The output raster with all masked out values == nodata or param