Skip to content

Commit 74d8318

Browse files
snowman2fmaussion
authored andcommitted
rasterio backend: added nodatavals attribute (#1740)
* Added nodatavals attribute Connected with issue #1736 * Replace None nodatavals with np.nan Fixes the type error with serialization: https://travis-ci.org/pydata/xarray/jobs/306107679 * Fix typo * Added nodatavals open_rasterio information. * added issue and author information * added tests for nodatavals attr * moved back to nan * added separate test for missing nodataval
1 parent f2ea7b6 commit 74d8318

File tree

3 files changed

+68
-0
lines changed

3 files changed

+68
-0
lines changed

doc/whats-new.rst

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -27,6 +27,8 @@ Documentation
2727

2828
Enhancements
2929
~~~~~~~~~~~~
30+
- Added nodatavals attribute to DataArray when using :py:func:`~xarray.open_rasterio`. (:issue:`1736`).
31+
By `Alan Snow <https://github.com/snowman2>`_.
3032

3133
- :py:func:`~plot.contourf()` learned to contour 2D variables that have both a
3234
1D co-ordinate (e.g. time) and a 2D co-ordinate (e.g. depth as a function of

xarray/backends/rasterio_.py

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -189,6 +189,10 @@ def open_rasterio(filename, chunks=None, cache=None, lock=None):
189189
# Affine transformation matrix (tuple of floats)
190190
# Describes coefficients mapping pixel coordinates to CRS
191191
attrs['transform'] = tuple(riods.transform)
192+
if hasattr(riods, 'nodatavals'):
193+
# The nodata values for the raster bands
194+
attrs['nodatavals'] = tuple([np.nan if nodataval is None else nodataval
195+
for nodataval in riods.nodatavals])
192196

193197
# Parse extra metadata from tags, if supported
194198
parsers = {'ENVI': _parse_envi}

xarray/tests/test_backends.py

Lines changed: 62 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -2170,6 +2170,68 @@ def test_serialization(self):
21702170
with xr.open_dataarray(tmp_nc_file) as ncds:
21712171
assert_identical(rioda, ncds)
21722172

2173+
@requires_scipy_or_netCDF4
2174+
def test_nodata(self):
2175+
import rasterio
2176+
from rasterio.transform import from_origin
2177+
2178+
# Create a geotiff file in utm proj
2179+
with create_tmp_file(suffix='.tif') as tmp_file:
2180+
# data
2181+
nx, ny, nz = 4, 3, 3
2182+
data = np.arange(nx*ny*nz,
2183+
dtype=rasterio.float32).reshape(nz, ny, nx)
2184+
transform = from_origin(5000, 80000, 1000, 2000.)
2185+
with rasterio.open(
2186+
tmp_file, 'w',
2187+
driver='GTiff', height=ny, width=nx, count=nz,
2188+
crs={'units': 'm', 'no_defs': True, 'ellps': 'WGS84',
2189+
'proj': 'utm', 'zone': 18},
2190+
transform=transform,
2191+
nodata=-9765,
2192+
dtype=rasterio.float32) as s:
2193+
s.write(data)
2194+
expected_nodatavals = [-9765, -9765, -9765]
2195+
with xr.open_rasterio(tmp_file) as rioda:
2196+
np.testing.assert_array_equal(rioda.attrs['nodatavals'],
2197+
expected_nodatavals)
2198+
with create_tmp_file(suffix='.nc') as tmp_nc_file:
2199+
rioda.to_netcdf(tmp_nc_file)
2200+
with xr.open_dataarray(tmp_nc_file) as ncds:
2201+
np.testing.assert_array_equal(ncds.attrs['nodatavals'],
2202+
expected_nodatavals)
2203+
2204+
@requires_scipy_or_netCDF4
2205+
def test_nodata_missing(self):
2206+
import rasterio
2207+
from rasterio.transform import from_origin
2208+
2209+
# Create a geotiff file in utm proj
2210+
with create_tmp_file(suffix='.tif') as tmp_file:
2211+
# data
2212+
nx, ny, nz = 4, 3, 3
2213+
data = np.arange(nx*ny*nz,
2214+
dtype=rasterio.float32).reshape(nz, ny, nx)
2215+
transform = from_origin(5000, 80000, 1000, 2000.)
2216+
with rasterio.open(
2217+
tmp_file, 'w',
2218+
driver='GTiff', height=ny, width=nx, count=nz,
2219+
crs={'units': 'm', 'no_defs': True, 'ellps': 'WGS84',
2220+
'proj': 'utm', 'zone': 18},
2221+
transform=transform,
2222+
dtype=rasterio.float32) as s:
2223+
s.write(data)
2224+
2225+
expected_nodatavals = [np.nan, np.nan, np.nan]
2226+
with xr.open_rasterio(tmp_file) as rioda:
2227+
np.testing.assert_array_equal(rioda.attrs['nodatavals'],
2228+
expected_nodatavals)
2229+
with create_tmp_file(suffix='.nc') as tmp_nc_file:
2230+
rioda.to_netcdf(tmp_nc_file)
2231+
with xr.open_dataarray(tmp_nc_file) as ncds:
2232+
np.testing.assert_array_equal(ncds.attrs['nodatavals'],
2233+
expected_nodatavals)
2234+
21732235
def test_utm(self):
21742236
import rasterio
21752237
from rasterio.transform import from_origin

0 commit comments

Comments
 (0)