Skip to content

Commit 274bd4b

Browse files
dtpcdcherian
andauthored
Fix open_rasterio() for WarpedVRT with specified src_crs (#4104)
* Test open_rasterio() support of WarpedVRT with specified src_crs * Pass additional WarpedVRT params when recreating in open_rasterio() * Add description to `whats-new.rst` * Update doc/whats-new.rst Co-authored-by: Deepak Cherian <[email protected]>
1 parent 09df5ca commit 274bd4b

File tree

3 files changed

+20
-2
lines changed

3 files changed

+20
-2
lines changed

doc/whats-new.rst

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -143,6 +143,8 @@ Bug fixes
143143
By `Mathias Hauser <https://github.com/mathause>`_.
144144
- Fix html repr in untrusted notebooks: fallback to plain text repr. (:pull:`4053`)
145145
By `Benoit Bovy <https://github.com/benbovy>`_.
146+
- Fix :py:func:`open_rasterio` for ``WarpedVRT`` with specified ``src_crs``. (:pull:`4104`)
147+
By `Dave Cole <https://github.com/dtpc>`_.
146148

147149
Documentation
148150
~~~~~~~~~~~~~

xarray/backends/rasterio_.py

Lines changed: 5 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -221,14 +221,17 @@ def open_rasterio(filename, parse_coordinates=None, chunks=None, cache=None, loc
221221
vrt = filename
222222
filename = vrt.src_dataset.name
223223
vrt_params = dict(
224+
src_crs=vrt.src_crs.to_string(),
224225
crs=vrt.crs.to_string(),
225226
resampling=vrt.resampling,
227+
tolerance=vrt.tolerance,
226228
src_nodata=vrt.src_nodata,
227229
nodata=vrt.nodata,
228-
tolerance=vrt.tolerance,
229-
transform=vrt.transform,
230230
width=vrt.width,
231231
height=vrt.height,
232+
src_transform=vrt.src_transform,
233+
transform=vrt.transform,
234+
dtype=vrt.working_dtype,
232235
warp_extras=vrt.warp_extras,
233236
)
234237

xarray/tests/test_backends.py

Lines changed: 13 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -4160,6 +4160,19 @@ def test_rasterio_vrt_with_transform_and_size(self):
41604160
assert actual_shape == expected_shape
41614161
assert actual_transform == expected_transform
41624162

4163+
def test_rasterio_vrt_with_src_crs(self):
4164+
# Test open_rasterio() support of WarpedVRT with specified src_crs
4165+
import rasterio
4166+
4167+
# create geotiff with no CRS and specify it manually
4168+
with create_tmp_geotiff(crs=None) as (tmp_file, expected):
4169+
src_crs = rasterio.crs.CRS({"init": "epsg:32618"})
4170+
with rasterio.open(tmp_file) as src:
4171+
assert src.crs is None
4172+
with rasterio.vrt.WarpedVRT(src, src_crs=src_crs) as vrt:
4173+
with xr.open_rasterio(vrt) as da:
4174+
assert da.crs == src_crs
4175+
41634176
@network
41644177
def test_rasterio_vrt_network(self):
41654178
# Make sure loading w/ rasterio give same results as xarray

0 commit comments

Comments
 (0)