mirror of
https://github.com/rasterio/rasterio.git
synced 2026-02-01 14:34:43 +00:00
Add support for dataset and band tags (#32)
GDAL calls these items "metadata" but I think "tags" suits rasterio better.
This commit is contained in:
parent
e0b93e02c6
commit
6b29b2e15c
@ -1,6 +1,10 @@
|
||||
Changes
|
||||
=======
|
||||
|
||||
Next
|
||||
------------------
|
||||
- Add support for dataset and band tags (#32).
|
||||
|
||||
0.5.1 (2014-02-02)
|
||||
------------------
|
||||
- Add mask option to shapes() function (#26).
|
||||
|
||||
@ -6,6 +6,9 @@ cdef extern from "cpl_conv.h":
|
||||
void CPLSetThreadLocalConfigOption (char *key, char *val)
|
||||
|
||||
cdef extern from "cpl_string.h":
|
||||
int CSLCount (char **papszStrList)
|
||||
char ** CSLAddNameValue (char **papszStrList, const char *pszName, const char *pszValue)
|
||||
int CSLFindName (char **papszStrList, const char *pszName)
|
||||
char ** CSLSetNameValue (char **list, char *name, char *val)
|
||||
void CSLDestroy (char **list)
|
||||
|
||||
@ -48,6 +51,11 @@ cdef extern from "gdal.h":
|
||||
const char * GDALGetDriverShortName(void *driver)
|
||||
const char * GDALGetDriverLongName(void *driver)
|
||||
|
||||
char** GDALGetMetadata (void *hObject, const char *pszDomain)
|
||||
int GDALSetMetadata (void *hObject, char **papszMD, const char *pszDomain)
|
||||
const char* GDALGetMetadataItem(void *hObject, const char *pszName, const char *pszDomain)
|
||||
int GDALSetMetadataItem (void *hObject, const char *pszName, const char *pszValue, const char *pszDomain)
|
||||
|
||||
cdef extern from "gdal_alg.h":
|
||||
int GDALPolygonize(void *src_band, void *mask_band, void *layer, int fidx, char **options, void *progress_func, void *progress_data)
|
||||
int GDALSieveFilter(void *src_band, void *mask_band, void *dst_band, int size, int connectivity, char **options, void *progress_func, void *progress_data)
|
||||
|
||||
@ -546,6 +546,41 @@ cdef class RasterReader:
|
||||
|
||||
return out
|
||||
|
||||
def tags(self, bidx=0, domain=None):
|
||||
"""Get the dataset's tags.
|
||||
|
||||
The optional domain and band arguments select a specific
|
||||
GDAL metadata domain and tags for a specific band.
|
||||
"""
|
||||
cdef char *item_c
|
||||
cdef void *hobj
|
||||
cdef const char *domain_c
|
||||
cdef char **papszStrList
|
||||
|
||||
if not self._hds:
|
||||
raise ValueError("can't read closed raster file")
|
||||
if bidx > 0:
|
||||
hobj = _gdal.GDALGetRasterBand(self._hds, bidx)
|
||||
if hobj is NULL:
|
||||
raise ValueError("NULL band")
|
||||
else:
|
||||
hobj = self._hds
|
||||
if domain:
|
||||
domain_b = domain.encode('utf-8')
|
||||
domain_c = domain_b
|
||||
else:
|
||||
domain_c = NULL
|
||||
papszStrList = _gdal.GDALGetMetadata(hobj, domain_c)
|
||||
num_items = _gdal.CSLCount(papszStrList)
|
||||
retval = {}
|
||||
for i in range(num_items):
|
||||
item_c = papszStrList[i]
|
||||
item_b = item_c
|
||||
item = item_b.decode('utf-8')
|
||||
key, value = item.split('=')
|
||||
retval[key] = value
|
||||
return retval
|
||||
|
||||
|
||||
cdef class RasterUpdater(RasterReader):
|
||||
# Read-write access to raster data and metadata.
|
||||
@ -774,3 +809,42 @@ cdef class RasterUpdater(RasterReader):
|
||||
raise ValueError("Invalid dtype")
|
||||
# TODO: handle errors (by retval).
|
||||
|
||||
def update_tags(self, bidx=0, domain=None, **kwargs):
|
||||
"""Update the dataset's tags from items in kwargs.
|
||||
|
||||
The optional domain and band arguments select a specific
|
||||
GDAL metadata domain and tags for a specific band.
|
||||
|
||||
Only keys and values that can be converted to strings are
|
||||
allowed.
|
||||
"""
|
||||
cdef char *key_c, *value_c
|
||||
cdef void *hobj
|
||||
cdef const char *domain_c
|
||||
|
||||
if not self._hds:
|
||||
raise ValueError("can't read closed raster file")
|
||||
if bidx > 0:
|
||||
hobj = _gdal.GDALGetRasterBand(self._hds, bidx)
|
||||
if hobj is NULL:
|
||||
raise ValueError("NULL band")
|
||||
else:
|
||||
hobj = self._hds
|
||||
if domain:
|
||||
domain_b = domain.encode('utf-8')
|
||||
domain_c = domain_b
|
||||
else:
|
||||
domain_c = NULL
|
||||
cdef char **papszStrList = _gdal.GDALGetMetadata(hobj, NULL)
|
||||
for key, value in kwargs.items():
|
||||
key_b = str(key).encode('utf-8')
|
||||
value_b = str(value).encode('utf-8')
|
||||
key_c = key_b
|
||||
value_c = value_b
|
||||
i = _gdal.CSLFindName(papszStrList, key_c)
|
||||
if i < 0:
|
||||
papszStrList = _gdal.CSLAddNameValue(papszStrList, key_c, value_c)
|
||||
else:
|
||||
papszStrList = _gdal.CSLSetNameValue(papszStrList, key_c, value_c)
|
||||
retval = _gdal.GDALSetMetadata(hobj, papszStrList, domain_c)
|
||||
|
||||
|
||||
29
rasterio/tests/test_tags.py
Normal file
29
rasterio/tests/test_tags.py
Normal file
@ -0,0 +1,29 @@
|
||||
|
||||
import pytest
|
||||
import rasterio
|
||||
|
||||
def test_tags_read():
|
||||
with rasterio.open('rasterio/tests/data/RGB.byte.tif') as src:
|
||||
assert src.tags() == {'AREA_OR_POINT': 'Area'}
|
||||
assert src.tags(domain='IMAGE_STRUCTURE') == {'INTERLEAVE': 'PIXEL'}
|
||||
assert 'STATISTICS_MAXIMUM' in src.tags(1)
|
||||
|
||||
def test_tags_update(tmpdir):
|
||||
tiffname = str(tmpdir.join('foo.tif'))
|
||||
with rasterio.open(
|
||||
tiffname,
|
||||
'w',
|
||||
driver='GTiff',
|
||||
count=1,
|
||||
dtype=rasterio.uint8,
|
||||
width=10,
|
||||
height=10) as dst:
|
||||
dst.update_tags(a='1', b='2')
|
||||
dst.update_tags(1, c=3)
|
||||
assert dst.tags() == {'a': '1', 'b': '2'}
|
||||
assert dst.tags(1) == {'c': '3'}
|
||||
|
||||
with rasterio.open(tiffname) as src:
|
||||
assert src.tags() == {'a': '1', 'b': '2'}
|
||||
assert src.tags(1) == {'c': '3'}
|
||||
|
||||
Loading…
x
Reference in New Issue
Block a user