-
Notifications
You must be signed in to change notification settings - Fork 228
Implement gmt xarray BackendEntrypoint #3919
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
Changes from 35 commits
b81aa06
2b98a65
f41a812
b71dcf5
60c6726
bf09a39
e622878
b9e3eb4
3a109ac
d764f78
98910f8
126560e
436f0fc
328b6ab
06aa39d
c79e50a
c0ccd4b
c5ebf9f
17fd9d7
1c8af3a
1d9f21a
ca9132f
9669bcf
f722b47
a952df9
0c6b8a3
3a256a6
64c50c9
e6ee6fa
5b90a8a
31d1999
6a4e1a4
91df0ee
c97bd82
2e5fb6b
b136e26
19c5252
3bb1f44
8984c1c
f2446fc
52c0b7e
143a375
acab8b2
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,74 @@ | ||
""" | ||
Tests for xarray 'gmt' backend engine. | ||
""" | ||
|
||
import re | ||
|
||
import numpy as np | ||
import numpy.testing as npt | ||
import pytest | ||
import xarray as xr | ||
from pygmt.enums import GridRegistration, GridType | ||
from pygmt.exceptions import GMTInvalidInput | ||
|
||
|
||
def test_xarray_backend_gmt_open_nc_grid(): | ||
""" | ||
Ensure that passing engine='gmt' to xarray.open_dataarray works for opening NetCDF | ||
grids. | ||
""" | ||
with xr.open_dataarray( | ||
"@static_earth_relief.nc", engine="gmt", raster_kind="grid" | ||
) as da: | ||
assert da.sizes == {"lat": 14, "lon": 8} | ||
assert da.dtype == "float32" | ||
assert da.gmt.registration == GridRegistration.PIXEL | ||
assert da.gmt.gtype == GridType.GEOGRAPHIC | ||
|
||
|
||
def test_xarray_backend_gmt_open_tif_image(): | ||
""" | ||
Ensure that passing engine='gmt' to xarray.open_dataarray works for opening GeoTIFF | ||
images. | ||
""" | ||
with xr.open_dataarray("@earth_day_01d", engine="gmt", raster_kind="image") as da: | ||
assert da.sizes == {"band": 3, "y": 180, "x": 360} | ||
assert da.dtype == "uint8" | ||
assert da.gmt.registration == GridRegistration.PIXEL | ||
assert da.gmt.gtype == GridType.GEOGRAPHIC | ||
|
||
|
||
def test_xarray_backend_gmt_load_grd_grid(): | ||
""" | ||
Ensure that passing engine='gmt' to xarray.load_dataarray works for loading GRD | ||
grids. | ||
""" | ||
da = xr.load_dataarray( | ||
"@earth_relief_20m_holes.grd", engine="gmt", raster_kind="grid" | ||
) | ||
# Ensure data is in memory. | ||
assert isinstance(da.data, np.ndarray) | ||
npt.assert_allclose(da.min(), -4929.5) | ||
assert da.sizes == {"lat": 31, "lon": 31} | ||
assert da.dtype == "float32" | ||
assert da.gmt.registration == GridRegistration.GRIDLINE | ||
assert da.gmt.gtype == GridType.GEOGRAPHIC | ||
|
||
|
||
def test_xarray_backend_gmt_read_invalid_kind(): | ||
""" | ||
Check that xarray.open_dataarray(..., engine="gmt") fails with missing or incorrect | ||
'raster_kind'. | ||
""" | ||
with pytest.raises( | ||
TypeError, | ||
match=re.escape( | ||
"GMTBackendEntrypoint.open_dataset() missing 1 required keyword-only argument: 'raster_kind'" | ||
), | ||
): | ||
xr.open_dataarray("nokind.nc", engine="gmt") | ||
|
||
with pytest.raises(GMTInvalidInput): | ||
xr.open_dataarray( | ||
filename_or_obj="invalid.tif", engine="gmt", raster_kind="invalid" | ||
) |
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Considering moving this file to There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Sounds good |
Original file line number | Diff line number | Diff line change | ||||||||||||||||
---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|---|
@@ -0,0 +1,106 @@ | ||||||||||||||||||
""" | ||||||||||||||||||
An xarray backend for reading raster grid/image files using the 'gmt' engine. | ||||||||||||||||||
""" | ||||||||||||||||||
|
||||||||||||||||||
from typing import Literal | ||||||||||||||||||
|
||||||||||||||||||
import xarray as xr | ||||||||||||||||||
from pygmt._typing import PathLike | ||||||||||||||||||
from pygmt.clib import Session | ||||||||||||||||||
from pygmt.exceptions import GMTInvalidInput | ||||||||||||||||||
from pygmt.helpers import build_arg_list | ||||||||||||||||||
from pygmt.src.which import which | ||||||||||||||||||
from xarray.backends import BackendEntrypoint | ||||||||||||||||||
|
||||||||||||||||||
|
||||||||||||||||||
class GMTBackendEntrypoint(BackendEntrypoint): | ||||||||||||||||||
""" | ||||||||||||||||||
Xarray backend to read raster grid/image files using 'gmt' engine. | ||||||||||||||||||
|
||||||||||||||||||
Internally, GMT uses the NetCDF C library to read NetCDF files, and GDAL for GeoTIFF | ||||||||||||||||||
and other raster formats. | ||||||||||||||||||
|
||||||||||||||||||
When using :py:func:`xarray.open_dataarray` or :py:func:`xarray.load_dataarray` with | ||||||||||||||||||
``engine="gmt"``, pass the ``raster_kind`` parameter that can be either: | ||||||||||||||||||
|
||||||||||||||||||
- ``"grid"`` - for reading single-band raster grids | ||||||||||||||||||
- ``"image"`` - for reading multi-band raster images | ||||||||||||||||||
|
||||||||||||||||||
Examples | ||||||||||||||||||
-------- | ||||||||||||||||||
Read a single-band NetCDF file using ``raster_kind="grid"`` | ||||||||||||||||||
|
||||||||||||||||||
>>> import pygmt | ||||||||||||||||||
>>> import xarray as xr | ||||||||||||||||||
>>> | ||||||||||||||||||
>>> da_grid = xr.open_dataarray( | ||||||||||||||||||
... "@static_earth_relief.nc", engine="gmt", raster_kind="grid" | ||||||||||||||||||
... ) | ||||||||||||||||||
>>> da_grid # doctest: +NORMALIZE_WHITESPACE, +ELLIPSIS | ||||||||||||||||||
<xarray.DataArray 'z' (lat: 14, lon: 8)>... | ||||||||||||||||||
[112 values with dtype=float32] | ||||||||||||||||||
Coordinates: | ||||||||||||||||||
* lat (lat) float64... -23.5 -22.5 -21.5 -20.5 ... -12.5 -11.5 -10.5 | ||||||||||||||||||
* lon (lon) float64... -54.5 -53.5 -52.5 -51.5 -50.5 -49.5 -48.5 -47.5 | ||||||||||||||||||
Attributes: | ||||||||||||||||||
Conventions: CF-1.7 | ||||||||||||||||||
title: Produced by grdcut | ||||||||||||||||||
history: grdcut @earth_relief_01d_p -R-55/-47/-24/-10 -Gstatic_eart... | ||||||||||||||||||
description: Reduced by Gaussian Cartesian filtering (111.2 km fullwidt... | ||||||||||||||||||
actual_range: [190. 981.] | ||||||||||||||||||
long_name: elevation (m) | ||||||||||||||||||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. A workaround for issue #3919 (comment)
Suggested change
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I put the ellipisis after the colon in 'Attributes:' line, and it seemed to remove the Only downside is the ellipsis shows as |
||||||||||||||||||
|
||||||||||||||||||
Read a multi-band GeoTIFF file using ``raster_kind="image"`` | ||||||||||||||||||
|
||||||||||||||||||
>>> da_image = xr.open_dataarray( | ||||||||||||||||||
... "@earth_night_01d", engine="gmt", raster_kind="image" | ||||||||||||||||||
... ) | ||||||||||||||||||
>>> da_image # doctest: +NORMALIZE_WHITESPACE, +ELLIPSIS | ||||||||||||||||||
<xarray.DataArray 'z' (band: 3, y: 180, x: 360)>... | ||||||||||||||||||
[194400 values with dtype=uint8] | ||||||||||||||||||
Coordinates: | ||||||||||||||||||
* y (y) float64... 89.5 88.5 87.5 86.5 ... -86.5 -87.5 -88.5 -89.5 | ||||||||||||||||||
* x (x) float64... -179.5 -178.5 -177.5 -176.5 ... 177.5 178.5 179.5 | ||||||||||||||||||
* band (band) uint8... 1 2 3 | ||||||||||||||||||
Attributes: | ||||||||||||||||||
long_name: z | ||||||||||||||||||
weiji14 marked this conversation as resolved.
Show resolved
Hide resolved
|
||||||||||||||||||
""" | ||||||||||||||||||
|
||||||||||||||||||
description = "Open raster (.grd, .nc or .tif) files in Xarray via GMT." | ||||||||||||||||||
open_dataset_parameters = ("filename_or_obj", "raster_kind") | ||||||||||||||||||
url = "https://pygmt.org/dev/api/generated/pygmt.GMTBackendEntrypoint.html" | ||||||||||||||||||
|
||||||||||||||||||
def open_dataset( # type: ignore[override] | ||||||||||||||||||
self, | ||||||||||||||||||
filename_or_obj: PathLike, | ||||||||||||||||||
*, | ||||||||||||||||||
drop_variables=None, # noqa: ARG002 | ||||||||||||||||||
raster_kind: Literal["grid", "image"], | ||||||||||||||||||
# other backend specific keyword arguments | ||||||||||||||||||
# `chunks` and `cache` DO NOT go here, they are handled by xarray | ||||||||||||||||||
) -> xr.Dataset: | ||||||||||||||||||
""" | ||||||||||||||||||
Backend open_dataset method used by Xarray in :py:func:`~xarray.open_dataset`. | ||||||||||||||||||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. With the fix in #3927, sphinx-autogen will generate a stub file for this method. The method documentation is not ideal, at least one "Parameters" section is needed. |
||||||||||||||||||
""" | ||||||||||||||||||
if raster_kind not in {"grid", "image"}: | ||||||||||||||||||
msg = f"Invalid raster kind: '{raster_kind}'. Valid values are 'grid' or 'image'." | ||||||||||||||||||
raise GMTInvalidInput(msg) | ||||||||||||||||||
|
||||||||||||||||||
with Session() as lib: | ||||||||||||||||||
with lib.virtualfile_out(kind=raster_kind) as voutfile: | ||||||||||||||||||
kwdict = {"T": {"grid": "g", "image": "i"}[raster_kind]} | ||||||||||||||||||
lib.call_module( | ||||||||||||||||||
module="read", | ||||||||||||||||||
args=[filename_or_obj, voutfile, *build_arg_list(kwdict)], | ||||||||||||||||||
) | ||||||||||||||||||
|
||||||||||||||||||
raster: xr.DataArray = lib.virtualfile_to_raster( | ||||||||||||||||||
vfname=voutfile, kind=raster_kind | ||||||||||||||||||
) | ||||||||||||||||||
# Add "source" encoding | ||||||||||||||||||
source = which(fname=filename_or_obj) | ||||||||||||||||||
raster.encoding["source"] = ( | ||||||||||||||||||
source[0] if isinstance(source, list) else source | ||||||||||||||||||
) | ||||||||||||||||||
_ = raster.gmt # Load GMTDataArray accessor information | ||||||||||||||||||
return raster.to_dataset() |
Uh oh!
There was an error while loading. Please reload this page.