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

load_tile_map: Add the new parameter 'crs' to set the CRS of the returned dataarray #3554

Merged
merged 10 commits into from
Dec 3, 2024

Conversation

seisman
Copy link
Member

@seisman seisman commented Oct 25, 2024

Description of proposed changes

This PR adds the crs parameter to the load_tile_map function so that users can change the CRS of the returned raster image without needing to know the details.

Address #2125 (comment).
Closes #3484.

Notes for maintainers:

As far as I know, there are multiple different ways to reproject raster images.

  1. contextily.warp_tiles: Internally calls rasterio.vrt.WarpedVRT
  2. rioxarray's rio.reproject: Internally calls rasterio.warp.reproject
  3. rasterio.vrt.WarpedVRT
  4. rasterio.warp.reproject
  5. gdal.Warp
  6. ...

Options 3-5 are complicated and I don't think we want to call them directly. Ideally, we should use option 1, so reprojection doesn't rely on an extra dependency rioxarray. However, as shown below, the result from contextily.warp_tiles doesn't work well with rio.to_raster.

Here are codes to compare the results with options 1 and 2:

import xarray as xr
import rioxarray
import contextily as cx
import numpy as np

def image2raster(image, extent, crs="EPSG:3857"):
    """
    Convert an image to a xarray.DataArray raster.

    Copied from pygmt.datasets.tile_map.load_tile_map, except that crs is an input parameter.
    """
    # Turn RGBA img from channel-last to channel-first and get 3-band RGB only
    _image = image.transpose(2, 0, 1)  # Change image from (H, W, C) to (C, H, W)
    rgb_image = _image[0:3, :, :]  # Get just RGB by dropping RGBA's alpha channel
    
    # Georeference RGB image into an xarray.DataArray
    left, right, bottom, top = extent
    dataarray = xr.DataArray(
        data=rgb_image,
        coords={
            "band": np.array(object=[1, 2, 3], dtype=np.uint8),  # Red, Green, Blue
            "y": np.linspace(start=top, stop=bottom, num=rgb_image.shape[1]),
            "x": np.linspace(start=left, stop=right, num=rgb_image.shape[2]),
        },
        dims=("band", "y", "x"),
    )
    dataarray = dataarray.rio.write_crs(input_crs=crs)
    return dataarray

# Image and extent of the global tile
image, extent = cx.bounds2img(w=-180, e=180, s=-90, n=90, ll=True, zoom=0)

# Option 1: Use contextily.warp_tiles
image_warp, extent_warp = cx.warp_tiles(image, extent, "OGC:CRS84")
dataarray_warp = image2raster(image_warp, extent_warp, crs="OGC:CRS84")
print("Using contextily.warp_tiles:", dataarray_warp.rio.bounds())
dataarray_warp.rio.to_raster("result_warp.tiff")

# Option 2: Use rio.reproject
dataarray = image2raster(image, extent, crs="EPSG:3857")
dataarray_reproject = dataarray.rio.reproject("OGC:CRS84")
print("Using rio.reproject:", dataarray_reproject.rio.bounds())
dataarray_reproject.rio.to_raster("result_reproject.tiff")

The outputs are:

Using contextily.warp_tiles: (-180.55157789333614, -85.9082903958547, 180.18036434850873, 85.66487761291295)
Using rio.reproject: (-180.0, -84.8307625666717, 180.00000000000006, 85.11165075221737)

As you can see, the longitude range is (-180.55157789333614, 180.18036434850873) with contextily.warp_tiles, which will lead to GMT errors like:

grdimage [ERROR]: Map region exceeds 360 degrees
grdimage [ERROR]: General map projection error

Two more notes:

Preview: https://pygmt-dev--3554.org.readthedocs.build/en/3554/api/generated/pygmt.datasets.load_tile_map.html

Reminders

  • Run make format and make check to make sure the code follows the style guide.
  • Add tests for new features or tests that would have caught the bug that you're fixing.
  • Add new public functions/methods/classes to doc/api/index.rst.
  • Write detailed docstrings for all functions/methods.
  • If wrapping a new module, open a 'Wrap new GMT module' issue and submit reasonably-sized PRs.
  • If adding new functionality, add an example to docstrings or tutorials.
  • Use underscores (not hyphens) in names of Python files and directories.

Slash Commands

You can write slash commands (/command) in the first line of a comment to perform
specific operations. Supported slash command is:

  • /format: automatically format and lint the code

