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

Wrap GMT's standard data type GMT_GRID for grids #2398

Merged
merged 15 commits into from
Apr 16, 2024
Merged
Show file tree
Hide file tree
Changes from 10 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
1 change: 1 addition & 0 deletions doc/api/index.rst
Original file line number Diff line number Diff line change
Expand Up @@ -292,6 +292,7 @@ Python objects to and from GMT virtual files:
clib.Session.virtualfile_in
clib.Session.virtualfile_out
clib.Session.virtualfile_to_dataset
clib.Session.virtualfile_to_raster

Low level access (these are mostly used by the :mod:`pygmt.clib` package):

Expand Down
71 changes: 70 additions & 1 deletion pygmt/clib/session.py
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@

import numpy as np
import pandas as pd
import xarray as xr
from packaging.version import Version
from pygmt.clib.conversion import (
array_to_datetime,
Expand Down Expand Up @@ -1739,7 +1740,9 @@
return c_inquire_virtualfile(self.session_pointer, vfname.encode())

def read_virtualfile(
self, vfname: str, kind: Literal["dataset", "grid", None] = None
self,
vfname: str,
kind: Literal["dataset", "grid", "image", "cube", None] = None,
):
"""
Read data from a virtual file and optionally cast into a GMT data container.
Expand Down Expand Up @@ -1798,6 +1801,8 @@
# _GMT_DATASET).
if kind is None: # Return the ctypes void pointer
return pointer
if kind in ["image", "cube"]:
raise NotImplementedError(f"kind={kind} is not supported yet.")

Check warning on line 1805 in pygmt/clib/session.py

View check run for this annotation

Codecov / codecov/patch

pygmt/clib/session.py#L1805

Added line #L1805 was not covered by tests
dtype = {"dataset": _GMT_DATASET, "grid": _GMT_GRID}[kind]
return ctp.cast(pointer, ctp.POINTER(dtype))

Expand Down Expand Up @@ -1946,6 +1951,70 @@
return result.to_numpy()
return result # pandas.DataFrame output

def virtualfile_to_raster(
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

In this method, kind can be grid, cube, image or None. If None, we will call Session.inquire_virtualfile to determine the family of the virtualfile.

Most GMT modules can write grids only, so we can set kind="grid" to avoid calling Session.inquire_virtualfile.

Here are some exceptions:

  1. grdmix writes images
  2. grdcut writes grids or images
  3. grdimage -A writes images [Unimplemented in PyGMT]
  4. greenspine writes grids or cubes
  5. grdinterpolate writes grids or cubes

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Most GMT modules can write grids only, so we can set kind="grid" to avoid calling Session.inquire_virtualfile.

Yes, I think we can set kind="grid", unless we know if GMT will allow passing GMT_IMAGE to certain grd* modules in the future. For example, maybe grdclip would allow GMT_IMAGE inputs at some point, and so we might want to set kind=None for future-proofing.

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

We need to do more refactoring if grdclip would support GMT_IMAGE in the future, since we need to check the kind of the input data (grid in, grid out or image in, image out), similar to what we have to do to grdcut (xref: #2398 (comment)).

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ok, I meant if upstream GMT were to support passing GMT_IMAGE into grdclip in the future, but not sure when that would happen. Ok to just default to kind='grid' for now, since that's the current behaviour I think.

self,
vfname: str,
kind: Literal["grid", "image", "cube", None] = None,
seisman marked this conversation as resolved.
Show resolved Hide resolved
outgrid: str | None = None,
) -> xr.DataArray | None:
"""
Output a raster data stored in a virtual file to an :class:`xarray.DataArray`
seisman marked this conversation as resolved.
Show resolved Hide resolved
object.

The raster data can be a grid, an image or a cube.

Parameters
----------
vfname
The virtual file name that stores the result grid.
seisman marked this conversation as resolved.
Show resolved Hide resolved
kind
Type of the raster data. Valid values are ``"grid"``, ``"image"``,
``"cube"`` or ``None``. If ``None``, will inquire the data type from the
virtual file name.
outgrid
Name of the output grid/image/cube. If specified, it means the raster data
was already saved into an actual file and will return ``None``.

Returns
-------
result
The result grid/image/cube. If ``outgrid`` is specified, return ``None``.

Examples
--------
>>> from pathlib import Path
>>> from pygmt.clib import Session
>>> from pygmt.helpers import GMTTempFile
>>> with Session() as lib:
... # file output
... with GMTTempFile(suffix=".nc") as tmpfile:
... outgrid = tmpfile.name
... with lib.virtualfile_out(kind="grid", fname=outgrid) as voutgrd:
... lib.call_module("read", f"@earth_relief_01d_g {voutgrd} -Tg")
... result = lib.virtualfile_to_raster(
... vfname=voutgrd, outgrid=outgrid
... )
... assert result == None
... assert Path(outgrid).stat().st_size > 0
...
... # xarray.DataArray output
... outgrid = None
... with lib.virtualfile_out(kind="grid", fname=outgrid) as voutgrd:
... lib.call_module("read", f"@earth_relief_01d_g {voutgrd} -Tg")
... result = lib.virtualfile_to_raster(vfname=voutgrd, outgrid=outgrid)
... assert isinstance(result, xr.DataArray)
"""
if outgrid is not None:
return None
if kind is None: # Inquire the data family from the virtualfile
family = self.inquire_virtualfile(vfname)
kind = { # type: ignore[assignment]
self["GMT_IS_GRID"]: "grid",
self["GMT_IS_IMAGE"]: "image",
self["GMT_IS_CUBE"]: "cube",
}[family]
return self.read_virtualfile(vfname, kind=kind).contents.to_dataarray()

def extract_region(self):
"""
Extract the WESN bounding box of the currently active figure.
Expand Down
192 changes: 191 additions & 1 deletion pygmt/datatypes/grid.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,197 @@
"""

import ctypes as ctp
from typing import ClassVar

import numpy as np
import xarray as xr
from pygmt.datatypes.header import _GMT_GRID_HEADER, gmt_grdfloat


class _GMT_GRID(ctp.Structure): # noqa: N801
pass
"""
GMT grid structure for holding a grid and its header.

The class is only meant for internal use and is not exposed to users. See the GMT
seisman marked this conversation as resolved.
Show resolved Hide resolved
source code gmt_resources.h for the original C structure definitions.

Examples
--------
>>> from pygmt.clib import Session
>>> with Session() as lib:
... with lib.virtualfile_out(kind="grid") as voutgrd:
... lib.call_module("read", f"@static_earth_relief.nc {voutgrd} -Tg")
... # Read the grid from the virtual file
... grid = lib.read_virtualfile(voutgrd, kind="grid").contents
... # The grid header
... header = grid.header.contents
... # Access the header properties
... print(header.n_rows, header.n_columns, header.registration)
... print(header.wesn[:], header.z_min, header.z_max, header.inc[:])
... print(header.z_scale_factor, header.z_add_offset)
... print(header.x_units, header.y_units, header.z_units)
... print(header.title)
... print(header.command)
... print(header.remark)
... print(header.nm, header.size, header.complex_mode)
... print(header.type, header.n_bands, header.mx, header.my)
... print(header.pad[:])
... print(header.mem_layout, header.nan_value, header.xy_off)
... # The x and y coordinates
... print(grid.x[: header.n_columns])
... print(grid.y[: header.n_rows])
... # The data array (with paddings)
... data = np.reshape(
... grid.data[: header.mx * header.my], (header.my, header.mx)
... )
... # The data array (without paddings)
... pad = header.pad[:]
... data = data[pad[2] : header.my - pad[3], pad[0] : header.mx - pad[1]]
... print(data)
14 8 1
[-55.0, -47.0, -24.0, -10.0] 190.0 981.0 [1.0, 1.0]
1.0 0.0
b'longitude [degrees_east]' b'latitude [degrees_north]' b'elevation (m)'
b'Produced by grdcut'
b'grdcut @earth_relief_01d_p -R-55/-47/-24/-10 -Gstatic_earth_relief.nc'
b'Reduced by Gaussian Cartesian filtering (111.2 km fullwidth) from ...'
112 216 0
18 1 12 18
[2, 2, 2, 2]
b'' nan 0.5
[-54.5, -53.5, -52.5, -51.5, -50.5, -49.5, -48.5, -47.5]
[-10.5, -11.5, -12.5, -13.5, -14.5, -15.5, ..., -22.5, -23.5]
[[347.5 331.5 309. 282. 190. 208. 299.5 348. ]
[349. 313. 325.5 247. 191. 225. 260. 452.5]
[345.5 320. 335. 292. 207.5 247. 325. 346.5]
[450.5 395.5 366. 248. 250. 354.5 550. 797.5]
[494.5 488.5 357. 254.5 286. 484.5 653.5 930. ]
[601. 526.5 535. 299. 398.5 645. 797.5 964. ]
[308. 595.5 555.5 556. 580. 770. 927. 920. ]
[521.5 682.5 796. 886. 571.5 638.5 739.5 881.5]
[310. 521.5 757. 570.5 538.5 524. 686.5 794. ]
[561.5 539. 446.5 481.5 439.5 553. 726.5 981. ]
[557. 435. 385.5 345.5 413.5 496. 519.5 833.5]
[373. 367.5 349. 352.5 419.5 428. 570. 667.5]
[383. 284.5 344.5 394. 491. 556.5 578.5 618.5]
[347.5 344.5 386. 640.5 617. 579. 646.5 671. ]]
"""

_fields_: ClassVar = [
# Pointer to full GMT header for grid
("header", ctp.POINTER(_GMT_GRID_HEADER)),
# Pointer to grid data
("data", ctp.POINTER(gmt_grdfloat)),
# Pointer to x coordinate vector
("x", ctp.POINTER(ctp.c_double)),
# Pointer to y coordinate vector
("y", ctp.POINTER(ctp.c_double)),
# Low-level information for GMT use only
("hidden", ctp.c_void_p),
]

def to_dataarray(self) -> xr.DataArray:
"""
Convert a _GMT_GRID object to a :class:`xarray.DataArray` object.

Returns
-------
dataarray
A :class:`xr.DataArray` object.

Examples
--------
>>> from pygmt.clib import Session
>>> with Session() as lib:
... with lib.virtualfile_out(kind="grid") as voutgrd:
... lib.call_module("read", f"@static_earth_relief.nc {voutgrd} -Tg")
... # Read the grid from the virtual file
... grid = lib.read_virtualfile(voutgrd, kind="grid")
... # Convert to xarray.DataArray and use it later
... da = grid.contents.to_dataarray()
>>> da # doctest: +NORMALIZE_WHITESPACE, +ELLIPSIS
<xarray.DataArray 'z' (lat: 14, lon: 8)>...
array([[347.5, 344.5, 386. , 640.5, 617. , 579. , 646.5, 671. ],
[383. , 284.5, 344.5, 394. , 491. , 556.5, 578.5, 618.5],
[373. , 367.5, 349. , 352.5, 419.5, 428. , 570. , 667.5],
[557. , 435. , 385.5, 345.5, 413.5, 496. , 519.5, 833.5],
[561.5, 539. , 446.5, 481.5, 439.5, 553. , 726.5, 981. ],
[310. , 521.5, 757. , 570.5, 538.5, 524. , 686.5, 794. ],
[521.5, 682.5, 796. , 886. , 571.5, 638.5, 739.5, 881.5],
[308. , 595.5, 555.5, 556. , 580. , 770. , 927. , 920. ],
[601. , 526.5, 535. , 299. , 398.5, 645. , 797.5, 964. ],
[494.5, 488.5, 357. , 254.5, 286. , 484.5, 653.5, 930. ],
[450.5, 395.5, 366. , 248. , 250. , 354.5, 550. , 797.5],
[345.5, 320. , 335. , 292. , 207.5, 247. , 325. , 346.5],
[349. , 313. , 325.5, 247. , 191. , 225. , 260. , 452.5],
[347.5, 331.5, 309. , 282. , 190. , 208. , 299.5, 348. ]])
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_ea...
description: Reduced by Gaussian Cartesian filtering (111.2 km fullwi...
long_name: elevation (m)
actual_range: [190. 981.]
>>> da.coords["lon"] # doctest: +NORMALIZE_WHITESPACE, +ELLIPSIS
<xarray.DataArray 'lon' (lon: 8)>...
array([-54.5, -53.5, -52.5, -51.5, -50.5, -49.5, -48.5, -47.5])
Coordinates:
* lon (lon) float64... -54.5 -53.5 -52.5 -51.5 -50.5 -49.5 -48.5 -47.5
Attributes:
long_name: longitude
units: degrees_east
standard_name: longitude
axis: X
actual_range: [-55. -47.]
>>> da.coords["lat"] # doctest: +NORMALIZE_WHITESPACE, +ELLIPSIS
<xarray.DataArray 'lat' (lat: 14)>...
array([-23.5, -22.5, -21.5, -20.5, -19.5, -18.5, -17.5, -16.5, -15.5, -14.5,
-13.5, -12.5, -11.5, -10.5])
Coordinates:
* lat (lat) float64... -23.5 -22.5 -21.5 -20.5 ... -12.5 -11.5 -10.5
Attributes:
long_name: latitude
units: degrees_north
standard_name: latitude
axis: Y
actual_range: [-24. -10.]
>>> da.gmt.registration, da.gmt.gtype
(1, 1)
"""
# The grid header
header = self.header.contents

# Get dimensions and their attributes from the header.
dims, dim_attrs = header.dims, header.dim_attrs
# The coordinates, given as a tuple of the form (dims, data, attrs)
coords = [
(dims[0], self.y[: header.n_rows], dim_attrs[0]),
(dims[1], self.x[: header.n_columns], dim_attrs[1]),
]

# The data array without paddings
pad = header.pad[:]
data = np.reshape(self.data[: header.mx * header.my], (header.my, header.mx))[
pad[2] : header.my - pad[3], pad[0] : header.mx - pad[1]
]

# Create the xarray.DataArray object
grid = xr.DataArray(
data, coords=coords, name=header.name, attrs=header.data_attrs
)

# Flip the coordinates and data if necessary so that coordinates are ascending.
# `grid.sortby(list(grid.dims))` sometimes causes crashes.
# The solution comes from https://github.com/pydata/xarray/discussions/6695.
for dim in grid.dims:
if grid[dim][0] > grid[dim][1]:
grid = grid.isel({dim: slice(None, None, -1)})
Comment on lines +188 to +193
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Could you elaborate a bit on this part? Sort ascending is probably needed for the x (longitude) dimension. But for the y (latitude) dimension, it is quite common to sort in descending order (y going from North Pole +90 to South Pole -90). What is the typical sorting order returned from GMT?

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think GMT requires ascending order for longitude/latitude. Similar codes are already used at (The codes below use sortby, but we can't use it here due to a mysterious reason (xref #2398 (comment)).)

if any(i < 0 for i in inc): # Sort grid when there are negative increments
inc = [abs(i) for i in inc]
grid = grid.sortby(variables=list(grid.dims), ascending=True)

I tried to find some evidence from the GMT source codes but failed because the codes are so complicated. Perhaps the most straightforward evidence is from following lines (https://github.com/GenericMappingTools/gmt/blob/e493ccc220d9e997abb655cb461142ed2e7cc777/src/gmt_nc.c#L754):

		/* Determine grid step */
		header->inc[GMT_X] = gmt_M_get_inc (GMT, dummy[0], dummy[1], header->n_columns, registration);
		if (gmt_M_is_dnan(header->inc[GMT_X]) || gmt_M_is_zero (header->inc[GMT_X])) header->inc[GMT_X] = 1.0;
		if (header->n_columns == 1) registration = GMT_GRID_PIXEL_REG;	/* The only way to have a grid like that */

		/* Extend x boundaries by half if we found pixel registration */
		if (registration == GMT_GRID_NODE_REG && header->registration == GMT_GRID_PIXEL_REG)
			header->wesn[XLO] = dummy[0] - header->inc[GMT_X] / 2.0,
			header->wesn[XHI] = dummy[1] + header->inc[GMT_X] / 2.0;
		else	/* Use as is */
			header->wesn[XLO] = dummy[0], header->wesn[XHI] = dummy[1];

in which dummy[0] and dummy[1] is the min/max value of the x vector. For pixel registration, GMT extends the x boundaries by half a grid. The codes only make sense when head->inc[GMT_X] is positive.

Copy link
Member

@weiji14 weiji14 Apr 15, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Ok, from what I can tell, this means GMT expects grids with coordinates sorted in ascending order. But here in the to_dataarray() method, we're going from GMT virtualfile to xarray.DataArray, so is this sorting still needed or can we return the grid as is? Maybe the sort ascending logic is only needed in virtualfile_from_grid (which is implemented in the dataarray_to_matrix conversion function already)?

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

But here in the to_dataarray() method, we're going from GMT virtualfile to xarray.DataArray, so is this sorting still needed or can we return the grid as is? Maybe the sort ascending logic is only needed in virtualfile_from_grid (which is implemented in the dataarray_to_matrix conversion function already)?

Hmmm, your comments make more sense. I tried to remove these sorting codes locally. Here is the output:

In [1]: import xarray as xr

In [2]: from pygmt import which

In [3]: path = which("@earth_relief_01d_g", download="a")

In [4]: grid1 = xr.load_dataarray(path)

In [6]: from pygmt.clib import Session
   ...: with Session() as lib:
   ...:     with lib.virtualfile_out(kind="grid") as voutgrd:
   ...:         lib.call_module("read", f"{path} {voutgrd} -Tg")
   ...:         grid = lib.read_virtualfile(voutgrd, kind="grid")
   ...:         grid2 = grid.contents.to_dataarray()
   ...: 

In [7]: grid1
Out[7]: 
<xarray.DataArray 'z' (lat: 181, lon: 361)> Size: 523kB
array([[ 2865. ,  2865. ,  2865. , ...,  2865. ,  2865. ,  2865. ],
       [ 3088. ,  3087.5,  3087. , ...,  3088.5,  3088. ,  3088. ],
       [ 3100.5,  3100.5,  3101. , ...,  3101.5,  3101. ,  3100.5],
       ...,
       [-3745.5, -3729. , -3722.5, ..., -3734. , -3742. , -3745.5],
       [-2940. , -2945. , -2951. , ..., -2895.5, -2921.5, -2940. ],
       [-3861. , -3861. , -3861. , ..., -3861. , -3861. , -3861. ]])
Coordinates:
  * lon      (lon) float64 3kB -180.0 -179.0 -178.0 -177.0 ... 178.0 179.0 180.0
  * lat      (lat) float64 1kB -90.0 -89.0 -88.0 -87.0 ... 87.0 88.0 89.0 90.0
Attributes:
    actual_range:  [-7174.  5350.]
    long_name:     elevation (m)

In [8]: grid2
Out[8]: 
<xarray.DataArray 'z' (lat: 181, lon: 361)> Size: 523kB
array([[-3861. , -3861. , -3861. , ..., -3861. , -3861. , -3861. ],
       [-2940. , -2945. , -2951. , ..., -2895.5, -2921.5, -2940. ],
       [-3745.5, -3729. , -3722.5, ..., -3734. , -3742. , -3745.5],
       ...,
       [ 3100.5,  3100.5,  3101. , ...,  3101.5,  3101. ,  3100.5],
       [ 3088. ,  3087.5,  3087. , ...,  3088.5,  3088. ,  3088. ],
       [ 2865. ,  2865. ,  2865. , ...,  2865. ,  2865. ,  2865. ]])
Coordinates:
  * lat      (lat) float64 1kB 90.0 89.0 88.0 87.0 ... -87.0 -88.0 -89.0 -90.0
  * lon      (lon) float64 3kB -180.0 -179.0 -178.0 -177.0 ... 178.0 179.0 180.0
Attributes:
    Conventions:   CF-1.7
    title:         SRTM15 Earth Relief v2.5.5 at 01 arc degree
    history:       
    description:   Reduced by Gaussian Cartesian filtering (314.5 km fullwidt...
    long_name:     elevation (m)
    actual_range:  [-7174.  5350.]

It seems xr.load_dataarray returns an xarray.DataArray with ascending lon/lat coodinates, but now GMT returns descending latitude. I think I was trying to make GMT_GRID.to_dataaray() return the exactly same xarray.DataArray object as xr.load_dataarray. But, maybe there are some misunderstandings here.

Copy link
Member

@weiji14 weiji14 Apr 16, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It seems xr.load_dataarray returns an xarray.DataArray with ascending lon/lat coodinates, but now GMT returns descending latitude. I think I was trying to make GMT_GRID.to_dataaray() return the exactly same xarray.DataArray object as xr.load_dataarray. But, maybe there are some misunderstandings here.

Hmm yeah, I don't think we need to match the behaviour of xr.load_dataarray here. I tried with rioxarray.open_rasterio and it gave descending latitude:

import pygmt
import rioxarray
import xarray as xr

path = pygmt.which("@earth_relief_01d_g", download="a")
grid2 = rioxarray.open_rasterio(path)
print(grid2)
# <xarray.DataArray 'z' (band: 1, y: 181, x: 361)> Size: 131kB
# [65341 values with dtype=int16]
# Coordinates:
#   * band         (band) int64 8B 1
#   * x            (x) float64 3kB -180.0 -179.0 -178.0 ... 178.0 179.0 180.0
#   * y            (y) float64 1kB 90.0 89.0 88.0 87.0 ... -87.0 -88.0 -89.0 -90.0
#     spatial_ref  int64 8B 0
# Attributes: (12/20)
#     lat#actual_range:   [-90.  90.]
#     lat#axis:           Y
#     lat#long_name:      latitude
#     lat#standard_name:  latitude
#     lat#units:          degrees_north
#     lon#actual_range:   [-180.  180.]
#     ...                 ...
#     title:              SRTM15 Earth Relief v2.5.5 at 01 arc degree
#     actual_range:       [-7174.  5350.]
#     long_name:          elevation (m)
#     scale_factor:       0.5
#     _FillValue:         -32768
#     add_offset:         0.0

so we could probably remove the sorting in to_dataarray(), assuming this doesn't break any of the unit tests.

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I don't think we need to match the behaviour of xr.load_dataarray here.

Matching the behavior of xr.load_dataarray makes it easier to test the output of the virtualfiles using xr.DataArray.indentical.

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I don't think we need to match the behaviour of xr.load_dataarray here.

Matching the behavior of xr.load_dataarray makes it easier to test the output of the virtualfiles using xr.DataArray.identical.

The returned attributes are not the same, so xr.testing.assert_identical will fail anyway.

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Actually, we're using xr.testing.assert_allclose in the tests. If we don't sort the coordinates ascendingly, then we have to update all these tests.

Copy link
Member

@weiji14 weiji14 Apr 16, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hmm yeah, I just tried running make test locally after commenting out the sorting code, and got 19 test failures. Most of these are from the expected_grid fixture returning a DataArray with y (latitude) sorted in ascending order:

================================================================================== short test summary info ===================================================================================
FAILED ../pygmt/datatypes/grid.py::pygmt.datatypes.grid._GMT_GRID.to_dataarray
FAILED ../pygmt/tests/test_dimfilter.py::test_dimfilter_no_outgrid - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdclip.py::test_grdclip_no_outgrid - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdfill.py::test_grdfill_dataarray_out - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdfill.py::test_grdfill_asymmetric_pad - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdfilter.py::test_grdfilter_dataarray_in_dataarray_out - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdgradient.py::test_grdgradient_no_outgrid - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdhisteq.py::test_equalize_grid_no_outgrid - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdlandmask.py::test_grdlandmask_no_outgrid - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdproject.py::test_grdproject_no_outgrid[M10c] - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdproject.py::test_grdproject_no_outgrid[EPSG:3395 +width=10] - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdproject.py::test_grdproject_no_outgrid[+proj=merc +ellps=WGS84 +units=m +width=10] - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_grdsample.py::test_grdsample_dataarray_out - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_surface.py::test_surface_input_file - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_surface.py::test_surface_input_data_array - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_surface.py::test_surface_input_xyz - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_triangulate.py::test_regular_grid_no_outgrid - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_xyz2grd.py::test_xyz2grd_input_array[array] - AssertionError: Left and right DataArray objects are not close
FAILED ../pygmt/tests/test_xyz2grd.py::test_xyz2grd_input_array[Dataset] - AssertionError: Left and right DataArray objects are not close
======================================================== 19 failed, 657 passed, 4 skipped, 5 xfailed, 2 xpassed, 5 warnings in 58.94s ========================================================

Let's just keep the ascending sort for now then.


# Set the gmt accesssor.
seisman marked this conversation as resolved.
Show resolved Hide resolved
# Must put at the end. The information get lost after specific grid operation.
grid.gmt.registration = header.registration
grid.gmt.gtype = header.gtype
return grid
10 changes: 6 additions & 4 deletions pygmt/helpers/decorators.py
Original file line number Diff line number Diff line change
Expand Up @@ -267,10 +267,12 @@
- ``file`` will save the result to the file specified by the ``outfile``
parameter.""",
"outgrid": """
outgrid : str or None
Name of the output netCDF grid file. For writing a specific grid
file format or applying basic data operations to the output grid,
see :gmt-docs:`gmt.html#grd-inout-full` for the available modifiers.""",
outgrid
Name of the output netCDF grid file. If not specified, will return an
:class:`xarray.DataArray` object. For writing a specific grid file format or
applying basic data operations to the output grid, see
:gmt-docs:`gmt.html#grd-inout-full` for the available modifiers.
""",
"panel": r"""
panel : bool, int, or list
[*row,col*\|\ *index*].
Expand Down
34 changes: 14 additions & 20 deletions pygmt/src/binstats.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,21 +3,13 @@
"""

from pygmt.clib import Session
from pygmt.helpers import (
GMTTempFile,
build_arg_string,
fmt_docstring,
kwargs_to_strings,
use_alias,
)
from pygmt.io import load_dataarray
from pygmt.helpers import build_arg_string, fmt_docstring, kwargs_to_strings, use_alias


@fmt_docstring
@use_alias(
C="statistic",
E="empty",
G="outgrid",
I="spacing",
N="normalize",
R="region",
Expand All @@ -31,7 +23,7 @@
r="registration",
)
@kwargs_to_strings(I="sequence", R="sequence", i="sequence_comma")
def binstats(data, **kwargs):
def binstats(data, outgrid: str | None = None, **kwargs):
r"""
Bin spatial data and determine statistics per bin.

Expand Down Expand Up @@ -110,13 +102,15 @@ def binstats(data, **kwargs):
- None if ``outgrid`` is set (grid output will be stored in file set by
``outgrid``)
"""
with GMTTempFile(suffix=".nc") as tmpfile:
with Session() as lib:
with lib.virtualfile_in(check_kind="vector", data=data) as vintbl:
if (outgrid := kwargs.get("G")) is None:
kwargs["G"] = outgrid = tmpfile.name # output to tmpfile
lib.call_module(
module="binstats", args=build_arg_string(kwargs, infile=vintbl)
)

return load_dataarray(outgrid) if outgrid == tmpfile.name else None
with Session() as lib:
with (
lib.virtualfile_in(check_kind="vector", data=data) as vintbl,
lib.virtualfile_out(kind="grid", fname=outgrid) as voutgrd,
):
kwargs["G"] = voutgrd
lib.call_module(
module="binstats", args=build_arg_string(kwargs, infile=vintbl)
)
return lib.virtualfile_to_raster(
vfname=voutgrd, kind="grid", outgrid=outgrid
)
Loading
Loading