Skip to content

Index

Reference - Exported functions

Geomorphometry.D8 Type

D8 Flow Direction method by Jenson (1988).

source
Geomorphometry.D8D Type
julia
D8D <: FlowDirectionConvention

D8 flow direction convention using power-of-2 encoding, clockwise from East:

julia
 32(↖) 64() 128(↗)
 16()  0(·)   1()
  8(↙)  4()   2(↘)

Axis convention: dim1 = x (East+), dim2 = y (North+).

julia
# 8 16  32    W
# 4  0  64  S   N
# 2  1 128    E

Values can be combined with bitwise OR to represent multiple flow directions (e.g., 1 | 2 | 4 = 7 for E+SE+S). Use decompose to extract individual directions.

source
Geomorphometry.DInf Type

DInf Flow Direction method by Tarboton (1997).

source
Geomorphometry.FD8 Type

FD8 Flow Direction method by Quinn (1991).

source
Geomorphometry.FlowDirection Type
julia
FlowDirection{C<:FlowDirectionConvention, T<:Integer} <: Integer

A flow direction value in convention C, stored as type T.

Constructors

  • FlowDirection{C}(v::Integer): Create from a raw direction value.

  • FlowDirection{C}(ci::CartesianIndex{2}): Create from a CartesianIndex offset.

Conversion

  • CartesianIndex(d::Direction): Convert to CartesianIndex offset.

  • convert(FlowDirection{C2}, d::FlowDirection{C1}): Convert between conventions.

  • Int(d::Direction): Get the raw integer value.

source
Geomorphometry.Horn Type

Third order finite difference estimator using all 8 neighbors by Horn, (1981).

source
Geomorphometry.LDD Type
julia
LDD <: FlowDirectionConvention

Local Drainage Direction (PCRaster) convention using 1-9 numpad encoding:

julia
7(↖) 8() 9(↗)
4() 5(·) 6()
1(↙) 2() 3(↘)

Axis convention: dim1 = x (East+), dim2 = y (North+). The table equals centered(reshape(1:9, 3, 3)):

julia
# 1 4 7    W
# 2 5 8  S   N
# 3 6 9    E
source
Geomorphometry.ZevenbergenThorne Type

Second order finite difference estimator using all 4 neighbors by Zevenbergen and Thorne, (1987).

source
Geomorphometry.aspect Method
julia
aspect(dem::Matrix{<:Real}, method=Horn())

Aspect is direction of slope.

source
Geomorphometry.bathymetric_position_index Function
julia
bathymetric_position_index(dem::AbstractMatrix{<:Real})

Bathymetric Position Index (BPI) (Lundblad et al., 2006), which is defined as the difference between a central pixel and the mean of the cells in an annulus around it.

source
Geomorphometry.depression_depth Method
julia
depression_depth(dem::AbstractMatrix; filled=filldepressions(dem))

Computes the depth of each cell below the filled surface.

Returns the difference between the depression-filled DEM and the original DEM, representing how deep each cell sits within a depression. Cells not in depressions will have a depth of zero.

This is useful for identifying potential cold air pooling zones and water retention areas.

source
Geomorphometry.depression_volume Method
julia
depression_volume(dem::AbstractMatrix; filled=filldepressions(dem), cellsize=cellsize(dem))

Computes the total volume of all depressions/basins in the DEM.

Returns the sum of depression depths multiplied by cell area.

source
Geomorphometry.drainage_potential Method
julia
drainage_potential(dem::AbstractMatrix; method=DInf(), cellsize=cellsize(dem))

Computes a drainage potential index indicating how well each cell drains.

High values indicate good drainage (steep slopes with low upstream accumulation). Low values indicate poor drainage (flat areas that accumulate flow from upslope).

This is useful as a proxy for cold air drainage - cells with high drainage potential will shed cold air downslope rather than pooling it.

Based on the relationship between slope and flow accumulation: drainage = sin(slope) / log(1 + accumulation)

source
Geomorphometry.entropy Method
julia
entropy(dem::AbstractMatrix{<:Real}; step=0.5)

Entropy calculates the Shannon entropy of the surrounding cells of a central cell. step is the bin size for the histogram.

source
Geomorphometry.filldepressions Function
julia
filldepressions(dem::AbstractMatrix, mask::Matrix{Bool})

Performs the Priority-Flood algorithm (Barnes et al., 2014) on the given digital elevation model (DEM) dem with an optional mask.

Arguments

  • dem::AbstractMatrix: A 2D array representing the digital elevation model (DEM).

  • mask::AbstractMatrix{Bool}: A 2D boolean array representing the mask. Cells with true values are considered in the computation, while cells with false values are ignored.

source
Geomorphometry.flowaccumulation Function
julia
flowaccumulation(dem::AbstractMatrix, closed::Matrix{Bool}, method::FlowDirectionMethod)

