Skip to content

Commit

Permalink
Use fillsinks_hybrid (#99)
Browse files Browse the repository at this point in the history
This adds a `hybrid` option to `GridObject.fillsinks` and `FlowObject` that dispatches to the hybrid morphological reconstruction. My experiments suggest that this is almost always faster than the sequential reconstruction, so `hybrid=True` is the new default, which chooses `fillsinks_hybrid`.

A test is added that ensures that the results from the two algorithms are the same.
  • Loading branch information
wkearn authored Dec 6, 2024
1 parent 4674305 commit b1505d9
Show file tree
Hide file tree
Showing 4 changed files with 46 additions and 6 deletions.
17 changes: 17 additions & 0 deletions src/lib/grid.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -32,6 +32,22 @@ void wrap_fillsinks(py::array_t<float> output, py::array_t<float> dem,
fillsinks(output_ptr, dem_ptr, bc_ptr, dims_ptr);
}

void wrap_fillsinks_hybrid(py::array_t<float> output,
py::array_t<ptrdiff_t> queue,
py::array_t<float> dem,
py::array_t<uint8_t> bc,
std::tuple<ptrdiff_t, ptrdiff_t> dims){

float *output_ptr = output.mutable_data();
float *dem_ptr = dem.mutable_data();
ptrdiff_t *queue_ptr = queue.mutable_data();
uint8_t *bc_ptr = bc.mutable_data();

std::array<ptrdiff_t, 2> dims_array = {std::get<0>(dims), std::get<1>(dims)};
ptrdiff_t *dims_ptr = dims_array.data();
fillsinks_hybrid(output_ptr, queue_ptr, dem_ptr, bc_ptr, dims_ptr);
}

// wrap_identifyflats:
// Parameters:
// output: A NumPy array to store the output, where flats, sill amd presills will be marked.
Expand Down Expand Up @@ -227,6 +243,7 @@ void wrap_gradient8(

PYBIND11_MODULE(_grid, m) {
m.def("fillsinks", &wrap_fillsinks);
m.def("fillsinks_hybrid", &wrap_fillsinks_hybrid);
m.def("identifyflats", &wrap_identifyflats);
m.def("excesstopography_fsm2d", &wrap_excesstopography_fsm2d);
m.def("excesstopography_fmm2d", &wrap_excesstopography_fmm2d);
Expand Down
15 changes: 12 additions & 3 deletions src/topotoolbox/flow_object.py
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,8 @@ class FlowObject():
"""

def __init__(self, grid: GridObject,
bc: np.ndarray | GridObject | None = None):
bc: np.ndarray | GridObject | None = None,
hybrid : bool = True):
"""The constructor for the FlowObject. Takes a GridObject as input,
computes flow direction information and saves them as an FlowObject.
Expand All @@ -30,6 +31,10 @@ def __init__(self, grid: GridObject,
indicate pixels that should be fixed to their values in the
original DEM and values of 0 indicate pixels that should be
filled.
hybrid: bool, optional
Should hybrid reconstruction algorithm be used to fill
sinks? Defaults to True. Hybrid reconstruction is faster
but requires additional memory be allocated for a queue.
Notes
-----
Expand Down Expand Up @@ -58,7 +63,11 @@ def __init__(self, grid: GridObject,
if isinstance(bc, GridObject):
bc = bc.z

_grid.fillsinks(filled_dem, dem, bc, dims)
queue = np.zeros_like(dem, dtype=np.int64, order='F')
if hybrid:
_grid.fillsinks_hybrid(filled_dem, queue, dem, bc, dims)
else:
_grid.fillsinks(filled_dem, dem, bc, dims)

if restore_nans:
dem[nans] = np.nan
Expand All @@ -73,7 +82,7 @@ def __init__(self, grid: GridObject,

dist = np.zeros_like(flats, dtype=np.float32, order='F')
prev = conncomps # prev: dtype=np.int64
heap = np.zeros_like(flats, dtype=np.int64, order='F')
heap = queue # heap: dtype=np.int64
back = np.zeros_like(flats, dtype=np.int64, order='F')
_grid.gwdt(dist, prev, costs, flats, heap, back, dims)

Expand Down
13 changes: 11 additions & 2 deletions src/topotoolbox/grid_object.py
Original file line number Diff line number Diff line change
Expand Up @@ -105,7 +105,8 @@ def reproject(self,
return dst

def fillsinks(self,
bc: 'np.ndarray | GridObject | None' = None) -> 'GridObject':
bc: 'np.ndarray | GridObject | None' = None,
hybrid: bool = True) -> 'GridObject':
"""Fill sinks in the digital elevation model (DEM).
Parameters
Expand All @@ -116,6 +117,10 @@ def fillsinks(self,
indicate pixels that should be fixed to their values in the
original DEM and values of 0 indicate pixels that should be
filled.
hybrid: bool, optional
Should hybrid reconstruction algorithm be used? Defaults to True. Hybrid
reconstruction is faster but requires additional memory be allocated
for a queue.
Returns
-------
Expand Down Expand Up @@ -146,7 +151,11 @@ def fillsinks(self,
if isinstance(bc, GridObject):
bc = bc.z

_grid.fillsinks(output, dem, bc, self.shape)
if hybrid:
queue = np.zeros_like(dem, dtype=np.int64)
_grid.fillsinks_hybrid(output, queue, dem, bc, self.shape)
else:
_grid.fillsinks(output, dem, bc, self.shape)

if restore_nans:
dem[nans] = np.nan
Expand Down
7 changes: 6 additions & 1 deletion tests/test_grid_object.py
Original file line number Diff line number Diff line change
Expand Up @@ -25,11 +25,16 @@ def test_fillsinks(square_dem, wide_dem, tall_dem):
# since grid is a fixture, it has to be assigned/called first
dem = grid
original_dem = dem.z.copy()
filled_dem = dem.fillsinks()
filled_dem = dem.fillsinks(hybrid=False)
filled_dem_hybrid = dem.fillsinks(hybrid=True)

# Ensure that DEM has not been modified by fillsinks
assert np.all(original_dem == dem.z)

# The sequential and hybrid reconstruction algorithms should
# produce the same result
assert np.all(filled_dem.z == filled_dem_hybrid.z)

# Loop over all cells of the DEM
for i in range(dem.shape[0]):
for j in range(dem.shape[1]):
Expand Down

0 comments on commit b1505d9

Please sign in to comment.