Reference
Geomorphometry.TPI
— MethodTPI(dem::Matrix{<:Real})
TPI stands for Topographic Position Index, which is defined as the difference between a central pixel and the mean of its surrounding cells (see Wilson et al 2007, Marine Geodesy 30:3-35).
Geomorphometry.TRI
— MethodTRI(dem::Matrix{<:Real})
TRI stands for Terrain Ruggedness Index, 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 difference between a central pixel and its surrounding cells. This is recommended for terrestrial use cases.
Geomorphometry.aspect
— Methodaspect(dem::Matrix{<:Real}, method=Horn())
Aspect is direction of slope
, as defined in Burrough, P. A., and McDonell, R. A., (1998, Principles of Geographical Information Systems).
Geomorphometry.curvature
— Methodcurvature(dem::Matrix{<:Real}; cellsize=1.0)
Curvature is derivative of slope
, so the second derivative of the DEM.
Geomorphometry.hillshade
— Methodhillshade(dem::Matrix{<:Real}; azimuth=315.0, zenith=45.0, cellsize=1.0)
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, Principles of Geographical Information Systems).
Geomorphometry.multihillshade
— Methodmultihillshade(dem::Matrix{<:Real}; cellsize=1.0)
multihillshade is the simulated illumination of a surface based on its slope
and aspect
. Like hillshade
, but now using multiple sources as defined in https://pubs.usgs.gov/of/1992/of92-422/of92-422.pdf, similar to GDALs -multidirectional.
Geomorphometry.opening!
— MethodApply the opening operation to A
with window size ω
.
Geomorphometry.opening
— MethodApply the opening operation to A
with window size ω
.
Geomorphometry.opening_circ!
— MethodApply the opening operation to A
with window size ω
.
Geomorphometry.pitremoval
— Methodpitremoval(dem::Matrix{<:Real})
Remove pits from a DEM Array if the center cell of a 3x3 patch is limit
lower or than the surrounding cells.
Geomorphometry.pmf
— MethodB, flags = pmf(A; ωₘ, slope, dhₘ, dh₀, cellsize, adjust, erode)
Applies the progressive morphological filter by Zhang et al. (2003) [zhang2003] to A
.
Output
B::Array{T,2}
Maximum allowable valuesflags::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]
Geomorphometry.pssm
— Methodimage = pssm(A; exaggeration=2.3, resolution=1.0)
Perceptually Shaded Slope Map by Pingel, Clarke. 2014 [pingel2014].
Output
image::Gray{T,2}
Grayscale image
Arguments
A::Array{Real,2}
Input Arrayexaggeration::Real=2.3
Factor to exaggerate elevationcellsize::Real=1.0
Size of cell to account for horizontal resolution if different from vertical resolution
Geomorphometry.roughness
— Methodroughness(dem::Matrix{<:Real})
Roughness is the largest inter-cell difference of a central pixel and its surrounding cell, as defined in Wilson et al (2007, Marine Geodesy 30:3-35).
Geomorphometry.skb
— Methodmask = skb(A; mean=mean(A))
Applies skewness balancing by Bartels e.a (2006) [bartels2006] to A
. Improved the performance by applying a binary search to find the threshold value.
Output
mask::BitMatrix
Mask of allowed values
Geomorphometry.skbr
— Methodmask = skbr(A; iterations=10)
Applies recursive skewness balancing by Bartels e.a (2006) [bartels2006] to A
. Applies skb
iterations
times to the object (non-terrain) mask, as to include more (sloped) terrain.
Output
mask::BitMatrix
Mask of allowed values
Geomorphometry.slope
— Methodslope(dem::Matrix{<:Real}; cellsize=1.0, method=Horn())
Slope is the rate of change between a cell and its neighbors as defined in Burrough, P. A., and McDonell, R. A., (1998, Principles of Geographical Information Systems).
Geomorphometry.smf
— MethodB = smf(A; ω, slope, dhₘ, dh₀, cellsize)
Applies the simple morphological filter by Pingel et al. (2013) [pingel2013] 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]
Geomorphometry.spread
— Methodspread(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.
Geomorphometry.spread
— Methodspread(points::Matrix{<:Real}, initial::Matrix{<:Real}, friction::Matrix{<:Real}; res=1, limit=Inf)
Total friction distance spread from points
as by Tomlin (1983) [tomlin1983]. This is also the method implemented by PCRaster.
Output
Array{Float64,2}
Total friction distance
Arguments
points::Matrix{<:Real}
Input Arrayinitial::Matrix{<:Real}
Factor to exaggerate elevationfriction::Matrix{<:Real}
Resolution of cell sizeres=1
Resolution or cell sizelimit=Inf
Initial fill value
Geomorphometry.spread2
— Methodspread2(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) [eastman1989]. This method should scale much better (linearly) than the [tomlin1983] 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 Arrayinitial::Matrix{<:Real}
Factor to exaggerate elevationfriction::Matrix{<:Real}
Resolution of cell sizeres=1
Resolution or cell sizelimit=Inf
Initial fill valueiterations=3
Number of pushbroom iterations
Geomorphometry.Horn
— TypeThird order finite difference estimator using all 8 neighbors (Horn, 1981).
Geomorphometry.MaximumDownwardGradient
— TypeMaximum Downward Gradient
Geomorphometry.ZevenbergenThorne
— TypeSecond order finite difference estimator using all 4 neighbors (Zevenbergen and Thorne, 1987).
- zhang2003Zhang, Keqi, Shu-Ching Chen, Dean Whitman, Mei-Ling Shyu, Jianhua Yan, and Chengcui Zhang. “A Progressive Morphological Filter for Removing Nonground Measurements from Airborne LIDAR Data.” IEEE Transactions on Geoscience and Remote Sensing 41, no. 4 (2003): 872–82. https://doi.org/10.1109/TGRS.2003.810682.
- pingel2014Pingel, Thomas, and Clarke, Keith. 2014. ‘Perceptually Shaded Slope Maps for the Visualization of Digital Surface Models’. Cartographica: The International Journal for Geographic Information and Geovisualization 49 (4): 225–40. https://doi.org/10/ggnthv.
- bartels2006Bartels, M., Hong Wei, and D.C. Mason. 2006. “DTM Generation from LIDAR Data Using Skewness Balancing.” In 18th International Conference on Pattern Recognition (ICPR’06), 1:566–69. https://doi.org/10/cwk4v2.
- bartels2006Bartels, M., Hong Wei, and D.C. Mason. 2006. “DTM Generation from LIDAR Data Using Skewness Balancing.” In 18th International Conference on Pattern Recognition (ICPR’06), 1:566–69. https://doi.org/10/cwk4v2.
- pingel2013Pingel, Thomas J., Keith C. Clarke, and William A. McBride. 2013. ‘An Improved Simple Morphological Filter for the Terrain Classification of Airborne LIDAR Data’. ISPRS Journal of Photogrammetry and Remote Sensing 77 (March): 21–30. https://doi.org/10.1016/j.isprsjprs.2012.12.002.
- tomlin1983Tomlin, Charles Dana. 1983. Digital Cartographic Modeling Techniques in Environmental Planning. Yale University.
- eastman1989Eastman, J. Ronald. 1989. ‘Pushbroom Algorithms for Calculating Distances in Raster Grids’. In Proceedings, Autocarto, 9:288–97.