|
4 | 4 |
|
5 | 5 | import numpy as np
|
6 | 6 | import xarray as xr
|
| 7 | +from scipy.spatial import KDTree |
7 | 8 |
|
8 | 9 |
|
9 | 10 | def ps_selection(
|
@@ -126,6 +127,126 @@ def ps_selection(
|
126 | 127 | return stm_masked
|
127 | 128 |
|
128 | 129 |
|
| 130 | +def network_stm_selection( |
| 131 | + stm: xr.Dataset, |
| 132 | + min_dist: int | float, |
| 133 | + include_index: list[int] = None, |
| 134 | + sortby_var: str = "pnt_nmad", |
| 135 | + crs: int | str = "radar", |
| 136 | + x_var: str = "azimuth", |
| 137 | + y_var: str = "range", |
| 138 | + azimuth_spacing: float = None, |
| 139 | + range_spacing: float = None, |
| 140 | +): |
| 141 | + """Select a Space-Time Matrix (STM) from a candidate STM for network processing. |
| 142 | +
|
| 143 | + The selection is based on two criteria: |
| 144 | + 1. A minimum distance between selected points. |
| 145 | + 2. A sorting metric to select better points. |
| 146 | +
|
| 147 | + The candidate STM will be sorted by the sorting metric. |
| 148 | + The selection will be performed iteratively, starting from the best point. |
| 149 | + In each iteration, the best point will be selected, and points within the minimum distance will be removed. |
| 150 | + The process will continue until no points are left in the candidate STM. |
| 151 | +
|
| 152 | + Parameters |
| 153 | + ---------- |
| 154 | + stm : xr.Dataset |
| 155 | + candidate Space-Time Matrix (STM). |
| 156 | + min_dist : int | float |
| 157 | + Minimum distance between selected points. |
| 158 | + include_index : list[int], optional |
| 159 | + Index of points in the candidate STM that must be included in the selection, by default None |
| 160 | + sortby_var : str, optional |
| 161 | + Sorting metric for selecting points, by default "pnt_nmad" |
| 162 | + crs : int | str, optional |
| 163 | + EPSG code of Coordinate Reference System of `x_var` and `y_var`, by default "radar". |
| 164 | + If crs is "radar", the distance will be calculated based on radar coordinates, and |
| 165 | + azimuth_spacing and range_spacing must be provided. |
| 166 | + x_var : str, optional |
| 167 | + Data variable name for x coordinate, by default "azimuth" |
| 168 | + y_var : str, optional |
| 169 | + Data variable name for y coordinate, by default "range" |
| 170 | + azimuth_spacing : float, optional |
| 171 | + Azimuth spacing, by default None. Required if crs is "radar". |
| 172 | + range_spacing : float, optional |
| 173 | + Range spacing, by default None. Required if crs is "radar". |
| 174 | +
|
| 175 | + Returns |
| 176 | + ------- |
| 177 | + xr.Dataset |
| 178 | + Selected network Space-Time Matrix (STM). |
| 179 | +
|
| 180 | + Raises |
| 181 | + ------ |
| 182 | + ValueError |
| 183 | + Raised when `azimuth_spacing` or `range_spacing` is not provided for radar coordinates. |
| 184 | + NotImplementedError |
| 185 | + Raised when an unsupported Coordinate Reference System is provided. |
| 186 | + """ |
| 187 | + match crs: |
| 188 | + case "radar": |
| 189 | + if (azimuth_spacing is None) or (range_spacing is None): |
| 190 | + raise ValueError("Azimuth and range spacing must be provided for radar coordinates.") |
| 191 | + case _: |
| 192 | + raise NotImplementedError |
| 193 | + |
| 194 | + # Get coordinates and sorting metric, load them into memory |
| 195 | + stm_select = None |
| 196 | + stm_remain = stm[[x_var, y_var, sortby_var]].compute() |
| 197 | + |
| 198 | + # Select the include_index if provided |
| 199 | + if include_index is not None: |
| 200 | + stm_select = stm_remain.isel(space=include_index) |
| 201 | + |
| 202 | + # Remove points within min_dist of the included points |
| 203 | + coords_include = np.column_stack( |
| 204 | + [stm_select["azimuth"].values * azimuth_spacing, stm_select["range"].values * range_spacing] |
| 205 | + ) |
| 206 | + coords_remain = np.column_stack( |
| 207 | + [stm_remain["azimuth"].values * azimuth_spacing, stm_remain["range"].values * range_spacing] |
| 208 | + ) |
| 209 | + idx_drop = _idx_within_distance(coords_include, coords_remain, min_dist) |
| 210 | + if idx_drop is not None: |
| 211 | + stm_remain = stm_remain.where(~(stm_remain["space"].isin(idx_drop)), drop=True) |
| 212 | + |
| 213 | + # Reorder the remaining points by the sorting metric |
| 214 | + stm_remain = stm_remain.sortby(sortby_var) |
| 215 | + |
| 216 | + # Build a list of the index of selected points |
| 217 | + if stm_select is None: |
| 218 | + space_idx_sel = [] |
| 219 | + else: |
| 220 | + space_idx_sel = stm_select["space"].values.tolist() |
| 221 | + |
| 222 | + while stm_remain.sizes["space"] > 0: |
| 223 | + # Select one point with best sorting metric |
| 224 | + stm_now = stm_remain.isel(space=0) |
| 225 | + |
| 226 | + # Append the selected point index |
| 227 | + space_idx_sel.append(stm_now["space"].values.tolist()) |
| 228 | + |
| 229 | + # Remove the selected point from the remaining points |
| 230 | + stm_remain = stm_remain.isel(space=slice(1, None)).copy() |
| 231 | + |
| 232 | + # Remove points in stm_remain within min_dist of stm_now |
| 233 | + coords_remain = np.column_stack( |
| 234 | + [stm_remain["azimuth"].values * azimuth_spacing, stm_remain["range"].values * range_spacing] |
| 235 | + ) |
| 236 | + coords_stmnow = np.column_stack( |
| 237 | + [stm_now["azimuth"].values * azimuth_spacing, stm_now["range"].values * range_spacing] |
| 238 | + ) |
| 239 | + idx_drop = _idx_within_distance(coords_stmnow, coords_remain, min_dist) |
| 240 | + if idx_drop is not None: |
| 241 | + stm_drop = stm_remain.isel(space=idx_drop) |
| 242 | + stm_remain = stm_remain.where(~(stm_remain["space"].isin(stm_drop["space"])), drop=True) |
| 243 | + |
| 244 | + # Get the selected points by space index from the original stm |
| 245 | + stm_out = stm.sel(space=space_idx_sel) |
| 246 | + |
| 247 | + return stm_out |
| 248 | + |
| 249 | + |
129 | 250 | def _nad_block(amp: xr.DataArray) -> xr.DataArray:
|
130 | 251 | """Compute Normalized Amplitude Dispersion (NAD) for a block of amplitude data.
|
131 | 252 |
|
@@ -170,3 +291,30 @@ def _nmad_block(amp: xr.DataArray) -> xr.DataArray:
|
170 | 291 | nmad = mad / (median_amplitude + np.finfo(amp.dtype).eps) # Normalized Median Absolute Deviation
|
171 | 292 |
|
172 | 293 | return nmad
|
| 294 | + |
| 295 | + |
| 296 | +def _idx_within_distance(coords_ref, coords_others, min_dist): |
| 297 | + """Get the index of points in coords_others that are within min_dist of coords_ref. |
| 298 | +
|
| 299 | + Parameters |
| 300 | + ---------- |
| 301 | + coords_ref : np.ndarray |
| 302 | + Coordinates of reference points. Shape (n, 2). |
| 303 | + coords_others : np.ndarray |
| 304 | + Coordinates of other points. Shape (m, 2). |
| 305 | + min_dist : int, float |
| 306 | + distance threshold. |
| 307 | +
|
| 308 | + Returns |
| 309 | + ------- |
| 310 | + np.ndarray |
| 311 | + Index of points in coords_others that are within `min_dist` of `coords_ref`. |
| 312 | + """ |
| 313 | + kd_ref = KDTree(coords_ref) |
| 314 | + kd_others = KDTree(coords_others) |
| 315 | + sdm = kd_ref.sparse_distance_matrix(kd_others, min_dist) |
| 316 | + if len(sdm) > 0: |
| 317 | + idx = np.array(list(sdm.keys()))[:, 1] |
| 318 | + return idx |
| 319 | + else: |
| 320 | + return None |
0 commit comments