@seisman seisman self-assigned this Oct 30, 2024
@seisman seisman force-pushed the load_tile_map/crs branch 3 times, most recently from c34a56f to d79af32 Compare November 27, 2024 06:01
@seisman seisman force-pushed the load_tile_map/crs branch 2 times, most recently from 259b0dd to a48f9df Compare November 27, 2024 06:26
@seisman seisman changed the title WIP: load_tile_map: Add the new parameter 'crs' to set the CRS of the returned dataarray load_tile_map: Add the new parameter 'crs' to set the CRS of the returned dataarray Nov 27, 2024
@seisman seisman marked this pull request as ready for review November 27, 2024 06:52
"""
kwargs = self._preprocess(**kwargs)

if not _HAS_RIOXARRAY:
Copy link
Member Author

Choose a reason for hiding this comment

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

The Figure.tilemap method no longer requires rioxarray explicitly. Instead, load_tile_map will raise the error.

raster = load_tile_map(
region=region,
zoom=zoom,
source=source,
lonlat=lonlat,
crs="OGC:CRS84" if lonlat is True else "EPSG:3857",
Copy link
Member Author

Choose a reason for hiding this comment

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

If lonlat=True, set the CRS to let load_tile_map do the raster reprojection for us.

@seisman seisman added enhancement Improving an existing feature needs review This PR has higher priority and needs review. labels Nov 27, 2024
@seisman seisman added this to the 0.14.0 milestone Nov 27, 2024
@seisman seisman removed the needs review This PR has higher priority and needs review. label Nov 27, 2024
@seisman seisman marked this pull request as draft November 27, 2024 14:23
@seisman seisman added the needs review This PR has higher priority and needs review. label Dec 2, 2024
@seisman seisman marked this pull request as ready for review December 2, 2024 03:37
Copy link
Member

@michaelgrund michaelgrund left a comment

Choose a reason for hiding this comment

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

Looks good to me!

Comment on lines +13 to +18
from rasterio.crs import CRS
from xyzservices import TileProvider

_HAS_CONTEXTILY = True
except ImportError:
CRS = None
Copy link
Member

Choose a reason for hiding this comment

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

Maybe put the rasterio imports in the rioxarray block below instead of here.

Copy link
Member Author

Choose a reason for hiding this comment

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

rasterio is a known dependency of contextily (xref: https://github.com/geopandas/contextily/blob/main/pyproject.toml), so it makes more sense to assume that they're installed together.

Copy link
Member

@weiji14 weiji14 Dec 3, 2024

Choose a reason for hiding this comment

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

Rasterio is also a dependency of rioxarray 😆. I would argue that rasterio is more tightly coupled to rioxarray, and since rasterio/rioxarray are both used for their projection system capabilities, it would make more sense to put those together. But ok too if you want to keep it here.

pygmt/datasets/tile_map.py Outdated Show resolved Hide resolved
@@ -128,6 +138,8 @@ def load_tile_map(
... raster.rio.crs.to_string()
'EPSG:3857'
"""
_default_crs = "EPSG:3857" # The default CRS returned from contextily
Copy link
Member

Choose a reason for hiding this comment

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

At https://contextily.readthedocs.io/en/latest/reference.html#contextily.bounds2img, there is this sentence:

IMPORTANT: tiles are assumed to be in the Spherical Mercator projection (EPSG:3857), unless the crs keyword is specified.

While tilemaps are usually served in EPSG:3857, I know that there are some tilemaps that can be served in a different projection system by default. I'm hesitant to use this default of EPSG:3857 here, unless we can be confident that contextily only supports EPSG:3857 tiles.

Copy link
Member Author

Choose a reason for hiding this comment

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

Looking at the contexily source codes, I think contexitly always assumes that CRS is EPSG:3857, at least when writing the image to raster files

https://github.com/geopandas/contextily/blob/eafeb3c993ab43356b65e7a29b9dd6df3ac76f00/contextily/tile.py#L358

Copy link
Member

Choose a reason for hiding this comment

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

That is for the warp_tiles function. I don't see anything about EPSG:3857 in the bounds2img function (though it's calling several functions and I haven't checked too closely).

Non-Web Mercator (EPSG:3857) XYZ tiles are not common, but they do exist, e.g. Openlayers allows configuring the projection system of the returned tiles:

Copy link
Member Author

Choose a reason for hiding this comment

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

Looking at the contexily source codes, I think contexitly always assumes that CRS is EPSG:3857, at least when writing the image to raster files

I meant this line:

https://github.com/geopandas/contextily/blob/eafeb3c993ab43356b65e7a29b9dd6df3ac76f00/contextily/tile.py#L175