Computes the flow accumulation of a digital elevation model (DEM) dem with an optional closed mask and a method for flow direction. Returns the flow accumulation and the flow direction (local drainage direction or ldd)

source
Geomorphometry.height_above_nearest_drainage Method
julia
height_above_nearest_drainage(dem::AbstractMatrix, stream_mask::AbstractMatrix{Bool})

Computes the Height Above Nearest Drainage (HAND, (Nobre et al., 2011)) of a digital elevation model (DEM) dem given a stream definition as a boolean stream_mask.

source
Geomorphometry.height_above_nearest_drainage Method
julia
height_above_nearest_drainage(dem::AbstractMatrix; method=DInf(), cellsize=cellsize(dem), threshold=1e10)

Compute Height Above Nearest Drainage (HAND, (Nobre et al., 2011)) of a digital elevation model (DEM) dem with an optional method for flow direction, a cellsize, and a flow accumulation threshold for stream definition.

source
Geomorphometry.hillshade Method
julia
hillshade(dem::Matrix{<:Real}; azimuth=315.0, zenith=45.0, cellsize=cellsize(dem))

hillshade is the simulated illumination of a surface based on its slope and aspect given a light source with azimuth and zenith angles in °, as defined in Burrough, P. A., and McDonell, R. A., (1998).

source
Geomorphometry.horizon_angle Method
julia
horizon_angle(dem; method=GridSweep(), cellsize=cellsize(dem))

Compute horizon angles for each cell in a DEM.

The horizon angle in a given direction is the maximum elevation angle from the cell to any point along that direction. Positive angles indicate terrain rising above the horizontal (sheltered), negative angles indicate terrain falling below (exposed).

Arguments

  • dem: Digital elevation model matrix

Keywords

  • directions: 16 by default, can be 4, 8, 16, 32, ...

  • cellsize: Cell size as (row_size, col_size) tuple

Returns

3D array of size (rows, cols, directions) with horizon angles in degrees.

source
Geomorphometry.multihillshade Method
julia
multihillshade(dem::AbstractMatrix{<:Real}; cellsize=cellsize(dem))

multihillshade is the simulated illumination of a surface based on its slope and aspect. Like hillshade, but now using multiple sources as defined in Mark, R.K. (1992), similar to GDALs -multidirectional.

source
Geomorphometry.percentile_elevation Method
julia
percentile_elevation(dem::AbstractMatrix{<:Real}; radius=1)

Computes the percentile rank of each cell's elevation relative to its neighborhood.

Returns a value between 0 and 1 for each cell, where:

  • 0 means the cell is lower than all neighbors (local minimum / valley bottom)

  • 1 means the cell is higher than all neighbors (local maximum / ridge top)

  • 0.5 means the cell is at the median elevation of its neighborhood

This metric is useful for identifying cold air pooling zones, which tend to occur in areas with low percentile elevation (valleys and depressions).

Based on Lundquist et al. (2008).

Arguments

  • dem::AbstractMatrix{<:Real}: Digital elevation model

  • radius::Int=1: Radius of the circular neighborhood in cells

Example

julia
dem = rand(100, 100) .* 1000  # Random terrain
pct = percentile_elevation(dem; radius=20)
valleys = pct .< 0.3  # Low-lying areas prone to cold air pooling
source
Geomorphometry.pitremoval Method
julia
pitremoval(dem::AbstractMatrix{<:Real})

Remove pits from a DEM Array if the center cell of a 3x3 patch is limit lower than the surrounding cells.

source
Geomorphometry.plan_curvature Method
julia
plan_curvature(dem::AbstractMatrix{<:Real}; cellsize = cellsize(dem), radius=1)

Calculate projected contour curvature (plan curvature) as defined by Minár et al., (2020).

source
Geomorphometry.profile_curvature Method
julia
profile_curvature(dem::AbstractMatrix{<:Real}; cellsize = cellsize(dem), radius=1)

Calculate normal slope line curvature (profile curvature) as defined by Minár et al., (2020).

source
Geomorphometry.progressive_morphological_filter Method
julia
B, flags = progressive_morphological_filter(A; ωₘ, slope, dhₘ, dh₀, cellsize, adjust, erode)

Applies the progressive morphological filter by Zhang (2003) to A.

Output

  • B::Array{T,2} Maximum allowable values

  • flags::Array{Float64,2} A sized array with window sizes if filtered, zero if not filtered.

Afterwards, one can retrieve the resulting mask for A by A .<= B or flags .== 0..

Arguments

  • A::Array{T,2} Input Array

  • ωₘ::Real=20. Maximum window size [m]

  • slope::Real=0.01 Terrain slope [m/m]

  • dhₘ::Real=2.5 Maximum elevation threshold [m]

  • dh₀::Real=0.2 Initial elevation threshold [m]

  • cellsize::Real=1. Cellsize in [m]

