diff --git a/.github/workflows/ci.yaml b/.github/workflows/ci.yaml index 6ae162b..c1b6307 100644 --- a/.github/workflows/ci.yaml +++ b/.github/workflows/ci.yaml @@ -19,7 +19,7 @@ jobs: debug-build: runs-on: ubuntu-latest env: - snapshot_version: v1.3.0 + snapshot_version: v1.4.0 steps: - uses: actions/checkout@v4 with: diff --git a/include/topotoolbox.h b/include/topotoolbox.h index ef110e1..6a1302a 100644 --- a/include/topotoolbox.h +++ b/include/topotoolbox.h @@ -876,6 +876,32 @@ void streamquad_trapz_f64(double *integral, double *integrand, ptrdiff_t *source, ptrdiff_t *target, float *weight, ptrdiff_t edge_count); +/** + @brief The anisotropic coefficient of variation (ACV) describes the general + geometry of the local land surface and can be used to destinguish elongated + from oval land forms. + + @param[out] output: The computed anisotropic of variation (ACV) + @parblock + A pointer to a `float` array of size `dims[0]` x `dims[1]` + @endparblock + + @param[in] dem: The input digital elevation model + @parblock + A pointer to a `float` array of size `dims[0]` x `dims[1]` + @endparblock + + @param[in] use_mp: If 1 then multiprocessing will be used to compute result + + @param[in] dims: The horizontal dimensions of the arrays + @parblock + A pointer to a `ptrdiff_t` array of size 2 + @endparblock + +*/ +TOPOTOOLBOX_API +void acv(float *output, float *dem, int use_mp, ptrdiff_t dims[2]); + /* Graphflood */ @@ -1028,4 +1054,5 @@ void graphflood_full(GF_FLOAT *Z, GF_FLOAT *hw, uint8_t *BCs, GF_FLOAT *Precipitations, GF_FLOAT *manning, GF_UINT *dim, GF_FLOAT dt, GF_FLOAT dx, bool SFD, bool D8, GF_UINT N_iterations, GF_FLOAT step); + #endif // TOPOTOOLBOX_H diff --git a/src/CMakeLists.txt b/src/CMakeLists.txt index 9ae69ad..e9b7c1d 100644 --- a/src/CMakeLists.txt +++ b/src/CMakeLists.txt @@ -7,6 +7,7 @@ add_library(topotoolbox identifyflats.c gwdt.c gradient8.c + acv.c morphology/reconstruct.c priority_queue.c excesstopography.c diff --git a/src/acv.c b/src/acv.c new file mode 100644 index 0000000..c80aa9e --- /dev/null +++ b/src/acv.c @@ -0,0 +1,265 @@ +#define TOPOTOOLBOX_BUILD + +#include +#include +#include + +#if TOPOTOOLBOX_OPENMP_VERSION > 0 +#include +#endif + +#include "topotoolbox.h" + +/* +The anisotropic coefficient of variation (ACV) describes the general +geometry of the local land surface and can be used to destinguish elongated +from oval land forms. + +output: Empty (or filled with zeros), same shape as dem +dem: DEM matrix, passed from python in order='C' +use_mp: If 1 then multiprocessing will be used to compute result +dims[2]: dimensions of DEM + +References: + Olaya, V. 2009: Basic land-surface parameters. In: Geomorphometry. + Concepts, Software, Applications, Hengl, T. & Reuter, H. I. (Eds.), + Elsevier, 33, 141-169. + +Author: + Theophil Bringezu (theophil.bringezu[at]uni-potsdam.de) + +Original MATLAB version by: + Wolfgang Schwanghart (w.schwanghart[at]geo.uni-potsdam.de) + +*/ + +TOPOTOOLBOX_API +void acv(float *output, float *dem, int use_mp, ptrdiff_t dims[2]) { + /* + filter_1 : + {{ 1, 0, 1, 0, 1}, + { 0, 0, 0, 0, 0}, + { 1, 0, 0, 0, -1}, + { 0, 0, 0, 0, 0}, + {-1, 0, -1, 0, -1}} + */ + float filter_1[25] = {1, 0, 1, 0, -1, 0, 0, 0, 0, 0, 1, 0, 0, + 0, -1, 0, 0, 0, 0, 0, 1, 0, -1, 0, -1}; + /* filter_2 + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {1, 0, 0, 0, -1}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}}, + + {{1, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, -1}}, + + {{0, 0, -1, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, -1, 0, 0}}, + + {{0, 0, 0, 0, -1}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {0, 0, 0, 0, 0}, + {1, 0, 0, 0, 0}}} + */ + float filter_2[4][25] = {{ + 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, + 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, + }, + {1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, + 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1}, + { + 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, + 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, + }, + {0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, + 0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0}}; + /* + filter_3 + {{{ 0, 0, 0}, + { 1, 0, -1}, + { 0, 0, 0}}, + + {{ 1, 0, 0}, + { 0, 0, 0}, + { 0, 0, -1}}, + + {{ 0, 1, 0}, + { 0, 0, 0}, + { 0, -1, 0}}, + + {{ 0, 0, 1}, + { 0, 0, 0}, + {-1, 0, 0}}}; + */ + float filter_3[4][9] = {{0, 1, 0, 0, 0, 0, 0, -1, 0}, + {1, 0, 0, 0, 0, 0, 0, 0, -1}, + {0, 0, 0, 1, 0, -1, 0, 0, 0}, + {0, 0, -1, 0, 0, 0, 1, 0, 0}}; + + /*{{21, 15, 10, 17, 23}, + {13, 5, 2, 7, 19}, + { 9, 1, 0, 4, 12}, + {14, 6, 3, 8, 20}, + {22, 16, 11, 18, 24}}; + search_order = {{row_shift, col_shift}, ...} + */ + ptrdiff_t search_order[25][2] = { + {0, 0}, {0, -1}, {-1, 0}, {1, 0}, {0, 1}, // 0 to 4 + {-1, -1}, {1, -1}, {-1, 1}, {1, 1}, // 5 to 8 + {0, -2}, {-2, 0}, {2, 0}, {0, 2}, // 9 to 12 + {-1, -2}, {1, -2}, {-2, -1}, {2, -1}, // 13 to 16 + {-2, 1}, {2, 1}, {-1, 2}, {1, 2}, // 17 to 20 + {-2, -2}, {2, -2}, {-2, 2}, {2, 2} // 21 to 24 + }; + + // ACV: +#if TOPOTOOLBOX_OPENMP_VERSION < 30 + ptrdiff_t col; +#pragma omp parallel for if (use_mp) + for (col = 0; col < dims[1]; col++) { + for (ptrdiff_t row = 0; row < dims[0]; row++) { +#else + ptrdiff_t col, row; +#pragma omp parallel for collapse(2) if (use_mp) + for (col = 0; col < dims[1]; col++) { + for (row = 0; row < dims[0]; row++) { +#endif + + // Catch NaN cells and skip them + ptrdiff_t index = col * dims[0] + row; + if (isnan(dem[index])) continue; + + float sum = 0.0f; + float dz_avg = 0.0f; + float anisotropic_cov = 0.0f; + + // filter_1 + for (ptrdiff_t k_col = -2; k_col <= 2; k_col++) { + for (ptrdiff_t k_row = -2; k_row <= 2; k_row++) { + // if filter cell is zero skip this filter cell (sum += 0.0f) + ptrdiff_t k_index = (k_col + 2) * 5 + (k_row + 2); + if (filter_1[k_index] == 0.0f) continue; + + ptrdiff_t true_row, true_col, true_index, search_pos; + int out_of_bounds; + // If out of bounds or isnan find nearest replacement value using + // Euclidean distance transform. (search_order[25][2]) + search_pos = 0; + while (true) { + true_col = col + k_col + search_order[search_pos][1]; + true_row = row + k_row + search_order[search_pos][0]; + out_of_bounds = (true_row < 0 || true_row >= dims[0] || + true_col < 0 || true_col >= dims[1]); + if (out_of_bounds) { + search_pos++; + continue; + } + true_index = true_col * dims[0] + true_row; + if (isnan(dem[true_index])) { + search_pos++; + continue; + } + // valid position found. While(true) loop will terminate because + // cell at `index` is valid. + break; + } + sum += filter_1[k_index] * dem[true_index]; + } + } + // dz_AVG = conv2(dem,k,'valid')/4; + dz_avg = sum / 4.0f; + + // filter_2 + for (ptrdiff_t n = 0; n < 4; n++) { + sum = 0.0f; + for (ptrdiff_t k_col = -2; k_col <= 2; k_col++) { + for (ptrdiff_t k_row = -2; k_row <= 2; k_row++) { + // if filter cell is zero skip this filter cell + ptrdiff_t k_index = (k_col + 2) * 5 + (k_row + 2); + if (filter_2[n][k_index] == 0.0f) continue; + + ptrdiff_t true_row, true_col, true_index, search_pos; + int out_of_bounds; + // If out of bounds or isnan find nearest replacement value using + // Euclidean distance transform. (search_order[25][2]) + search_pos = 0; + while (true) { + true_col = col + k_col + search_order[search_pos][1]; + true_row = row + k_row + search_order[search_pos][0]; + out_of_bounds = (true_row < 0 || true_row >= dims[0] || + true_col < 0 || true_col >= dims[1]); + if (out_of_bounds) { + search_pos++; + continue; + } + true_index = true_col * dims[0] + true_row; + if (isnan(dem[true_index])) { + search_pos++; + continue; + } + // valid position found. While(true) loop will terminate because + // cell at `index` is valid. + break; + } + sum += filter_2[n][k_index] * dem[true_index]; + } + } + // ACV = ACV + (conv2(dem,F{r},'valid') - dz_AVG).^2; + anisotropic_cov += powf(sum - dz_avg, 2.0f); + } + + // filter_3 + for (ptrdiff_t n = 0; n < 4; n++) { + sum = 0.0f; + for (ptrdiff_t k_col = -1; k_col <= 1; k_col++) { + for (ptrdiff_t k_row = -1; k_row <= 1; k_row++) { + // if filter cell is zero skip this filter cell + ptrdiff_t k_index = (k_col + 1) * 3 + (k_row + 1); + if (filter_3[n][k_index] == 0.0f) continue; + + ptrdiff_t true_row, true_col, true_index, search_pos; + int out_of_bounds; + // If out of bounds or isnan find nearest replacement value using + // Euclidean distance transform. (search_order[25][2]) + search_pos = 0; + while (true) { + true_col = col + k_col + search_order[search_pos][1]; + true_row = row + k_row + search_order[search_pos][0]; + out_of_bounds = (true_row < 0 || true_row >= dims[0] || + true_col < 0 || true_col >= dims[1]); + if (out_of_bounds) { + search_pos++; + continue; + } + true_index = true_col * dims[0] + true_row; + if (isnan(dem[true_index])) { + search_pos++; + continue; + } + // valid position found. While(true) loop will terminate because + // cell at `index` is valid. + break; + } + sum += filter_3[n][k_index] * dem[true_index]; + } + } + // ACV = ACV + (conv2(dem,F{r},'valid') - dz_AVG).^2; + anisotropic_cov += powf(sum - dz_avg, 2.0f); + } + // dz_AVG = max(abs(dz_AVG),0.001); + dz_avg = fmaxf(0.001f, fabsf(dz_avg)); + + // C = log(1 + sqrt(ACV./8)./dz_AVG); + output[index] = logf(1.0f + sqrtf(anisotropic_cov / 8.0f) / dz_avg); + } + } +} diff --git a/test/random_dem.cpp b/test/random_dem.cpp index b258853..f243b4c 100644 --- a/test/random_dem.cpp +++ b/test/random_dem.cpp @@ -353,6 +353,10 @@ int32_t test_gradient8_mp(float *gradient, float *gradient_mp, } return 0; } +/* + +*/ +int32_t test_acv() { return 0; } /* Flow direction should point downstream or across flats */ diff --git a/test/snapshot.cpp b/test/snapshot.cpp index 0dfb75d..8cb66e9 100644 --- a/test/snapshot.cpp +++ b/test/snapshot.cpp @@ -70,12 +70,14 @@ struct SnapshotData { std::vector filled_dem; std::vector flats; std::vector sills; + std::vector acv; // Output arrays std::vector test_dem; std::vector test_filled_dem; std::vector test_flats; std::vector test_sills; + std::vector test_acv; // Intermediate arrays std::vector bc; @@ -111,6 +113,12 @@ struct SnapshotData { assert(dims[0] == dims_check[0] && dims[1] == dims_check[1]); } + if (exists(snapshot_path / "acv.tif")) { + load_data_from_file(snapshot_path / "acv.tif", acv, + dims_check); + assert(dims[0] == dims_check[0] && dims[1] == dims_check[1]); + } + // Allocate and resize output and intermediate arrays if (dem.size() > 0) { if (filled_dem.size() > 0) { @@ -121,10 +129,13 @@ struct SnapshotData { if (flats.size() > 0 && sills.size() > 0) { test_flats.resize(dims[0] * dims[1]); } + if (acv.size() > 0) { + test_acv.resize(dims[0] * dims[1]); + } } } - int test_fillsinks() { + int run_fillsinks() { // Initialize bcs for (ptrdiff_t j = 0; j < dims[1]; j++) { for (ptrdiff_t i = 0; i < dims[0]; i++) { @@ -161,7 +172,7 @@ struct SnapshotData { return 0; } - int test_identifyflats() { + int run_identifyflats() { // identifyflats // // Use the snapshot filled DEM rather than the generated one in @@ -188,11 +199,28 @@ struct SnapshotData { return 0; } + int run_acv() { + // acv + tt::acv(test_acv.data(), dem.data(), 0, dims.data()); + + for (ptrdiff_t j = 0; j < dims[1]; j++) { + for (ptrdiff_t i = 0; i < dims[0]; i++) { + if (test_acv[j * dims[0] + i] != acv[j * dims[0] + i]) { + write_data_to_file( + path / "test_acv.tif", path / "acv.tif", test_acv, dims); + return -1; + } + } + } + + return 0; + } + int runtests() { int result = 0; // fillsinks if (test_filled_dem.size() > 0) { - if (test_fillsinks() < 0) { + if (run_fillsinks() < 0) { result = -1; std::cout << "[FAILURE] (fillsinks) " << path << std::endl; @@ -202,7 +230,7 @@ struct SnapshotData { } if (flats.size() > 0 && sills.size() > 0) { - if (test_identifyflats() < 0) { + if (run_identifyflats() < 0) { result = -1; std::cout << "[FAILURE] (identifyflats) " << path << std::endl; @@ -211,6 +239,16 @@ struct SnapshotData { } } + if (acv.size() > 0) { + if (run_acv() < 0) { + result = -1; + + std::cout << "[FAILURE] (acv) " << path << std::endl; + } else { + std::cout << "[SUCCESS] (acv) " << path << std::endl; + } + } + return result; } }; @@ -242,4 +280,4 @@ int main(int argc, char* argv[]) { } return result; -} +} \ No newline at end of file