In bounds2raster, it first calls bounds2img to get the image and extent (can be any CRS), but when writing to a raster image, bounds2raster assumes the image is EPSG:3857.

Looking at the xyzservices (https://github.com/geopandas/xyzservices/blob/c847f312d2e67a0a45ceae119fdcdf7cf60858b6/provider_sources/xyzservices-providers.json#L50), I can also see some providers that has the crs attribute, but it's unclear how the crs attribute is handled in contexitly.

Copy link
Member

Choose a reason for hiding this comment

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

I tried using xyzservices.providers.MapTiler.Basic4326 which comes in EPSG:4326:

import contextily
import xyzservices
import matplotlib.pyplot as plt

source = xyzservices.providers.MapTiler.Basic4326(
    key="BwTZdryxoypHc4BnCZ4"  # add an 'o' to the end
)
img, bounds = contextily.bounds2img(w=0, s=-80, e=150, n=80, ll=True, source=source)
# img.shape # (4096, 2048, 4)
plt.imshow(img[:, :, :3])

produces

temp

The extent doesn't seem to be correct (it should be plotting from Greenwich to Asia), possibly a bug in contextily hardcoding to EPSG:3857. But it does show that non-Web Mercator XYZ tiles are indeed possible.

Copy link
Member

Choose a reason for hiding this comment

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

Found geopandas/contextily#119, and had a closer look at the bounds2img code, specifically these lines - https://github.com/geopandas/contextily/blob/00ca0318033e69e29a164ce571d9c932a057ecc6/contextily/tile.py#L270-L272. It seems like contextily is using mercantile which only supports Spherical Mercator (EPSG:3857), see mapbox/mercantile#109 (comment). Unless contextily switches to something like https://github.com/developmentseed/morecantile in the future, then it would only 'work' with EPSG:3857.

Copy link
Member Author

@seisman seisman Dec 3, 2024

Choose a reason for hiding this comment

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

The extent doesn't seem to be correct (it should be plotting from Greenwich to Asia), possibly a bug in contextily hardcoding to EPSG:3857.

https://github.com/geopandas/contextily/blob/eafeb3c993ab43356b65e7a29b9dd6df3ac76f00/contextily/tile.py#L256

When ll=True, bounds2img calls the private _sm2ll function to convert bounds in Spherical Mercator coordinates to lon/lat. I believe "EPSG:3857" is hardcoded/assumed in many places in the contextily package.

Copy link
Member

Choose a reason for hiding this comment

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

How about checking source.crs (an attribute from an instance of the xyzservices.TileProvider class), and if it is present, use that as the default crs. If not, fallback to EPSG:3857.

Copy link
Member Author

@seisman seisman Dec 3, 2024

Choose a reason for hiding this comment

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

Sounds good. Done in 6733e19.

_default_crs = getattr(source, "crs", "EPSG:3857")

Please note that this line of code works for TileProvide, str and None. str and None don't have crs, so EPSG:3857 is returned.

Edit: renamed _default_crs to _source_crs in c1d55b0.

Please note that, currently, the returned xarray.DataArray is still "EPSG:3857" even if the source CRS is not (i.e. we don't the reprojection). This may not be ideal.

pygmt/datasets/tile_map.py Outdated Show resolved Hide resolved
@@ -128,6 +138,8 @@ def load_tile_map(
... raster.rio.crs.to_string()
'EPSG:3857'
"""
_default_crs = "EPSG:3857" # The default CRS returned from contextily
Copy link
Member

Choose a reason for hiding this comment

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

That is for the warp_tiles function. I don't see anything about EPSG:3857 in the bounds2img function (though it's calling several functions and I haven't checked too closely).

Non-Web Mercator (EPSG:3857) XYZ tiles are not common, but they do exist, e.g. Openlayers allows configuring the projection system of the returned tiles:

@weiji14 weiji14 added final review call This PR requires final review and approval from a second reviewer and removed needs review This PR has higher priority and needs review. labels Dec 3, 2024
@seisman seisman merged commit 1478fa2 into main Dec 3, 2024
22 of 23 checks passed
@seisman seisman deleted the load_tile_map/crs branch December 3, 2024 07:28
@seisman seisman removed the final review call This PR requires final review and approval from a second reviewer label Dec 3, 2024
seisman added a commit that referenced this pull request Dec 3, 2024
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
enhancement Improving an existing feature
Projects
None yet
Development

Successfully merging this pull request may close these issues.

load_map_tiles: Allow a new parameter to control the projection of the returned DataArray
3 participants