Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Convert dfsu2d to UGRID netcdf #693

Draft
wants to merge 5 commits into
base: main
Choose a base branch
from
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion .github/workflows/full_test.yml
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@ jobs:
strategy:
matrix:
os: [ubuntu-latest, windows-latest]
python-version: [3.8, "3.12"]
python-version: [3.9, "3.12"]

steps:
- uses: actions/checkout@v3
Expand Down
2 changes: 1 addition & 1 deletion .github/workflows/legacy_test.yml
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@ jobs:
- name: Set up Python
uses: actions/setup-python@v4
with:
python-version: 3.8
python-version: 3.9

- name: Install older dependencies
run: |
Expand Down
48 changes: 48 additions & 0 deletions mikeio/dataset/_dataarray.py
Original file line number Diff line number Diff line change
Expand Up @@ -1873,6 +1873,54 @@ def to_xarray(self) -> "xarray.DataArray":
)
return xr_da

def to_xugrid(self):
import xugrid as xu
import xarray as xr

assert isinstance(
self.geometry, GeometryFM2D
), "Only flexible mesh 2D data is supported"
# face_code_connectivity is slightly different than the element_table

g = self.geometry

if g.is_tri_only:
conn = np.zeros(shape=(self.geometry.n_elements, 3), dtype=np.int32)
for i in range(g.n_elements):
conn[i, 0:3] = g.element_table[i]
else:
conn = np.zeros(shape=(self.geometry.n_elements, 4), dtype=np.int32)
conn[:] = -1

for i in range(g.n_elements):
if len(g.element_table[i]) == 3:
conn[i, 0:3] = g.element_table[i]
else:
conn[i] = g.element_table[i]

xn = g.node_coordinates[:, 0]
yn = g.node_coordinates[:, 1]

grid = xu.Ugrid2d(xn, yn, fill_value=-1, face_node_connectivity=conn)

if "time" not in self.dims:
da = xr.DataArray(
name=self.name,
data=self.to_numpy(),
coords={},
dims=[grid.face_dimension],
)
else:
da = xr.DataArray(
name=self.name,
data=self.to_numpy(),
coords={"time": self.time},
dims=["time", grid.face_dimension],
)

uda = xu.UgridDataArray(da, grid)
return uda

# ===============================================

def __repr__(self) -> str:
Expand Down
9 changes: 9 additions & 0 deletions mikeio/dataset/_dataset.py
Original file line number Diff line number Diff line change
Expand Up @@ -1898,6 +1898,15 @@ def to_xarray(self) -> "xarray.Dataset":
data = {da.name: da.to_xarray() for da in self}
return xarray.Dataset(data)

def to_xugrid(self):
import xugrid as xu

data = [da.to_xugrid() for da in self]
ds = xu.UgridDataset(grids=data[0].grid)
for da in data:
ds[da.name] = da
return ds

# ===============================================

def __repr__(self) -> str:
Expand Down
3 changes: 3 additions & 0 deletions mypy.ini
Original file line number Diff line number Diff line change
Expand Up @@ -33,4 +33,7 @@ ignore_missing_imports = True
[mypy-yaml.*]
ignore_missing_imports = True

[mypy-xugrid.*]
ignore_missing_imports = True


2 changes: 1 addition & 1 deletion pyproject.toml
Original file line number Diff line number Diff line change
Expand Up @@ -58,7 +58,7 @@ dev = ["pytest",
"mypy==1.6.1",
]

test = ["pytest", "pytest-cov", "xarray","mypy==1.6.1","shapely","pyproj"]
test = ["pytest", "pytest-cov", "xarray","mypy==1.6.1","shapely","pyproj","xugrid"]

notebooks= [
"nbformat",
Expand Down
26 changes: 26 additions & 0 deletions tests/test_dfsu.py
Original file line number Diff line number Diff line change
Expand Up @@ -970,3 +970,29 @@ def test_writing_non_equdistant_dfsu_is_not_possible(tmp_path):
with pytest.raises(ValueError, match="equidistant"):
fp = tmp_path / "not_gonna_work.dfsu"
dss.to_dfs(fp)


def test_convert_dfsu2d_to_xugrid(tmp_path):

import xugrid

ds = mikeio.read("tests/testdata/NorthSea_HD_and_windspeed.dfsu")
ds.sel(y=55, x=0).isel(time=-1).to_numpy() == pytest.approx(0.25443718)
xu_ds = ds.to_xugrid()

xu_ds.ugrid.sel(y=55, x=0)["Surface elevation"].isel(
time=-1
).to_numpy() == pytest.approx(0.25443718)
assert xu_ds.time[-1] == ds.time[-1]

fp = tmp_path / "north_sea.nc"
xu_ds.ugrid.to_netcdf(fp)

xu_ds2 = xugrid.open_dataset(fp)
assert xu_ds2.time[-1] == ds.time[-1]


def test_mixed_mesh_time_invariant_to_xugrid():
ds = mikeio.read("tests/testdata/FakeLake.dfsu", time=0)
xu_ds = ds.to_xugrid()
assert len(xu_ds.mesh2d_nFaces) == 1011
Loading