-
Notifications
You must be signed in to change notification settings - Fork 222
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
Add function to load Blue Marble dataset #2235
Changes from 34 commits
9bde7f7
d8c8759
c6bd235
093fc26
2f29458
1c99f21
e839ef7
aa8a474
cd4cab6
fc060ea
c3017d5
c3e2e66
a2a18a9
7c38754
5c65097
767d7c7
a8600d3
971ec11
6eabbf4
4ef65b0
54184fb
e6c2964
19a7755
01d1b80
0ca210a
28d668d
a6a08e3
f13d12b
471c6e8
e886408
5b6bb42
755275e
f7e381e
9a5e36b
a0fa29b
70f6041
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,97 @@ | ||
""" | ||
Function to download the NASA Blue Marble image datasets from the GMT data server, and | ||
load as :class:`xarray.DataArray`. | ||
|
||
The images are available in various resolutions. | ||
""" | ||
|
||
from collections.abc import Sequence | ||
from typing import Literal | ||
|
||
import xarray as xr | ||
from pygmt.datasets.load_remote_dataset import _load_remote_dataset | ||
|
||
__doctest_skip__ = ["load_blue_marble"] | ||
|
||
|
||
def load_blue_marble( | ||
resolution: Literal[ | ||
"01d", | ||
"30m", | ||
"20m", | ||
"15m", | ||
"10m", | ||
"06m", | ||
"05m", | ||
"04m", | ||
"03m", | ||
"02m", | ||
"01m", | ||
"30s", | ||
] = "01d", | ||
region: Sequence[float] | str | None = None, | ||
) -> xr.DataArray: | ||
r""" | ||
Load NASA Blue Marble images in various resolutions. | ||
|
||
.. figure:: https://www.generic-mapping-tools.org/remote-datasets/_images/GMT_earth_daynight.jpg | ||
:width: 80% | ||
:align: center | ||
|
||
Earth day/night dataset. | ||
|
||
The images are downloaded to a user data directory (usually | ||
``~/.gmt/server/earth/earth_day/``) the first time you invoke this function. | ||
Afterwards, it will load the image from the data directory. So you'll need an | ||
internet connection the first time around. | ||
|
||
These images can also be accessed by passing in the file name | ||
**@earth_day**\_\ *res* to any image processing function or plotting method. *res* | ||
is the image resolution (see below). | ||
|
||
Refer to :gmt-datasets:`earth-daynight.html` for more details about available | ||
datasets, including version information and references. | ||
|
||
Parameters | ||
---------- | ||
resolution | ||
The image resolution. The suffix ``d``, ``m``, and ``s`` stand for arc-degree, | ||
arc-minute, and arc-second. | ||
|
||
region | ||
The subregion of the image to load, in the form of a sequence [*xmin*, *xmax*, | ||
*ymin*, *ymax*]. | ||
|
||
Returns | ||
------- | ||
image | ||
The NASA Blue Marble image. Coordinates are latitude and longitude in degrees. | ||
|
||
Note | ||
---- | ||
The registration and coordinate system type of the returned | ||
:class:`xarray.DataArray` image can be accessed via the GMT accessors (i.e., | ||
``image.gmt.registration`` and ``image.gmt.gtype`` respectively). However, these | ||
properties may be lost after specific image operations (such as slicing) and will | ||
need to be manually set before passing the image to any PyGMT data processing or | ||
plotting functions. Refer to :class:`pygmt.GMTDataArrayAccessor` for detailed | ||
explanations and workarounds. | ||
|
||
Examples | ||
-------- | ||
|
||
>>> from pygmt.datasets import load_blue_marble | ||
>>> # load the default image (pixel-registered 1 arc-degree image) | ||
>>> image = load_blue_marble() | ||
""" | ||
image = _load_remote_dataset( | ||
name="earth_day", | ||
prefix="earth_day", | ||
resolution=resolution, | ||
region=region, | ||
registration="pixel", | ||
) | ||
# If rioxarray is installed, set the coordinate reference system | ||
if hasattr(image, "rio"): | ||
image = image.rio.write_crs(input_crs="OGC:CRS84") | ||
return image |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,41 @@ | ||
""" | ||
Test basic functionality for loading Blue Marble datasets. | ||
""" | ||
|
||
import numpy as np | ||
import numpy.testing as npt | ||
from pygmt.datasets import load_blue_marble | ||
|
||
|
||
def test_blue_marble_01d(): | ||
""" | ||
Test some properties of the Blue Marble 01d data. | ||
""" | ||
data = load_blue_marble(resolution="01d") | ||
assert data.name == "z" | ||
assert data.long_name == "blue_marble" | ||
assert data.attrs["horizontal_datum"] == "WGS84" | ||
assert data.attrs["description"] == "NASA Day Images" | ||
assert data.shape == (3, 180, 360) | ||
assert data.dtype == "uint8" | ||
assert data.gmt.registration == 1 | ||
assert data.gmt.gtype == 1 | ||
npt.assert_allclose(data.y, np.arange(89.5, -90.5, -1)) | ||
npt.assert_allclose(data.x, np.arange(-179.5, 180.5, 1)) | ||
npt.assert_allclose(data.min(), 10, atol=1) | ||
npt.assert_allclose(data.max(), 255, atol=1) | ||
|
||
|
||
def test_blue_marble_01d_with_region(): | ||
""" | ||
Test loading low-resolution Blue Marble with 'region'. | ||
""" | ||
data = load_blue_marble(resolution="01d", region=[-10, 10, -5, 5]) | ||
assert data.shape == (3, 10, 20) | ||
assert data.dtype == "uint8" | ||
assert data.gmt.registration == 1 | ||
assert data.gmt.gtype == 1 | ||
npt.assert_allclose(data.y, np.arange(4.5, -5.5, -1)) | ||
npt.assert_allclose(data.x, np.arange(-9.5, 10.5, 1)) | ||
npt.assert_allclose(data.min(), 10, atol=1) | ||
npt.assert_allclose(data.max(), 77, atol=1) |
Original file line number | Diff line number | Diff line change | ||
---|---|---|---|---|
|
@@ -3,7 +3,8 @@ | |||
""" | ||||
|
||||
import pytest | ||||
from pygmt import Figure, which | ||||
from pygmt import Figure | ||||
from pygmt.datasets import load_blue_marble | ||||
|
||||
rioxarray = pytest.importorskip("rioxarray") | ||||
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. Can this line be removed?
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. Not yet, because we are still using |
||||
|
||||
|
@@ -14,12 +15,9 @@ def fixture_xr_image(): | |||
Load the image data from Blue Marble as an xarray.DataArray with shape {"band": 3, | ||||
"y": 180, "x": 360}. | ||||
""" | ||||
geotiff = which(fname="@earth_day_01d_p", download="c") | ||||
with rioxarray.open_rasterio(filename=geotiff) as rda: | ||||
if len(rda.band) == 3: | ||||
xr_image = rda.load() | ||||
assert xr_image.sizes == {"band": 3, "y": 180, "x": 360} | ||||
return xr_image | ||||
xr_image = load_blue_marble(resolution="01d") | ||||
assert xr_image.sizes == {"band": 3, "y": 180, "x": 360} | ||||
return xr_image | ||||
|
||||
|
||||
@pytest.mark.mpl_image_compare | ||||
|
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
In
GMT_IMAGE.to_dataarray()
, perhaps we should parse theheader->ProjRefPROJ4
and set the correct CRS to the 3-band xarray.DataArray. If done, then we probably don't need to set the CRS here.There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Yeah, I was thinking of parsing the projection information from the header when you mentioned this
ProjRefPROJ4
field at #3128 (comment). But ideally we'll need to handle PROJ4/WKT/EPSG:pygmt/pygmt/datatypes/header.py
Lines 138 to 143 in f97c3a4
Something to consider for a separate PR, because we might need to use
pyproj
for this.