|
| 1 | +import numpy as np |
| 2 | + |
| 3 | + |
| 4 | +class SpatialHash: |
| 5 | + """Custom data structure that is used for performing grid searches using Spatial Hashing. This class constructs an overlying |
| 6 | + uniformly spaced rectilinear grid, called the "hash grid" on top parcels.xgrid.XGrid. It is particularly useful for grid searching |
| 7 | + on curvilinear grids. Faces in the Xgrid are related to the cells in the hash grid by determining the hash cells the bounding box |
| 8 | + of the unstructured face cells overlap with. |
| 9 | +
|
| 10 | + Parameters |
| 11 | + ---------- |
| 12 | + grid : parcels.xgrid.XGrid |
| 13 | + Source grid used to construct the hash grid and hash table |
| 14 | + reconstruct : bool, default=False |
| 15 | + If true, reconstructs the spatial hash |
| 16 | +
|
| 17 | + Note |
| 18 | + ---- |
| 19 | + Does not currently support queries on periodic elements. |
| 20 | + """ |
| 21 | + |
| 22 | + def __init__( |
| 23 | + self, |
| 24 | + grid, |
| 25 | + reconstruct=False, |
| 26 | + ): |
| 27 | + # TODO : Enforce grid to be an instance of parcels.xgrid.XGrid |
| 28 | + # Currently, this is not done due to circular import with parcels.xgrid |
| 29 | + |
| 30 | + self._source_grid = grid |
| 31 | + self.reconstruct = reconstruct |
| 32 | + |
| 33 | + # Hash grid size |
| 34 | + self._dh = self._hash_cell_size() |
| 35 | + |
| 36 | + # Lower left corner of the hash grid |
| 37 | + lon_min = self._source_grid.lon.min() |
| 38 | + lat_min = self._source_grid.lat.min() |
| 39 | + lon_max = self._source_grid.lon.max() |
| 40 | + lat_max = self._source_grid.lat.max() |
| 41 | + |
| 42 | + # Get corner vertices of each face |
| 43 | + self._xbound = np.stack( |
| 44 | + ( |
| 45 | + self._source_grid.lon[:-1, :-1], |
| 46 | + self._source_grid.lon[:-1, 1:], |
| 47 | + self._source_grid.lon[1:, 1:], |
| 48 | + self._source_grid.lon[1:, :-1], |
| 49 | + ), |
| 50 | + axis=-1, |
| 51 | + ) |
| 52 | + self._ybound = np.stack( |
| 53 | + ( |
| 54 | + self._source_grid.lat[:-1, :-1], |
| 55 | + self._source_grid.lat[:-1, 1:], |
| 56 | + self._source_grid.lat[1:, 1:], |
| 57 | + self._source_grid.lat[1:, :-1], |
| 58 | + ), |
| 59 | + axis=-1, |
| 60 | + ) |
| 61 | + |
| 62 | + self._xmin = lon_min - self._dh |
| 63 | + self._ymin = lat_min - self._dh |
| 64 | + self._xmax = lon_max + self._dh |
| 65 | + self._ymax = lat_max + self._dh |
| 66 | + |
| 67 | + # Number of x points in the hash grid; used for |
| 68 | + # array flattening |
| 69 | + Lx = self._xmax - self._xmin |
| 70 | + Ly = self._ymax - self._ymin |
| 71 | + self._nx = int(np.ceil(Lx / self._dh)) |
| 72 | + self._ny = int(np.ceil(Ly / self._dh)) |
| 73 | + |
| 74 | + # Generate the mapping from the hash indices to unstructured grid elements |
| 75 | + self._face_hash_table = None |
| 76 | + self._face_hash_table = self._initialize_face_hash_table() |
| 77 | + |
| 78 | + def _hash_cell_size(self): |
| 79 | + """Computes the size of the hash cells from the source grid. |
| 80 | + The hash cell size is set to 1/2 of the square root of the median cell area |
| 81 | + """ |
| 82 | + return np.sqrt(np.median(_planar_quad_area(self._source_grid.lat, self._source_grid.lon))) * 0.5 |
| 83 | + |
| 84 | + def _hash_index2d(self, coords): |
| 85 | + """Computes the 2-d hash index (i,j) for the location (x,y), where x and y is the same units |
| 86 | + as the source grid coordinates |
| 87 | + """ |
| 88 | + # Wrap longitude to [-180, 180] |
| 89 | + if self._source_grid.mesh == "spherical": |
| 90 | + lon = (coords[:, 1] + 180.0) % (360.0) - 180.0 |
| 91 | + else: |
| 92 | + lon = coords[:, 1] |
| 93 | + i = ((lon - self._xmin) / self._dh).astype(np.int32) |
| 94 | + j = ((coords[:, 0] - self._ymin) / self._dh).astype(np.int32) |
| 95 | + return j, i |
| 96 | + |
| 97 | + def _hash_index(self, coords): |
| 98 | + """Computes the flattened hash index for the location (x,y), where x and y are given in spherical |
| 99 | + coordinates (in degrees). The single dimensioned hash index orders the flat index with all of the |
| 100 | + i-points first and then all the j-points. |
| 101 | + """ |
| 102 | + j, i = self._hash_index2d(coords) |
| 103 | + return i + self._nx * j |
| 104 | + |
| 105 | + def _grid_ji_for_eid(self, eid): |
| 106 | + """Returns the (i,j) grid coordinates for the given element id (eid)""" |
| 107 | + j = eid // (self._source_grid.xdim) |
| 108 | + i = eid - j * (self._source_grid.xdim) |
| 109 | + return j, i |
| 110 | + |
| 111 | + def _initialize_face_hash_table(self): |
| 112 | + """Create a mapping that relates unstructured grid faces to hash indices by determining |
| 113 | + which faces overlap with which hash cells |
| 114 | + """ |
| 115 | + if self._face_hash_table is None or self.reconstruct: |
| 116 | + index_to_face = [[] for i in range(self._nx * self._ny)] |
| 117 | + # Get the bounds of each curvilinear faces |
| 118 | + lat_bounds, lon_bounds = _curvilinear_grid_facebounds( |
| 119 | + self._source_grid.lat, |
| 120 | + self._source_grid.lon, |
| 121 | + ) |
| 122 | + coords = np.stack( |
| 123 | + ( |
| 124 | + lat_bounds[:, :, 0].flatten(), |
| 125 | + lon_bounds[:, :, 0].flatten(), |
| 126 | + ), |
| 127 | + axis=-1, |
| 128 | + ) |
| 129 | + yi1, xi1 = self._hash_index2d(coords) |
| 130 | + coords = np.stack( |
| 131 | + ( |
| 132 | + lat_bounds[:, :, 1].flatten(), |
| 133 | + lon_bounds[:, :, 1].flatten(), |
| 134 | + ), |
| 135 | + axis=-1, |
| 136 | + ) |
| 137 | + yi2, xi2 = self._hash_index2d(coords) |
| 138 | + nface = (self._source_grid.xdim) * (self._source_grid.ydim) |
| 139 | + for eid in range(nface): |
| 140 | + for j in range(yi1[eid], yi2[eid] + 1): |
| 141 | + if xi1[eid] <= xi2[eid]: |
| 142 | + # Normal case, no wrap |
| 143 | + for i in range(xi1[eid], xi2[eid] + 1): |
| 144 | + index_to_face[(i % self._nx) + self._nx * j].append(eid) |
| 145 | + else: |
| 146 | + # Wrap-around case |
| 147 | + for i in range(xi1[eid], self._nx): |
| 148 | + index_to_face[(i % self._nx) + self._nx * j].append(eid) |
| 149 | + for i in range(0, xi2[eid] + 1): |
| 150 | + index_to_face[(i % self._nx) + self._nx * j].append(eid) |
| 151 | + return index_to_face |
| 152 | + |
| 153 | + def query( |
| 154 | + self, |
| 155 | + coords, |
| 156 | + tol=1e-6, |
| 157 | + ): |
| 158 | + """Queries the hash table. |
| 159 | +
|
| 160 | + Parameters |
| 161 | + ---------- |
| 162 | + coords : array_like |
| 163 | + coordinate pairs in degrees (lat, lon) to query. |
| 164 | +
|
| 165 | +
|
| 166 | + Returns |
| 167 | + ------- |
| 168 | + faces : ndarray of shape (coords.shape[0]), dtype=np.int32 |
| 169 | + Face id's in the self._source_grid where each coords element is found. When a coords element is not found, the |
| 170 | + corresponding array entry in faces is set to -1. |
| 171 | + """ |
| 172 | + num_coords = coords.shape[0] |
| 173 | + |
| 174 | + # Preallocate results |
| 175 | + faces = np.full((num_coords, 2), -1, dtype=np.int32) |
| 176 | + |
| 177 | + # Get the list of candidate faces for each coordinate |
| 178 | + candidate_faces = [self._face_hash_table[pid] for pid in self._hash_index(coords)] |
| 179 | + |
| 180 | + for i, (coord, candidates) in enumerate(zip(coords, candidate_faces, strict=False)): |
| 181 | + for face_id in candidates: |
| 182 | + yi, xi = self._grid_ji_for_eid(face_id) |
| 183 | + nodes = np.stack( |
| 184 | + ( |
| 185 | + self._ybound[yi, xi, :], |
| 186 | + self._xbound[yi, xi, :], |
| 187 | + ), |
| 188 | + axis=-1, |
| 189 | + ) |
| 190 | + |
| 191 | + bcoord = np.asarray(_barycentric_coordinates(nodes, coord)) |
| 192 | + err = abs(np.dot(bcoord, nodes[:, 0]) - coord[0]) + abs(np.dot(bcoord, nodes[:, 1]) - coord[1]) |
| 193 | + if (bcoord >= 0).all() and err < tol: |
| 194 | + faces[i, :] = [yi, xi] |
| 195 | + break |
| 196 | + |
| 197 | + return faces |
| 198 | + |
| 199 | + |
| 200 | +def _triangle_area(A, B, C): |
| 201 | + """Compute the area of a triangle given by three points.""" |
| 202 | + d1 = B - A |
| 203 | + d2 = C - A |
| 204 | + d3 = np.cross(d1, d2) |
| 205 | + return 0.5 * np.linalg.norm(d3) |
| 206 | + |
| 207 | + |
| 208 | +def _barycentric_coordinates(nodes, point, min_area=1e-8): |
| 209 | + """ |
| 210 | + Compute the barycentric coordinates of a point P inside a convex polygon using area-based weights. |
| 211 | + So that this method generalizes to n-sided polygons, we use the Waschpress points as the generalized |
| 212 | + barycentric coordinates, which is only valid for convex polygons. |
| 213 | +
|
| 214 | + Parameters |
| 215 | + ---------- |
| 216 | + nodes : numpy.ndarray |
| 217 | + Spherical coordinates (lat,lon) of each corner node of a face |
| 218 | + point : numpy.ndarray |
| 219 | + Spherical coordinates (lat,lon) of the point |
| 220 | +
|
| 221 | + Returns |
| 222 | + ------- |
| 223 | + numpy.ndarray |
| 224 | + Barycentric coordinates corresponding to each vertex. |
| 225 | +
|
| 226 | + """ |
| 227 | + n = len(nodes) |
| 228 | + sum_wi = 0 |
| 229 | + w = [] |
| 230 | + |
| 231 | + for i in range(0, n): |
| 232 | + vim1 = nodes[i - 1] |
| 233 | + vi = nodes[i] |
| 234 | + vi1 = nodes[(i + 1) % n] |
| 235 | + a0 = _triangle_area(vim1, vi, vi1) |
| 236 | + a1 = max(_triangle_area(point, vim1, vi), min_area) |
| 237 | + a2 = max(_triangle_area(point, vi, vi1), min_area) |
| 238 | + sum_wi += a0 / (a1 * a2) |
| 239 | + w.append(a0 / (a1 * a2)) |
| 240 | + barycentric_coords = [w_i / sum_wi for w_i in w] |
| 241 | + |
| 242 | + return barycentric_coords |
| 243 | + |
| 244 | + |
| 245 | +def _planar_quad_area(lat, lon): |
| 246 | + """Computes the area of each quadrilateral face in a curvilinear grid. |
| 247 | + The lon and lat arrays are assumed to be 2D arrays of points with dimensions (n_y, n_x). |
| 248 | + The area is computed using the Shoelace formula. |
| 249 | + This method is only used during hashgrid construction to determine the size of the hash cells. |
| 250 | +
|
| 251 | + Parameters |
| 252 | + ---------- |
| 253 | + lon : np.ndarray |
| 254 | + 2D array of shape (n_y, n_x) containing the longitude of each corner node of the curvilinear grid. |
| 255 | + lat : np.ndarray |
| 256 | + 2D array of shape (n_y, n_x) containing the latitude of each corner node of the curvilinear grid. |
| 257 | +
|
| 258 | + Returns |
| 259 | + ------- |
| 260 | + area : np.ndarray |
| 261 | + 2D array of shape (n_y-1, n_x-1) containing the area of each quadrilateral face in the curvilinear grid. |
| 262 | + """ |
| 263 | + x0 = lon[:-1, :-1] |
| 264 | + x1 = lon[:-1, 1:] |
| 265 | + x2 = lon[1:, 1:] |
| 266 | + x3 = lon[1:, :-1] |
| 267 | + |
| 268 | + y0 = lat[:-1, :-1] |
| 269 | + y1 = lat[:-1, 1:] |
| 270 | + y2 = lat[1:, 1:] |
| 271 | + y3 = lat[1:, :-1] |
| 272 | + |
| 273 | + # Shoelace formula: 0.5 * |sum(x_i*y_{i+1} - x_{i+1}*y_i)| |
| 274 | + area = 0.5 * np.abs(x0 * y1 + x1 * y2 + x2 * y3 + x3 * y0 - y0 * x1 - y1 * x2 - y2 * x3 - y3 * x0) |
| 275 | + return area |
| 276 | + |
| 277 | + |
| 278 | +def _curvilinear_grid_facebounds(lat, lon): |
| 279 | + """Computes the bounds of each curvilinear face in the grid. |
| 280 | + The lon and lat arrays are assumed to be 2D arrays of points with dimensions (n_y, n_x). |
| 281 | + The bounds are for faces whose corner node vertices are defined by lat,lon. |
| 282 | + Face(yi,xi) is surrounding by points (yi,xi), (yi,xi+1), (yi+1,xi+1), (yi+1,xi). |
| 283 | + This method is only used during hashgrid construction to determine which curvilinear |
| 284 | + faces overlap with which hash cells. |
| 285 | +
|
| 286 | + Parameters |
| 287 | + ---------- |
| 288 | + lon : np.ndarray |
| 289 | + 2D array of shape (n_y, n_x) containing the longitude of each corner node of the curvilinear grid. |
| 290 | + lat : np.ndarray |
| 291 | + 2D array of shape (n_y, n_x) containing the latitude of each corner node of the curvilinear grid. |
| 292 | +
|
| 293 | + Returns |
| 294 | + ------- |
| 295 | + xbounds : np.ndarray |
| 296 | + Array of shape (n_y-1, n_x-1, 2) containing the bounds of each face in the x-direction. |
| 297 | + ybounds : np.ndarray |
| 298 | + Array of shape (n_y-1, n_x-1, 2) containing the bounds of each face in the y-direction. |
| 299 | + """ |
| 300 | + xf = np.stack((lon[:-1, :-1], lon[:-1, 1:], lon[1:, 1:], lon[1:, :-1]), axis=-1) |
| 301 | + xf_low = xf.min(axis=-1) |
| 302 | + xf_high = xf.max(axis=-1) |
| 303 | + xbounds = np.stack([xf_low, xf_high], axis=-1) |
| 304 | + |
| 305 | + yf = np.stack((lat[:-1, :-1], lat[:-1, 1:], lat[1:, 1:], lat[1:, :-1]), axis=-1) |
| 306 | + yf_low = yf.min(axis=-1) |
| 307 | + yf_high = yf.max(axis=-1) |
| 308 | + ybounds = np.stack([yf_low, yf_high], axis=-1) |
| 309 | + |
| 310 | + return ybounds, xbounds |
0 commit comments