source
Geomorphometry.pssm Method
julia
image = pssm(dem; exaggeration=2.3, resolution=1.0)

Perceptually Shaded Slope Map by Pingel and Clarke (2014).

Output

  • image::Gray{T,2} Grayscale image

Arguments

  • A::Array{Real,2} Input Array

  • exaggeration::Real=2.3 Factor to exaggerate elevation

  • cellsize::Real=1.0 Size of cell to account for horizontal resolution if different from vertical resolution

source
Geomorphometry.roughness Function
julia
roughness(dem::AbstractMatrix{<:Real})

Roughness is the largest inter-cell difference of a central pixel and its surrounding cell, as defined in Wilson et al. (2007).

source
Geomorphometry.roughness_index_elevation Function
julia
roughness_index_elevation(dem::AbstractMatrix{<:Real})

Roughness Index Elevation (RIE), which quantifies the standard deviation of residual topography (Cavalli et al., 2008)

source
Geomorphometry.rugosity Method
julia
rugosity(dem::AbstractMatrix{<:Real})

Compute the rugosity of a DEM, which is the ratio between the surface area divided by the planimetric area.

Jenness 2019 https://onlinelibrary.wiley.com/doi/abs/10.2193/0091-7648(2004)032[0829%3ACLSAFD]2.0.CO%3B2

source
Geomorphometry.simple_morphological_filter Method
julia
B = simple_morphological_filter(A; ω, slope, dhₘ, dh₀, cellsize)

Applies the simple morphological filter by Pingel et al. (2013) to A.

Output

  • B::Array{Float64,2} A filtered version of A

Arguments

  • A::Array{T,2} Input Array

  • ω::Float64=18. Maximum window size [m]

  • slope::Float64=0.01 Terrain slope [m/m]

  • cellsize::Float64=1. Cellsize in [m]

source
Geomorphometry.skewness_balancing Method
julia
mask = skewness_balancing(A; mean=mean(A))

Applies skewness balancing by Bartels et al. (2006) to A. Improved the performance by applying a binary search to find the threshold value.

Output

  • mask::BitMatrix Mask of allowed values
source
Geomorphometry.sky_view_factor Method
julia
sky_view_factor(dem; directions=16, cellsize=cellsize(dem))

Compute the Sky View Factor (SVF) for each cell in a DEM.

SVF is the fraction of the sky hemisphere visible from each point, ranging from 0 (fully obstructed) to 1 (full sky visible). It is computed as the mean of cos²(horizon_angle) across all directions.

Arguments

  • dem: Digital elevation model matrix

Keywords

  • directions: Number of directions (default: 16)

  • cellsize: Cell size as (row_size, col_size) tuple

Returns

A matrix of SVF values in the range [0, 1].

source
Geomorphometry.slope Method
julia
slope(dem::Matrix{<:Real}; cellsize=cellsize(dem), method=Horn(), exaggeration=1.0)

Slope is the rate of change between a cell and its neighbors.

source
Geomorphometry.spread Method
julia
spread(points::Matrix{<:Real}, initial::Matrix{<:Real}, friction::Real; distance=Euclidean(), res=1.0)
spread(points::Matrix{<:Real}, initial::Real, friction::Real; distance=Euclidean(), res=1.0)

Optimized (and more accurate) function based on the same friction everywhere.

When the friction is the same everywhere, there's no need for searching the shortest cost path, as one can just take a direct line to the input points.

The calculated cost is more accurate, as there's no 'zigzag' from cell center to cell center.

source
Geomorphometry.spread Method
julia
spread(::Eastman, points::Matrix{<:Real}, initial::Matrix{<:Real}, friction::Matrix{<:Real}; res=1, limit=Inf, iterations=3)

Pushbroom method for friction costs as discussed by Eastman (1989). This method should scale better (linearly) than the Tomlin (1983) method, but can require more iterations than set by default (3) in the case of maze-like, uncrossable obstacles.

Output

  • Array{Float64,2} Total friction distance

Arguments

  • points::Matrix{<:Real} Input Array

  • initial::Matrix{<:Real} Factor to exaggerate elevation

  • friction::Matrix{<:Real} Resolution of cell size

  • res=1 Resolution or cell size

  • limit=Inf Initial fill value

source
Geomorphometry.spread Method
julia
spread(points::Matrix{<:Real}, initial::Matrix{<:Real}, friction::Matrix{<:Real}; cellsize=(1,1), limit=Inf, method=Tomlin())

Total friction distance spread from points from initial with friction. By default uses Tomlin, see SpreadMethod for other algorithms.

