Skip to content

Commit ede4c44

Browse files
authored
Add sky-view factor function (#962) (#972)
* Add sky-view factor function (#962) New xrspatial.sky_view_factor module with all four backends: NumPy (@ngjit), CuPy (@cuda.jit), Dask+NumPy, and Dask+CuPy. * Add tests for sky_view_factor, fix edge cell handling (#962) 25 tests covering correctness, NaN handling, range validation, cross-backend equivalence (NumPy/Dask/CuPy/Dask+CuPy), and Dataset support. Edge cells now compute with truncated rays instead of returning NaN, keeping numpy and dask paths consistent. * Add docs and user guide notebook for sky_view_factor (#962) - Add Sky-View Factor entry to docs/source/reference/surface.rst - Add examples/user_guide/22_Sky_View_Factor.ipynb with parameter exploration and hillshade comparison * Add sky_view_factor to README feature matrix (#962)
1 parent ec31dc6 commit ede4c44

File tree

6 files changed

+601
-0
lines changed

6 files changed

+601
-0
lines changed

README.md

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -267,6 +267,7 @@ In the GIS world, rasters are used for representing continuous phenomena (e.g. e
267267
| [Curvature](xrspatial/curvature.py) | Measures rate of slope change (concavity/convexity) at each cell | ✅️ |✅️ |✅️ | ✅️ |
268268
| [Hillshade](xrspatial/hillshade.py) | Simulates terrain illumination from a given sun angle and azimuth | ✅️ | ✅️ | ✅️ | ✅️ |
269269
| [Roughness](xrspatial/terrain_metrics.py) | Computes local relief as max minus min elevation in a 3×3 window | ✅️ | ✅️ | ✅️ | ✅️ |
270+
| [Sky-View Factor](xrspatial/sky_view_factor.py) | Measures the fraction of visible sky hemisphere at each cell | ✅️ | ✅️ | ✅️ | ✅️ |
270271
| [Slope](xrspatial/slope.py) | Computes terrain gradient steepness at each cell in degrees | ✅️ | ✅️ | ✅️ | ✅️ |
271272
| [Terrain Generation](xrspatial/terrain.py) | Generates synthetic terrain from fBm or ridged fractal noise with optional domain warping, Worley blending, and hydraulic erosion | ✅️ | ✅️ | ✅️ | ✅️ |
272273
| [TPI](xrspatial/terrain_metrics.py) | Computes Topographic Position Index (center minus mean of neighbors) | ✅️ | ✅️ | ✅️ | ✅️ |

docs/source/reference/surface.rst

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -39,6 +39,13 @@ Terrain Generation
3939

4040
xrspatial.terrain.generate_terrain
4141

42+
Sky-View Factor
43+
===============
44+
.. autosummary::
45+
:toctree: _autosummary
46+
47+
xrspatial.sky_view_factor.sky_view_factor
48+
4249
Viewshed
4350
========
4451
.. autosummary::
Lines changed: 101 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,101 @@
1+
{
2+
"cells": [
3+
{
4+
"cell_type": "markdown",
5+
"id": "0fnki2gcn3md",
6+
"source": "# Sky-View Factor\n\nSky-view factor (SVF) measures the fraction of the sky hemisphere visible from each cell in a digital elevation model. Values range from 0 (fully obstructed, e.g. a deep canyon) to 1 (flat open terrain with no horizon obstruction).\n\nSVF is used in:\n- **LiDAR archaeology** to reveal subtle terrain features hidden under canopy\n- **Urban heat island studies** to quantify street canyon geometry\n- **Solar energy modeling** as a proxy for diffuse sky irradiance\n- **Terrain visualization** as an illumination-independent alternative to hillshade\n\nThe algorithm casts rays at evenly spaced azimuths from each cell and records the maximum elevation angle to the horizon along each ray. SVF is then:\n\n$$\\text{SVF} = 1 - \\frac{1}{N}\\sum_{d=1}^{N} \\sin(\\theta_d)$$\n\nwhere $\\theta_d$ is the maximum horizon angle in direction $d$.",
7+
"metadata": {}
8+
},
9+
{
10+
"cell_type": "code",
11+
"id": "050wobw30odo",
12+
"source": "import numpy as np\nimport xarray as xr\nimport matplotlib.pyplot as plt\n\nfrom xrspatial import sky_view_factor",
13+
"metadata": {},
14+
"execution_count": null,
15+
"outputs": []
16+
},
17+
{
18+
"cell_type": "markdown",
19+
"id": "9w2idcivtrw",
20+
"source": "## Generate synthetic terrain\n\nWe'll create a terrain surface with valleys, ridges, and a flat plateau to show how SVF responds to different landforms.",
21+
"metadata": {}
22+
},
23+
{
24+
"cell_type": "code",
25+
"id": "yubj4xlrxd",
26+
"source": "# Create a 200x200 terrain with ridges and valleys\nrows, cols = 200, 200\ny = np.arange(rows, dtype=np.float64)\nx = np.arange(cols, dtype=np.float64)\nY, X = np.meshgrid(y, x, indexing='ij')\n\n# Sinusoidal ridges + a central peak\nterrain = (\n 40 * np.sin(X / 15) * np.cos(Y / 20)\n + 80 * np.exp(-((Y - 100)**2 + (X - 100)**2) / (2 * 30**2))\n + 20 * np.sin(Y / 10)\n + 200\n)\n\ndem = xr.DataArray(\n terrain,\n dims=['y', 'x'],\n coords={'y': np.arange(rows, dtype=float), 'x': np.arange(cols, dtype=float)},\n attrs={'res': (1.0, 1.0)},\n)\n\nfig, ax = plt.subplots(figsize=(8, 7))\ndem.plot(ax=ax, cmap='terrain', cbar_kwargs={'label': 'Elevation'})\nax.set_title('Synthetic DEM')\nax.set_aspect('equal')\nplt.tight_layout()\nplt.show()",
27+
"metadata": {},
28+
"execution_count": null,
29+
"outputs": []
30+
},
31+
{
32+
"cell_type": "markdown",
33+
"id": "35xi1fqq831",
34+
"source": "## Compute SVF with default parameters\n\nThe default uses `max_radius=10` cells and `n_directions=16` azimuth rays.",
35+
"metadata": {}
36+
},
37+
{
38+
"cell_type": "code",
39+
"id": "mrt8dfxrni",
40+
"source": "svf = sky_view_factor(dem, max_radius=10, n_directions=16)\n\nfig, ax = plt.subplots(figsize=(8, 7))\nsvf.plot(ax=ax, cmap='gray', vmin=0, vmax=1,\n cbar_kwargs={'label': 'Sky-View Factor'})\nax.set_title('SVF (max_radius=10, n_directions=16)')\nax.set_aspect('equal')\nplt.tight_layout()\nplt.show()",
41+
"metadata": {},
42+
"execution_count": null,
43+
"outputs": []
44+
},
45+
{
46+
"cell_type": "markdown",
47+
"id": "ot8twzbdqgq",
48+
"source": "## Effect of `max_radius`\n\nA larger search radius captures more distant obstructions. This matters in terrain with broad features like wide valleys or distant ridgelines.",
49+
"metadata": {}
50+
},
51+
{
52+
"cell_type": "code",
53+
"id": "pq6q4dvd2dd",
54+
"source": "radii = [5, 15, 30]\nfig, axes = plt.subplots(1, 3, figsize=(18, 5))\n\nfor ax, r in zip(axes, radii):\n result = sky_view_factor(dem, max_radius=r, n_directions=16)\n result.plot(ax=ax, cmap='gray', vmin=0, vmax=1, add_colorbar=False)\n ax.set_title(f'max_radius={r}')\n ax.set_aspect('equal')\n\nfig.suptitle('SVF at different search radii', fontsize=14, y=1.02)\nplt.tight_layout()\nplt.show()",
55+
"metadata": {},
56+
"execution_count": null,
57+
"outputs": []
58+
},
59+
{
60+
"cell_type": "markdown",
61+
"id": "l1rcmgtk8",
62+
"source": "## Effect of `n_directions`\n\nMore azimuth directions give a smoother result at the cost of longer computation. For most uses, 16 directions is a good balance.",
63+
"metadata": {}
64+
},
65+
{
66+
"cell_type": "code",
67+
"id": "b33ppidyce",
68+
"source": "directions = [4, 16, 64]\nfig, axes = plt.subplots(1, 3, figsize=(18, 5))\n\nfor ax, n in zip(axes, directions):\n result = sky_view_factor(dem, max_radius=15, n_directions=n)\n result.plot(ax=ax, cmap='gray', vmin=0, vmax=1, add_colorbar=False)\n ax.set_title(f'n_directions={n}')\n ax.set_aspect('equal')\n\nfig.suptitle('SVF at different direction counts', fontsize=14, y=1.02)\nplt.tight_layout()\nplt.show()",
69+
"metadata": {},
70+
"execution_count": null,
71+
"outputs": []
72+
},
73+
{
74+
"cell_type": "markdown",
75+
"id": "n7hl0m8uhmq",
76+
"source": "## Comparison: SVF vs Hillshade\n\nHillshade depends on a specific sun angle, which creates directional bias. SVF provides omnidirectional illumination information, making features visible regardless of orientation.",
77+
"metadata": {}
78+
},
79+
{
80+
"cell_type": "code",
81+
"id": "v8qltgseoba",
82+
"source": "from xrspatial import hillshade\n\nhs = hillshade(dem, azimuth=315, angle_altitude=45)\nsvf_result = sky_view_factor(dem, max_radius=15, n_directions=16)\n\nfig, axes = plt.subplots(1, 2, figsize=(14, 6))\n\nhs.plot(ax=axes[0], cmap='gray', add_colorbar=False)\naxes[0].set_title('Hillshade (azimuth=315)')\naxes[0].set_aspect('equal')\n\nsvf_result.plot(ax=axes[1], cmap='gray', vmin=0, vmax=1, add_colorbar=False)\naxes[1].set_title('Sky-View Factor')\naxes[1].set_aspect('equal')\n\nplt.tight_layout()\nplt.show()",
83+
"metadata": {},
84+
"execution_count": null,
85+
"outputs": []
86+
}
87+
],
88+
"metadata": {
89+
"kernelspec": {
90+
"display_name": "Python 3",
91+
"language": "python",
92+
"name": "python3"
93+
},
94+
"language_info": {
95+
"name": "python",
96+
"version": "3.10.0"
97+
}
98+
},
99+
"nbformat": 4,
100+
"nbformat_minor": 5
101+
}

xrspatial/__init__.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -73,6 +73,7 @@
7373
from xrspatial.snap_pour_point import snap_pour_point # noqa
7474
from xrspatial.stream_link import stream_link # noqa
7575
from xrspatial.stream_order import stream_order # noqa
76+
from xrspatial.sky_view_factor import sky_view_factor # noqa
7677
from xrspatial.slope import slope # noqa
7778
from xrspatial.surface_distance import surface_allocation # noqa
7879
from xrspatial.surface_distance import surface_direction # noqa

xrspatial/sky_view_factor.py

Lines changed: 253 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,253 @@
1+
"""Sky-view factor (SVF) for digital elevation models.
2+
3+
SVF measures the fraction of the sky hemisphere visible from each cell,
4+
ranging from 0 (fully obstructed) to 1 (flat open terrain). For each
5+
cell, rays are cast at *n_directions* evenly spaced azimuths out to
6+
*max_radius* cells. Along each ray the maximum elevation angle to the
7+
horizon is recorded and SVF is computed as:
8+
9+
SVF(i,j) = 1 - mean(sin(max_horizon_angle(d)) for d in directions)
10+
11+
References
12+
----------
13+
Zakek, K., Ostir, K., Kokalj, Z. (2011). Sky-View Factor as a Relief
14+
Visualization Technique. Remote Sensing, 3(2), 398-415.
15+
16+
Yokoyama, R., Sidle, R.C., Noguchi, S. (2002). Application of
17+
Topographic Shading to Terrain Visualization. International Journal
18+
of Geographical Information Science, 16(5), 489-502.
19+
"""
20+
from __future__ import annotations
21+
22+
from functools import partial
23+
from math import atan2 as _atan2, cos as _cos, pi as _pi, sin as _sin, sqrt as _sqrt
24+
from typing import Optional
25+
26+
import numpy as np
27+
import xarray as xr
28+
from numba import cuda
29+
30+
try:
31+
import cupy
32+
except ImportError:
33+
class cupy:
34+
ndarray = False
35+
36+
try:
37+
import dask.array as da
38+
except ImportError:
39+
da = None
40+
41+
from xrspatial.dataset_support import supports_dataset
42+
from xrspatial.utils import (
43+
ArrayTypeFunctionMapping,
44+
_boundary_to_dask,
45+
_validate_raster,
46+
_validate_scalar,
47+
cuda_args,
48+
ngjit,
49+
)
50+
51+
52+
# ---------------------------------------------------------------------------
53+
# CPU kernel
54+
# ---------------------------------------------------------------------------
55+
56+
@ngjit
57+
def _svf_cpu(data, max_radius, n_directions):
58+
"""Compute SVF over an entire 2-D array on the CPU."""
59+
rows, cols = data.shape
60+
out = np.empty((rows, cols), dtype=np.float64)
61+
out[:] = np.nan
62+
63+
for y in range(rows):
64+
for x in range(cols):
65+
center = data[y, x]
66+
if center != center: # NaN check
67+
continue
68+
69+
svf_sum = 0.0
70+
for d in range(n_directions):
71+
angle = 2.0 * _pi * d / n_directions
72+
dx = _cos(angle)
73+
dy = _sin(angle)
74+
75+
max_elev_angle = 0.0
76+
for r in range(1, max_radius + 1):
77+
sx = x + int(round(r * dx))
78+
sy = y + int(round(r * dy))
79+
if sx < 0 or sx >= cols or sy < 0 or sy >= rows:
80+
break
81+
elev = data[sy, sx]
82+
if elev != elev: # NaN
83+
break
84+
dz = elev - center
85+
dist = float(r)
86+
elev_angle = _atan2(dz, dist)
87+
if elev_angle > max_elev_angle:
88+
max_elev_angle = elev_angle
89+
90+
svf_sum += _sin(max_elev_angle)
91+
92+
out[y, x] = 1.0 - svf_sum / n_directions
93+
return out
94+
95+
96+
# ---------------------------------------------------------------------------
97+
# GPU kernels
98+
# ---------------------------------------------------------------------------
99+
100+
@cuda.jit
101+
def _svf_gpu(data, out, max_radius, n_directions):
102+
"""CUDA global kernel: one thread per cell."""
103+
y, x = cuda.grid(2)
104+
rows, cols = data.shape
105+
106+
if y >= rows or x >= cols:
107+
return
108+
109+
center = data[y, x]
110+
if center != center: # NaN
111+
return
112+
113+
svf_sum = 0.0
114+
pi2 = 2.0 * 3.141592653589793
115+
116+
for d in range(n_directions):
117+
angle = pi2 * d / n_directions
118+
dx = _cos(angle)
119+
dy = _sin(angle)
120+
121+
max_elev_angle = 0.0
122+
for r in range(1, max_radius + 1):
123+
sx = x + int(round(r * dx))
124+
sy = y + int(round(r * dy))
125+
if sx < 0 or sx >= cols or sy < 0 or sy >= rows:
126+
break
127+
elev = data[sy, sx]
128+
if elev != elev: # NaN
129+
break
130+
dz = elev - center
131+
dist = float(r)
132+
elev_angle = _atan2(dz, dist)
133+
if elev_angle > max_elev_angle:
134+
max_elev_angle = elev_angle
135+
136+
svf_sum += _sin(max_elev_angle)
137+
138+
out[y, x] = 1.0 - svf_sum / n_directions
139+
140+
141+
# ---------------------------------------------------------------------------
142+
# Backend wrappers
143+
# ---------------------------------------------------------------------------
144+
145+
def _run_numpy(data, max_radius, n_directions):
146+
data = data.astype(np.float64)
147+
return _svf_cpu(data, max_radius, n_directions)
148+
149+
150+
def _run_cupy(data, max_radius, n_directions):
151+
data = data.astype(cupy.float64)
152+
out = cupy.full(data.shape, cupy.nan, dtype=cupy.float64)
153+
griddim, blockdim = cuda_args(data.shape)
154+
_svf_gpu[griddim, blockdim](data, out, max_radius, n_directions)
155+
return out
156+
157+
158+
def _run_dask_numpy(data, max_radius, n_directions):
159+
data = data.astype(np.float64)
160+
_func = partial(_svf_cpu, max_radius=max_radius, n_directions=n_directions)
161+
out = data.map_overlap(
162+
_func,
163+
depth=(max_radius, max_radius),
164+
boundary=np.nan,
165+
meta=np.array(()),
166+
)
167+
return out
168+
169+
170+
def _run_dask_cupy(data, max_radius, n_directions):
171+
data = data.astype(cupy.float64)
172+
_func = partial(_run_cupy, max_radius=max_radius, n_directions=n_directions)
173+
out = data.map_overlap(
174+
_func,
175+
depth=(max_radius, max_radius),
176+
boundary=cupy.nan,
177+
meta=cupy.array(()),
178+
)
179+
return out
180+
181+
182+
# ---------------------------------------------------------------------------
183+
# Public API
184+
# ---------------------------------------------------------------------------
185+
186+
@supports_dataset
187+
def sky_view_factor(
188+
agg: xr.DataArray,
189+
max_radius: int = 10,
190+
n_directions: int = 16,
191+
name: Optional[str] = 'sky_view_factor',
192+
) -> xr.DataArray:
193+
"""Compute the sky-view factor for each cell of a DEM.
194+
195+
SVF measures the fraction of the sky hemisphere visible from each
196+
cell on a scale from 0 (fully obstructed) to 1 (flat open terrain).
197+
Rays are cast at *n_directions* evenly spaced azimuths out to
198+
*max_radius* cells, and the maximum elevation angle along each
199+
ray determines the horizon obstruction.
200+
201+
Parameters
202+
----------
203+
agg : xarray.DataArray or xr.Dataset
204+
2D NumPy, CuPy, NumPy-backed Dask, or CuPy-backed Dask
205+
xarray DataArray of elevation values.
206+
If a Dataset is passed, the operation is applied to each
207+
data variable independently.
208+
max_radius : int, default=10
209+
Maximum search distance in cells along each ray direction.
210+
Cells within *max_radius* of the raster edge will be NaN.
211+
n_directions : int, default=16
212+
Number of azimuth directions to sample, evenly spaced
213+
around 360 degrees.
214+
name : str, default='sky_view_factor'
215+
Name of the output DataArray.
216+
217+
Returns
218+
-------
219+
xarray.DataArray or xr.Dataset
220+
2D array of SVF values in [0, 1] with the same shape, coords,
221+
dims, and attrs as the input. Edge cells use truncated rays
222+
that stop at the raster boundary.
223+
224+
References
225+
----------
226+
Zakek, K., Ostir, K., Kokalj, Z. (2011). Sky-View Factor as a
227+
Relief Visualization Technique. Remote Sensing, 3(2), 398-415.
228+
229+
Examples
230+
--------
231+
>>> from xrspatial import sky_view_factor
232+
>>> svf = sky_view_factor(dem, max_radius=100, n_directions=16)
233+
"""
234+
_validate_raster(agg, func_name='sky_view_factor', name='agg')
235+
_validate_scalar(max_radius, func_name='sky_view_factor',
236+
name='max_radius', dtype=int, min_val=1)
237+
_validate_scalar(n_directions, func_name='sky_view_factor',
238+
name='n_directions', dtype=int, min_val=1)
239+
240+
mapper = ArrayTypeFunctionMapping(
241+
numpy_func=_run_numpy,
242+
cupy_func=_run_cupy,
243+
dask_func=_run_dask_numpy,
244+
dask_cupy_func=_run_dask_cupy,
245+
)
246+
out = mapper(agg)(agg.data, max_radius, n_directions)
247+
return xr.DataArray(
248+
out,
249+
name=name,
250+
coords=agg.coords,
251+
dims=agg.dims,
252+
attrs=agg.attrs,
253+
)

0 commit comments

Comments
 (0)