Skip to content

Commit d88a2b3

Browse files
authored
Merge pull request #1408 from UXARRAY/rajeeja/fix_issue1406
Fix .sel() wrapper when a slice indexer and method are used together
2 parents d42921f + 261ee90 commit d88a2b3

File tree

6 files changed

+93
-4
lines changed

6 files changed

+93
-4
lines changed

test/core/test_api.py

Lines changed: 7 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -36,17 +36,21 @@ def test_open_dataset(gridpath, datasetpath, mesh_constants):
3636
nt.assert_equal(len(uxds_var2_ne30.uxgrid._ds.data_vars), mesh_constants['DATAVARS_outCSne30'])
3737
nt.assert_equal(uxds_var2_ne30.source_datasets, str(data_path))
3838

39-
def test_open_mf_dataset(gridpath, test_data_dir, mesh_constants):
39+
def test_open_mf_dataset(gridpath, datasetpath, mesh_constants):
4040
"""Loads multiple datasets with their grid topology file using
4141
uxarray's open_dataset call."""
4242

4343
grid_path = gridpath("ugrid", "outCSne30", "outCSne30.ug")
44-
dsfiles_mf_ne30 = str(test_data_dir) + "/ugrid/outCSne30/outCSne30_*.nc"
44+
dsfiles_mf_ne30 = datasetpath(
45+
"ugrid",
46+
"outCSne30",
47+
["outCSne30_var2.nc", "outCSne30_vortex.nc"],
48+
)
4549
uxds_mf_ne30 = ux.open_mfdataset(grid_path, dsfiles_mf_ne30)
4650

4751
nt.assert_equal(uxds_mf_ne30.uxgrid.node_lon.size, mesh_constants['NNODES_outCSne30'])
4852
nt.assert_equal(len(uxds_mf_ne30.uxgrid._ds.data_vars), mesh_constants['DATAVARS_outCSne30'])
49-
nt.assert_equal(uxds_mf_ne30.source_datasets, dsfiles_mf_ne30)
53+
nt.assert_equal(uxds_mf_ne30.source_datasets, str(dsfiles_mf_ne30))
5054

5155
def test_open_grid(gridpath, mesh_constants):
5256
"""Loads only a grid topology file using uxarray's open_grid call."""

test/core/test_dataset.py

Lines changed: 32 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -67,3 +67,35 @@ def test_get_dual(gridpath, datasetpath):
6767

6868
assert isinstance(dual, UxDataset)
6969
assert len(uxds.data_vars) == len(dual.data_vars)
70+
71+
72+
def _load_sel_subset(gridpath, datasetpath):
73+
grid_file = gridpath("ugrid", "outCSne30", "outCSne30.ug")
74+
data_file = datasetpath("ugrid", "outCSne30", "outCSne30_sel_timeseries.nc")
75+
uxds = ux.open_dataset(grid_file, data_file)
76+
base_time = np.datetime64("2018-04-28T00:00:00")
77+
offsets = np.arange(uxds.sizes["time"], dtype="timedelta64[h]")
78+
uxds = uxds.assign_coords(time=(base_time + offsets).astype("datetime64[ns]"))
79+
return uxds
80+
81+
82+
def test_sel_time_slice(gridpath, datasetpath):
83+
uxds = _load_sel_subset(gridpath, datasetpath)
84+
85+
times = uxds["time"].values
86+
sliced = uxds.sel(time=slice(times[0], times[2]))
87+
88+
assert sliced.dims["time"] == 3
89+
np.testing.assert_array_equal(sliced["time"].values, times[:3])
90+
91+
92+
def test_sel_method_forwarded(gridpath, datasetpath):
93+
uxds = _load_sel_subset(gridpath, datasetpath)
94+
95+
target = np.datetime64("2018-04-28T02:20:00")
96+
nearest = uxds.sel(time=target, method="nearest")
97+
98+
np.testing.assert_array_equal(
99+
nearest["time"].values,
100+
np.array(uxds["time"].values[2], dtype="datetime64[ns]"),
101+
)
Lines changed: 11 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,11 @@
1+
# outCSne30 UGRID assets
2+
3+
This folder houses the NE30 cubed-sphere mesh (5400 faces) and a few companion datasets used throughout the test suite.
4+
5+
## Files
6+
7+
- `outCSne30.ug` – canonical NE30 UGRID mesh with 5400 faces (nMesh2_face), max 4 vertices per face, and 5402 nodes. Provides `Mesh2_face_nodes`, `Mesh2_node_x`, and `Mesh2_node_y`.
8+
- `outCSne30_var2.nc` – single-variable dataset (`var2(ncol=5400)`) used for basic integration/zonal workflows.
9+
- `outCSne30_vortex.nc` – single-variable dataset (`psi(ncol=5400)`) containing the barotropic vortex test case.
10+
- `outCSne30_sel_timeseries.nc` – synthetic 6-step hourly time series (`psi(time=6, ncol=5400)`) where each time slice ramps with longitude; input for the `.sel()` regression tests. Regenerate via `python generate_sel_timeseries.py`.
11+
- `generate_sel_timeseries.py` – helper script that rewrites the synthetic time series using `outCSne30_var2.nc` as a template for face count.
Lines changed: 36 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,36 @@
1+
"""Generate the synthetic NE30 time series used in sel() regression tests.
2+
3+
Run from repo root: ``python test/meshfiles/ugrid/outCSne30/generate_sel_timeseries.py``
4+
"""
5+
6+
from pathlib import Path
7+
8+
import numpy as np
9+
import xarray as xr
10+
11+
TIME_STEPS = 6
12+
13+
def main() -> None:
14+
base = Path(__file__).parent
15+
16+
data_src = base / "outCSne30_var2.nc"
17+
template = xr.open_dataset(data_src)
18+
n_face = template.dims["ncol"]
19+
template.close()
20+
21+
times = np.arange(TIME_STEPS, dtype=np.int64)
22+
time_vals = (np.datetime64("2018-04-28T00:00:00") + times * np.timedelta64(1, "h"))
23+
faces = np.arange(n_face, dtype=np.float32)
24+
field = (times[:, None].astype(np.float32) + faces[None, :] / 100.0)
25+
26+
ds = xr.Dataset(
27+
{"psi": (("time", "ncol"), field)},
28+
coords={"time": time_vals},
29+
attrs={"source": "Synthetic field for sel regression"},
30+
)
31+
data_path = base / "outCSne30_sel_timeseries.nc"
32+
ds.to_netcdf(data_path)
33+
34+
35+
if __name__ == "__main__":
36+
main()
135 KB
Binary file not shown.

uxarray/core/dataset.py

Lines changed: 7 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -690,7 +690,13 @@ def sel(
690690
self, indexers=None, method=None, tolerance=None, drop=False, **indexers_kwargs
691691
):
692692
return UxDataset(
693-
self.to_xarray().sel(indexers, tolerance, drop, **indexers_kwargs),
693+
self.to_xarray().sel(
694+
indexers=indexers,
695+
method=method,
696+
tolerance=tolerance,
697+
drop=drop,
698+
**indexers_kwargs,
699+
),
694700
uxgrid=self.uxgrid,
695701
)
696702

0 commit comments

Comments
 (0)