source
Geomorphometry.spread Method
julia
spread(::Tomlin, points::Matrix{<:Real}, initial::Matrix{<:Real}, friction::Matrix{<:Real}; res=1, limit=Inf, method=Tomlin())

Total friction distance spread from points as described by Tomlin (1983). This is also the method implemented by PCRaster.

Output

  • Array{Float64,2} Total friction distance

Arguments

  • points::Matrix{<:Real} Input Array

  • initial::Matrix{<:Real} Initial values of the result

  • friction::Matrix{<:Real} Resolution of cell size

  • res=1 Resolution or cell size

  • limit=Inf Initial fill value

source
Geomorphometry.stream_power_index Method
julia
stream_power_index(dem::AbstractMatrix; method=DInf(), cellsize=cellsize(dem))

Computes the Stream Power Index (SPI) of a digital elevation model (DEM) dem with an optional method for flow direction and a cellsize.

source
Geomorphometry.tangential_curvature Method
julia
tangential_curvature(dem::AbstractMatrix{<:Real}; cellsize = cellsize(dem), radius=1)

Calculate normal contour curvature (tangential curvature) as defined by Minár et al., (2020).

source
Geomorphometry.terrain_ruggedness_index Method
julia
terrain_ruggedness_index(dem::AbstractMatrix{<:Real})

Terrain Ruggedness Index (TRI), which measures the difference between a central pixel and its surrounding cells. This algorithm uses the square root of the sum of the square of the absolute difference between a central pixel and its surrounding cells. This is recommended for terrestrial use cases.

source
Geomorphometry.topographic_position_index Function
julia
topographic_position_index(dem::AbstractMatrix{<:Real})

Topographic Position Index (TPI), which is defined as the difference between a central pixel and the mean of its surrounding cells, as defined in Wilson et al. (2007).

source
Geomorphometry.topographic_wetness_index Method
julia
topographic_wetness_index(dem::AbstractMatrix; method=DInf(), cellsize=cellsize(dem))

Computes the Topographic Wetness Index (TWI) of a digital elevation model (DEM) dem with an optional method for flow direction and a cellsize.

source

Reference - Internal functions

Geomorphometry.FlowDirectionConvention Type
julia
FlowDirectionConvention

Abstract type for flow direction encoding conventions. Subtypes define how direction integers map to neighbor offsets. See LDD.

source
Geomorphometry.MaximumDownwardGradient Type

Maximum Downward Gradient

source
Geomorphometry.SpreadMethod Type

Spread algorithms.

source
Geomorphometry._orient Method
julia
_orient(ci::CartesianIndex{2}, cellsize)

Convert a pixel CartesianIndex offset to the table convention (dim1=East+, dim2=North+). Accounts for the sign of cellsize: a GeoTIFF typically has negative cellsize[2] (+dim2 = South), so the second component is flipped. This function is its own inverse.

source
Geomorphometry.cellsize Method
julia
cellsize(dem)

Return an Tuple with the x and y length of each cell of the dem. Set them negatively to flip the image.

source
Geomorphometry.convention Method

Return the FlowDirectionConvention type of a direction.

source
Geomorphometry.decompose Method
julia
decompose(d::Direction)

Decompose a direction into a tuple of individual single-direction values. For LDD, returns a 1-tuple. For D8D, extracts each set bit.

Examples

julia
decompose(FlowDirection{D8D}(7))   # (→, ↘, ↓) → E + SE + S
decompose(FlowDirection{LDD}(9))   # (↘,) → SE
source
Geomorphometry.entropy! Method
julia
entropy!(dem::AbstractMatrix{<:Real})

Entropy calculates the Shannon entropy of the surrounding cells of a central cell.

source
Geomorphometry.ismulti Method

Whether a convention supports encoding multiple directions in a single value.

source
Geomorphometry.ispit Method

Whether this direction represents a pit (no outflow).

source
Geomorphometry.issingle Method

Whether this direction encodes exactly one flow direction (or pit).

source
Geomorphometry.ndirections Method

Number of flow directions encoded (0 for pit).

source
Geomorphometry.opening! Method

Apply the opening operation to A with window size ω.

source
Geomorphometry.opening_circ! Method

Apply the opening operation to A with window size ω.

source
Geomorphometry.prominence Method
julia
prominence(dem::AbstractMatrix{<:Real})

Prominence calculates the number of cells that are lower or equal than the central cell. Thus, 8 is a local maximum (peak), while 0 is a local minimum (pit).

source
Geomorphometry.round_odd Method
julia
round_odd(x)

Rounds x to the nearest odd number.

source
Geomorphometry.skbr Method
julia
mask = skbr(A; iterations=10)

Applies recursive skewness balancing by Bartels et al. (2010) to A. Applies skewness_balancing iterations times to the object (non-terrain) mask, as to include more (sloped) terrain.

Output

  • mask::BitMatrix Mask of allowed values
source