darts_preprocessing.engineering.arcticdem
¶
Computation of ArcticDEM derived products.
Distance
dataclass
¶
Convenience class to represent a distance in pixels and meters.
__repr__
¶
parse
classmethod
¶
parse(
v: int | float | str, res: float
) -> darts_preprocessing.engineering.arcticdem.Distance
Parse a distance from a string or numeric value.
If the input is a string, it can be in the format of "10px" or "10m". If it is a numeric value, it is interpreted as meters and converted to pixels based on the resolution.
Parameters:
Raises:
-
ValueError–If the input distance is negative.
-
ValueError–If the input distance is not a valid string format.
-
TypeError–If the input distance is not a string, int, or float.
Returns:
-
Distance(darts_preprocessing.engineering.arcticdem.Distance) –The parsed distance in pixels and meters.
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
calculate_aspect
¶
Calculate aspect (compass direction) of the terrain surface from an ArcticDEM Dataset.
Aspect indicates the downslope direction of the maximum rate of change in elevation.
Parameters:
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
aspect (float32): Aspect in degrees clockwise from north [0-360], or -1 for flat areas.
-
long_name: "Aspect"
- units: "degrees"
- description: Compass direction of slope
- source: "ArcticDEM"
-
Note
Aspect values:
- 0° or 360°: North-facing
- 90°: East-facing
- 180°: South-facing
- 270°: West-facing
- -1: Flat (no dominant direction)
Example
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
calculate_curvature
¶
Calculate curvature of the terrain surface from an ArcticDEM Dataset.
Curvature measures the rate of change of slope, indicating terrain convexity or concavity.
Parameters:
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
curvature (float32): Curvature values.
-
long_name: "Curvature"
- description: Rate of change of slope
- source: "ArcticDEM"
-
Note
Curvature interpretation:
- Positive values: Convex terrain (hills, ridges)
- Negative values: Concave terrain (valleys, depressions)
- Near zero: Planar terrain
Example
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
calculate_dissection_index
¶
calculate_dissection_index(
arcticdem_ds: xarray.Dataset, neighborhood_size: int
) -> xarray.Dataset
Calculate the Dissection Index (DI) from an ArcticDEM Dataset.
DI measures the degree to which a landscape has been cut by valleys and ravines. Values range from 0 (smooth, undissected) to 1 (highly dissected).
Parameters:
-
arcticdem_ds(xarray.Dataset) –Dataset containing: - dem (float32): Digital Elevation Model
-
neighborhood_size(int) –Neighborhood window size for DI calculation. Can be specified as string with units (e.g., "100m" or "10px").
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
di (float32): Dissection Index [0-1].
-
long_name: "Dissection Index"
- description: Documents neighborhood size used
- source: "ArcticDEM"
-
Note
The dissection index quantifies landscape dissection by comparing elevation ranges within the neighborhood window. Higher values indicate more deeply incised terrain with greater vertical relief.
The neighborhood_size parameter is converted to pixels based on DEM resolution.
Example
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
calculate_hillshade
¶
calculate_hillshade(
arcticdem_ds: xarray.Dataset,
azimuth: int = 225,
angle_altitude: int = 25,
) -> xarray.Dataset
Calculate hillshade of the terrain surface from an ArcticDEM Dataset.
Hillshade simulates illumination of terrain from a specified sun position, useful for visualization and terrain analysis.
Parameters:
-
arcticdem_ds(xarray.Dataset) –Dataset containing: - dem (float32): Digital Elevation Model
-
azimuth(int, default:225) –Light source azimuth in degrees clockwise from north [0-360]. Defaults to 225 (southwest).
-
angle_altitude(int, default:25) –Light source altitude angle in degrees above horizon [0-90]. Defaults to 25.
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
hillshade (float32): Illumination values [0-255], where 0 is shadow and 255 is fully lit.
-
long_name: "Hillshade"
- description: Documents azimuth and angle_altitude used
- source: "ArcticDEM"
-
Note
Common azimuth/altitude combinations:
- 315°/45°: Classic northwest illumination (default for many GIS applications)
- 225°/25°: Southwest with low sun (better for visualizing subtle features)
The hillshade calculation accounts for both slope and aspect of the terrain.
Example
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
calculate_slope
¶
Calculate slope of the terrain surface from an ArcticDEM Dataset.
Slope represents the rate of change of elevation, indicating terrain steepness.
Parameters:
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
slope (float32): Slope in degrees [0-90].
-
long_name: "Slope"
- units: "degrees"
- source: "ArcticDEM"
-
Note
Slope is calculated using finite difference methods on the DEM. Values approaching 90° indicate near-vertical terrain.
Example
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
calculate_terrain_ruggedness_index
¶
calculate_terrain_ruggedness_index(
arcticdem_ds: xarray.Dataset, neighborhood_size: int
) -> xarray.Dataset
Calculate the Terrain Ruggedness Index (TRI) from an ArcticDEM Dataset.
TRI quantifies topographic heterogeneity by measuring elevation differences between a cell and its surrounding cells. Higher values indicate more rugged terrain.
Parameters:
-
arcticdem_ds(xarray.Dataset) –Dataset containing: - dem (float32): Digital Elevation Model
-
neighborhood_size(int) –Neighborhood window size for TRI calculation. Can be specified as string with units (e.g., "100m" or "10px").
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
tri (float32): Terrain Ruggedness Index in meters.
-
long_name: "Terrain Ruggedness Index"
- units: "m"
- description: Documents kernel size used
- source: "ArcticDEM"
-
Note
TRI methodology from Riley et al (1999):
- Measures elevation difference from center cell to 8 surrounding cells
- Squares and averages these differences
- Takes square root for final TRI value
The neighborhood_size parameter controls the kernel size. A square kernel is used, with the actual size rounded to the nearest pixel based on DEM resolution.
References
Riley, S.J., DeGloria, S.D., Elliot, R., 1999. A Terrain Ruggedness Index That Quantifies Topographic Heterogeneity. Intermountain Journal of Sciences, vol. 5, No. 1-4, pp. 23-27.
Example
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 | |
calculate_topographic_position_index
¶
calculate_topographic_position_index(
arcticdem_ds: xarray.Dataset,
outer_radius: int,
inner_radius: int,
) -> xarray.Dataset
Calculate the Topographic Position Index (TPI) from an ArcticDEM Dataset.
TPI measures the relative topographic position of a point by comparing its elevation to the mean elevation of the surrounding neighborhood. Positive values indicate higher positions (ridges), negative values indicate lower positions (valleys).
Parameters:
-
arcticdem_ds(xarray.Dataset) –The ArcticDEM Dataset containing the 'dem' variable (float32).
-
outer_radius(int) –The outer radius of the neighborhood in meters. Can also be specified as string with units (e.g., "100m" or "10px").
-
inner_radius(int) –The inner radius of the annulus kernel in meters. If > 0, creates an annulus (ring) instead of a circle. Set to 0 for a circular kernel. Can also be specified as string with units (e.g., "50m" or "5px").
Returns:
-
xarray.Dataset–xr.Dataset: The input Dataset with a new data variable added:
-
xarray.Dataset–-
tpi (float32): Topographic Position Index values.
-
long_name: "Topographic Position Index (TPI)"
- description: Details about the kernel used
-
Note
Kernel shape combinations:
- inner_radius=0: Circular kernel comparing each cell to all neighbors within outer_radius
- inner_radius>0: Annulus kernel comparing each cell to neighbors in a ring between inner_radius and outer_radius. Useful for multi-scale terrain analysis.
The actual radii used are rounded to the nearest pixel based on the DEM resolution.
Example
Calculate TPI with circular and annulus kernels:
from darts_preprocessing import calculate_topographic_position_index
# Circular kernel (100m radius)
arcticdem_with_tpi = calculate_topographic_position_index(
arcticdem_ds=arcticdem,
outer_radius=100,
inner_radius=0
)
# Annulus kernel (50-100m ring)
arcticdem_multi_scale = calculate_topographic_position_index(
arcticdem_ds=arcticdem,
outer_radius=100,
inner_radius=50
)
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 | |
calculate_vector_ruggedness_measure
¶
calculate_vector_ruggedness_measure(
arcticdem_ds: xarray.Dataset, neighborhood_size: int
) -> xarray.Dataset
Calculate the Vector Ruggedness Measure (VRM) from an ArcticDEM Dataset.
VRM quantifies terrain ruggedness using vector analysis of slope and aspect, providing a measure independent of absolute elevation. Values range from 0 (smooth) to 1 (rugged).
Parameters:
-
arcticdem_ds(xarray.Dataset) –Dataset containing: - dem (float32): Digital Elevation Model - slope (float32): Slope in degrees (will be calculated if not present) - aspect (float32): Aspect in degrees (will be calculated if not present)
-
neighborhood_size(int) –Neighborhood window size for VRM calculation. Can be specified as string with units (e.g., "100m" or "10px").
Returns:
-
xarray.Dataset–xr.Dataset: Input Dataset with new data variable added:
-
xarray.Dataset–-
vrm (float32): Vector Ruggedness Measure [0-1].
-
long_name: "Vector Ruggedness Measure"
- description: Documents neighborhood size used
- source: "ArcticDEM"
-
Note
VRM calculation:
- Converts slope and aspect to 3D unit vectors (x, y, z components)
- Sums vectors in the neighborhood window
- Calculates magnitude of resultant vector
- VRM = 1 - resultant magnitude
Flat areas (aspect = -1) are handled by setting aspect to 0.
Requires slope and aspect to be already calculated on the dataset.
References
Sappington, J.M., K.M. Longshore, and D.B. Thomson. 2007. Quantifying Landscape Ruggedness for Animal Habitat Analysis: A case Study Using Bighorn Sheep in the Mojave Desert. Journal of Wildlife Management. 71(5): 1419-1426.
Example
from darts_preprocessing import (
calculate_slope, calculate_aspect,
calculate_vector_ruggedness_measure
)
# VRM requires slope and aspect
arcticdem = calculate_slope(arcticdem)
arcticdem = calculate_aspect(arcticdem)
arcticdem_with_vrm = calculate_vector_ruggedness_measure(
arcticdem_ds=arcticdem,
neighborhood_size=100
)
Source code in darts-preprocessing/src/darts_preprocessing/engineering/arcticdem.py
455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 | |
dissection_index
¶
dissection_index(
agg: xarray.DataArray,
window_size: int = 3,
name: str | None = "dissection_index",
) -> xarray.DataArray
Compute the dissection index of a 2D array.
Parameters:
-
agg(xarray.DataArray) –The input data array.
-
window_size(int, default:3) –The size of the window to use for the computation. Defaults to 3.
-
name(str | None, default:'dissection_index') –The name of the output data array. Defaults to "dissection_index".
Returns: