API Documentation

Raster Class

spatialist.raster.subset_tolerance = 0

this parameter can be set to increase the pixel tolerance in percent when subsetting Raster objects with the extent of other spatial objects.

Examples

Coordinates are in EPSG:32632, pixel resolution of the image to be subsetted is 90 m:
(subsetting extent)
{‘xmin’: 534093.341, ‘xmax’: 830103.341, ‘ymin’: 5030609.645, ‘ymax’: 5250929.645}
subset_tolerance = 0
{‘xmin’: 534003.341, ‘xmax’: 830103.341, ‘ymin’: 5030519.645, ‘ymax’: 5250929.645}
subset_tolerance = 0.02
{‘xmin’: 534093.341, ‘xmax’: 830103.341, ‘ymin’: 5030609.645, ‘ymax’: 5250929.645}
class spatialist.raster.Raster(filename, list_separate=True)[source]

This is intended as a raster meta information handler with options for reading and writing raster data in a convenient manner by simplifying the numerous options provided by the GDAL python binding. Several methods are provided along with this class to directly modify the raster object in memory or directly write a newly created file to disk (without modifying the raster object itself). Upon initializing a Raster object, only metadata is loaded. The actual data can be, for example, loaded to memory by calling methods matrix() or load().

Parameters:
  • filename (str or list or gdal.Dataset) – the raster file(s)/object to read
  • list_separate (bool) – treat a list of files as separate layers or otherwise as a single layer? The former is intended for single layers of a stack, the latter for tiles of a mosaic.
__getitem__(index)[source]

subset the object by slices or vector geometry. If slices are provided, one slice for each raster dimension needs to be defined. I.e., if the raster object contains several image bands, three slices are necessary. Integer slices are treated as pixel coordinates and float slices as map coordinates. If a Vector geometry is defined, it is internally projected to the raster CRS if necessary, its extent derived, and the extent converted to raster pixel slices, which are then used for subsetting.

Parameters:index (tuple of slice or Vector) – the subsetting indices to be used
Returns:a new raster object referenced through an in-memory GDAL VRT file
Return type:Raster

Examples

>>> filename = 'test'
>>> with Raster(filename) as ras:
>>>     print(ras)
class      : spatialist Raster object
dimensions : 2908, 2069, 115 (rows, cols, bands)
resolution : 20.0, -20.0 (x, y)
extent     : 713315.198, 754695.198, 4068985.595, 4127145.595 (xmin, xmax, ymin, ymax)
coord. ref.: +proj=utm +zone=29 +datum=WGS84 +units=m +no_defs
data source: test
>>>
>>>
>>> xmin = 0
>>> xmax = 100
>>> ymin = 4068985.595
>>> ymax = 4088985.595
>>> with Raster(filename)[ymin:ymax, xmin:xmax, :] as ras:
>>>     print(ras)
class      : spatialist Raster object
dimensions : 1000, 100, 115 (rows, cols, bands)
resolution : 20.0, -20.0 (x, y)
extent     : 713315.198, 715315.198, 4068985.595, 4088985.595 (xmin, xmax, ymin, ymax)
coord. ref.: +proj=utm +zone=29 +datum=WGS84 +units=m +no_defs
data source: /tmp/tmpk5weyhhq.vrt
>>>
>>>
>>> ext = {'xmin': 713315.198, 'xmax': 715315.198, 'ymin': ymin, 'ymax': ymax}
>>>
>>> with bbox(ext, crs=32629) as vec:
>>>     with Raster(filename)[vec] as ras:
>>>         print(ras)
class      : spatialist Raster object
dimensions : 1000, 100, 115 (rows, cols, bands)
resolution : 20.0, -20.0 (x, y)
extent     : 713315.198, 715315.198, 4068985.595, 4088985.595 (xmin, xmax, ymin, ymax)
coord. ref.: +proj=utm +zone=29 +datum=WGS84 +units=m +no_defs
data source: /tmp/tmps4rc9o09.vrt
allstats(approximate=False)[source]

Compute some basic raster statistics

Parameters:approximate (bool) – approximate statistics from overviews or a subset of all tiles?
Returns:a list with a dictionary of statistics for each band. Keys: min, max, mean, sdev. See gdal.Band.ComputeStatistics.
Return type:list of dicts
array()[source]

read all raster bands into a numpy ndarray

Returns:the array containing all raster data
Return type:numpy.ndarray
assign(array, band)[source]

assign an array to an existing Raster object

Parameters:
  • array (numpy.ndarray) – the array to be assigned to the Raster object
  • band (int) – the index of the band to assign to
bandnames
Returns:the names of the bands
Return type:list
bands
Returns:the number of image bands
Return type:int
bbox(outname=None, driver='ESRI Shapefile', overwrite=True, source='image')[source]
Parameters:
  • outname (str or None) – the name of the file to write; If None, the bounding box is returned as Vector object
  • driver (str) – The file format to write
  • overwrite (bool) – overwrite an already existing file?
  • source ({'image', 'gcp'}) – get the bounding box of either the image or the ground control points
Returns:

the bounding box vector object

Return type:

Vector or None

close()[source]

closes the GDAL raster file connection

cols
Returns:the number of image columns
Return type:int
coord_img2map(x=None, y=None)[source]

convert image pixel coordinates to map coordinates in the raster CRS. Either x, y or both must be defined.

Parameters:
Returns:

the converted coordinate for either x, y or both

Return type:

float or tuple

coord_map2img(x=None, y=None)[source]

convert map coordinates in the raster CRS to image pixel coordinates. Either x, y or both must be defined.

Parameters:
Returns:

the converted coordinate for either x, y or both

Return type:

int or tuple

dim
Returns:(rows, columns, bands)
Return type:tuple
driver
Returns:a GDAL raster driver object.
Return type:gdal.Driver
dtype
Returns:the data type description; e.g. Float32
Return type:str
epsg
Returns:the CRS EPSG code
Return type:int
extent
Returns:the extent of the image
Return type:dict
extract(px, py, radius=1, nodata=None)[source]

extract weighted average of pixels intersecting with a defined radius to a point.

Parameters:
  • px (int or float) – the x coordinate in units of the Raster SRS
  • py (int or float) – the y coordinate in units of the Raster SRS
  • radius (int or float) – the radius around the point to extract pixel values from; defined as multiples of the pixel resolution
  • nodata (int) – a value to ignore from the computations; If None, the nodata value of the Raster object is used
Returns:

the the weighted average of all pixels within the defined radius

Return type:

int or float

files
Returns:a list of all absolute names of files associated with this raster data set
Return type:list of str
format
Returns:the name of the image format
Return type:str
geo

General image geo information.

Returns:a dictionary with keys xmin, xmax, xres, rotation_x, ymin, ymax, yres, rotation_y
Return type:dict
geogcs
Returns:an identifier of the geographic coordinate system
Return type:str or None
is_valid()[source]

Check image integrity. Tries to compute the checksum for each raster layer and returns False if this fails. See this forum entry: How to check if image is valid?.

Returns:is the file valid?
Return type:bool
layers()[source]
Returns:a list containing a gdal.Band object for each image band
Return type:list of gdal.Band
load()[source]

load all raster data to internal memory arrays. This shortens the read time of other methods like matrix().

matrix(band=1, mask_nan=True)[source]

read a raster band (subset) into a numpy ndarray

Parameters:
  • band (int) – the band to read the matrix from; 1-based indexing
  • mask_nan (bool) – convert nodata values to numpy.nan? As numpy.nan requires at least float values, any integer array is cast to float32.
Returns:

the matrix (subset) of the selected band

Return type:

numpy.ndarray

nodata
Returns:the raster nodata value(s)
Return type:float or list
proj4
Returns:the CRS PROJ4 description
Return type:str
proj4args
Returns:the PROJ4 string arguments as a dictionary
Return type:dict
projcs
Returns:an identifier of the projected coordinate system; If the CRS is not projected None is returned
Return type:str or None
projection
Returns:the CRS Well Known Text (WKT) description
Return type:str
res

the raster resolution in x and y direction

Returns:(xres, yres)
Return type:tuple
rescale(fun)[source]

perform raster computations with custom functions and assign them to the existing raster object in memory

Parameters:fun (function) – the custom function to compute on the data

Examples

>>> with Raster('filename') as ras:
>>>     ras.rescale(lambda x: 10 * x)
rows
Returns:the number of image rows
Return type:int
srs
Returns:the spatial reference system of the data set.
Return type:osr.SpatialReference
write(outname, dtype='default', format='ENVI', nodata='default', compress_tif=False, overwrite=False, cmap=None)[source]

write the raster object to a file.

Parameters:
  • outname (str) – the file to be written
  • dtype (str) – the data type of the written file; data type notations of GDAL (e.g. Float32) and numpy (e.g. int8) are supported.
  • format (str) – the file format; e.g. ‘GTiff’
  • nodata (int or float) – the nodata value to write to the file
  • compress_tif (bool) – if the format is GeoTiff, compress the written file?
  • overwrite (bool) – overwrite an already existing file?
  • cmap (gdal.ColorTable) – a color map to apply to each band. Can for example be created with function cmap_mpl2gdal().

Raster Tools

spatialist.raster.stack(srcfiles, dstfile, resampling, targetres, dstnodata, srcnodata=None, shapefile=None, layernames=None, sortfun=None, separate=False, overwrite=False, compress=True, cores=4, pbar=False)[source]

function for mosaicking, resampling and stacking of multiple raster files into a 3D data cube

Parameters:
  • srcfiles (list) – a list of file names or a list of lists; each sub-list is treated as a task to mosaic its containing files
  • dstfile (str) – the destination file or a directory (if separate is True)
  • resampling ({near, bilinear, cubic, cubicspline, lanczos, average, mode, max, min, med, Q1, Q3}) – the resampling method; see documentation of gdalwarp.
  • targetres (tuple or list) – two entries for x and y spatial resolution in units of the source CRS
  • srcnodata (int, float or None) – the nodata value of the source files; if left at the default (None), the nodata values are read from the files
  • dstnodata (int or float) – the nodata value of the destination file(s)
  • shapefile (str, Vector or None) – a shapefile for defining the spatial extent of the destination files
  • layernames (list) – the names of the output layers; if None, the basenames of the input files are used; overrides sortfun
  • sortfun (function) – a function for sorting the input files; not used if layernames is not None. This is first used for sorting the items in each sub-list of srcfiles; the basename of the first item in a sub-list will then be used as the name for the mosaic of this group. After mosaicing, the function is again used for sorting the names in the final output (only relevant if separate is False)
  • separate (bool) – should the files be written to a single raster stack (ENVI format) or separate files (GTiff format)?
  • overwrite (bool) – overwrite the file if it already exists?
  • compress (bool) – compress the geotiff files?
  • cores (int) – the number of CPU threads to use
  • pbar (bool) – add a progressbar? This is currently only used if separate==False
Raises:

RuntimeError

Notes

This function does not reproject any raster files. Thus, the CRS must be the same for all input raster files. This is checked prior to executing gdalwarp. In case a shapefile is defined, it is internally reprojected to the raster CRS prior to retrieving its extent.

Examples

from pyroSAR.ancillary import groupbyTime, find_datasets, seconds
from spatialist.raster import stack

# find pyroSAR files by metadata attributes
archive_s1 = '/.../sentinel1/GRD/processed'
scenes_s1 = find_datasets(archive_s1, sensor=('S1A', 'S1B'), acquisition_mode='IW')

# group images by acquisition time
groups = groupbyTime(images=scenes_s1, function=seconds, time=30)

# mosaic individual groups and stack the mosaics to a single ENVI file
# only files overlapping with the shapefile are selected and resampled to its extent
stack(srcfiles=groups, dstfile='stack', resampling='bilinear', targetres=(20, 20),
      srcnodata=-99, dstnodata=-99, shapefile='site.shp', separate=False)
spatialist.raster.rasterize(vectorobject, reference, outname=None, burn_values=1, expressions=None, nodata=0, append=False)[source]

rasterize a vector object

Parameters:
  • vectorobject (Vector) – the vector object to be rasterized
  • reference (Raster) – a reference Raster object to retrieve geo information and extent from
  • outname (str or None) – the name of the GeoTiff output file; if None, an in-memory object of type Raster is returned and parameter outname is ignored
  • burn_values (int or list) – the values to be written to the raster file
  • expressions (list) – SQL expressions to filter the vector object by attributes
  • nodata (int) – the nodata value of the target raster file
  • append (bool) – if the output file already exists, update this file with new rasterized values? If True and the output file exists, parameters reference and nodata are ignored.
Returns:

if outname is None, a raster object pointing to an in-memory dataset else None

Return type:

Raster or None

Example

>>> from spatialist import Vector, Raster, rasterize
>>> outname1 = 'target1.tif'
>>> outname2 = 'target2.tif'
>>> with Vector('source.shp') as vec:
>>>     with Raster('reference.tif') as ref:
>>>         burn_values = [1, 2]
>>>         expressions = ['ATTRIBUTE=1', 'ATTRIBUTE=2']
>>>         rasterize(vec, reference, outname1, burn_values, expressions)
>>>         expressions = ["ATTRIBUTE2='a'", "ATTRIBUTE2='b'"]
>>>         rasterize(vec, reference, outname2, burn_values, expressions)

Vector Class

class spatialist.vector.Vector(filename=None, driver=None)[source]

This is intended as a vector meta information handler with options for reading and writing vector data in a convenient manner by simplifying the numerous options provided by the OGR python binding.

Parameters:
  • filename (str or None) –

    the vector file to read; if filename is None, a new in-memory Vector object is created. In this case driver is overridden and set to ‘Memory’. The following file extensions are auto-detected:

    • .geojson (GeoJSON)
    • .gpkg (GPKG)
    • .shp (ESRI Shapefile)
  • driver (str) – the vector file format; needs to be defined if the format cannot be auto-detected from the filename extension
__getitem__(expression)[source]

subset the vector object by index or attribute.

Parameters:expression (int or str) – the key or expression to be used for subsetting. See ogr.Layer.SetAttributeFilter for details on the expression syntax.
Returns:a vector object matching the specified criteria
Return type:Vector

Examples

Assuming we have a shapefile called testsites.shp, which has an attribute sitename, we can subset individual sites and write them to new files like so:

>>> from spatialist import Vector
>>> filename = 'testsites.shp'
>>> with Vector(filename)["sitename='site1'"] as site1:
>>>     site1.write('site1.shp')
addfeature(geometry, fields=None)[source]

add a feature to the vector object from a geometry

Parameters:
  • geometry (ogr.Geometry) – the geometry to add as a feature
  • fields (dict or None) – the field names and values to assign to the new feature
addfield(name, type, width=10)[source]

add a field to the vector layer

Parameters:
  • name (str) – the field name
  • type (int) – the OGR Field Type (OFT), e.g. ogr.OFTString. See Module ogr.
  • width (int) – the width of the new field (only for ogr.OFTString fields)
addlayer(name, srs, geomType)[source]

add a layer to the vector layer

Parameters:
addvector(vec)[source]

add a vector object to the layer of the current Vector object

Parameters:
  • vec (Vector) – the vector object to add
  • merge (bool) – merge overlapping polygons?
bbox(outname=None, driver=None, overwrite=True)[source]

create a bounding box from the extent of the Vector object

Parameters:
  • outname (str or None) – the name of the vector file to be written; if None, a Vector object is returned
  • driver (str) – the name of the file format to write
  • overwrite (bool) – overwrite an already existing file?
Returns:

if outname is None, the bounding box Vector object

Return type:

Vector or None

close()[source]

closes the OGR vector file connection

convert2wkt(set3D=True)[source]

export the geometry of each feature as a wkt string

Parameters:set3D (bool) – keep the third (height) dimension?
extent

the extent of the vector object

Returns:a dictionary with keys xmin, xmax, ymin, ymax
Return type:dict
fieldDefs
Returns:the field definition for each field of the Vector object
Return type:list of ogr.FieldDefn
fieldnames
Returns:the names of the fields
Return type:list of str
geomType
Returns:the layer geometry type
Return type:int
geomTypes
Returns:the geometry type of each feature
Return type:list
getArea()[source]
Returns:the area of the vector geometries
Return type:float
getFeatureByAttribute(fieldname, attribute)[source]

get features by field attribute

Parameters:
  • fieldname (str) – the name of the queried field
  • attribute (int or str) – the field value of interest
Returns:

the feature(s) matching the search query

Return type:

list of ogr.Feature or ogr.Feature

getFeatureByIndex(index)[source]

get features by numerical (positional) index

Parameters:index (int) – the queried index
Returns:the requested feature
Return type:ogr.Feature
getProjection(type)[source]

get the CRS of the Vector object. See spatialist.auxil.crsConvert().

Parameters:type (str) – the type of projection required.
Returns:the output CRS
Return type:int, str or osr.SpatialReference
getUniqueAttributes(fieldname)[source]
Parameters:fieldname (str) – the name of the field of interest
Returns:the unique attributes of the field
Return type:list of str or int
getfeatures()[source]
Returns:a list of cloned features
Return type:list of ogr.Feature
init_features()[source]

delete all in-memory features

init_layer()[source]

initialize a layer object

layerdef
Returns:the layer’s feature definition
Return type:ogr.FeatureDefn
layername
Returns:the name of the layer
Return type:str
load()[source]

load all feature into memory

nfeatures
Returns:the number of features
Return type:int
nfields
Returns:the number of fields
Return type:int
nlayers
Returns:the number of layers
Return type:int
proj4
Returns:the CRS in PRO4 format
Return type:str
reproject(projection)[source]

in-memory reprojection

Parameters:projection (int, str, osr.SpatialReference) – the target CRS. See spatialist.auxil.crsConvert().
setCRS(crs)[source]

directly reset the spatial reference system of the vector object. This is not going to reproject the Vector object, see reproject() instead.

Parameters:crs (int, str, osr.SpatialReference) – the input CRS

Example

>>> site = Vector('shape.shp')
>>> site.setCRS('+proj=longlat +ellps=WGS84 +datum=WGS84 +no_defs ')
srs
Returns:the geometry’s spatial reference system
Return type:osr.SpatialReference
write(outfile, driver=None, overwrite=True)[source]

write the Vector object to a file

Parameters:
  • outfile

    the name of the file to write; the following extensions are automatically detected for determining the format driver:

    • .geojson (GeoJSON)
    • .gpkg (GPKG)
    • .shp (ESRI Shapefile)
  • driver (str) – the output file format; needs to be defined if the format cannot be auto-detected from the filename extension
  • overwrite (bool) – overwrite an already existing file?

Vector Tools

spatialist.vector.intersect(obj1, obj2)[source]

intersect two Vector objects

Parameters:
  • obj1 (Vector) – the first vector object; this object is reprojected to the CRS of obj2 if necessary
  • obj2 (Vector) – the second vector object
Returns:

the intersect of obj1 and obj2 if both intersect and None otherwise

Return type:

Vector or None

spatialist.vector.bbox(coordinates, crs, outname=None, driver=None, overwrite=True)[source]

create a bounding box vector object or shapefile from coordinates and coordinate reference system. The CRS can be in either WKT, EPSG or PROJ4 format

Parameters:
  • coordinates (dict) – a dictionary containing numerical variables with keys xmin, xmax, ymin and ymax
  • crs (int, str, osr.SpatialReference) – the CRS of the coordinates. See crsConvert() for options.
  • outname (str) – the file to write to. If None, the bounding box is returned as Vector object
  • driver (str) –
    the output file format; needs to be defined if the format cannot
    be auto-detected from the filename extension
  • overwrite (bool) – overwrite an existing file?
Returns:

the bounding box Vector object

Return type:

Vector or None

spatialist.vector.feature2vector(feature, ref, layername=None)[source]

create a Vector object from ogr features

Parameters:
  • feature (list of ogr.Feature or ogr.Feature) – a single feature or a list of features
  • ref (Vector) – a reference Vector object to retrieve geo information from
  • layername (str or None) – the name of the output layer; retrieved from ref if None
Returns:

the new Vector object

Return type:

Vector

spatialist.vector.wkt2vector(wkt, srs, layername='wkt')[source]

convert a well-known text string geometry to a Vector object

Parameters:
Returns:

the vector representation

Return type:

Vector

Examples

>>> from spatialist.vector import wkt2vector
>>> wkt = 'POLYGON ((0. 0., 0. 1., 1. 1., 1. 0., 0. 0.))'
>>> with wkt2vector(wkt, srs=4326) as vec:
>>>     print(vec.getArea())
1.0

General Spatial Tools

spatialist.auxil.cmap_mpl2gdal(mplcolor, values)[source]

convert a matplotlib color table to a GDAL representation.

Parameters:
  • mplcolor (str) – a color table code
  • values (list) – the integer data values for which to retrieve colors
Returns:

the color table in GDAL format

Return type:

gdal.ColorTable

Notes

This function is currently only developed for handling discrete integer data values in an 8 Bit file. Colors are thus scaled between 0 and 255.

Examples

>>> from osgeo import gdal
>>> from spatialist.auxil import cmap_mpl2gdal
>>> values = list(range(0, 100))
>>> cmap = cmap_mpl2gdal(mplcolor='YlGnBu', values=values)
>>> print(isinstance(cmap, gdal.ColorTable))
True
spatialist.auxil.coordinate_reproject(x, y, s_crs, t_crs)[source]

reproject a coordinate from one CRS to another

Parameters:
Returns:

Return type:

tuple

spatialist.auxil.crsConvert(crsIn, crsOut)[source]

convert between different types of spatial references

Parameters:
  • crsIn (int, str, osr.SpatialReference) – the input CRS
  • crsOut ({'wkt', 'proj4', 'epsg', 'osr', 'opengis' or 'prettyWkt'}) – the output CRS type
Returns:

the output CRS

Return type:

int, str, osr.SpatialReference

Examples

convert an integer EPSG code to PROJ4:

>>> crsConvert(4326, 'proj4')
'+proj=longlat +datum=WGS84 +no_defs '

convert a PROJ4 string to an opengis URL:

>>> crsConvert('+proj=longlat +datum=WGS84 +no_defs ', 'opengis')
'http://www.opengis.net/def/crs/EPSG/0/4326'

convert the opengis URL back to EPSG:

>>> crsConvert('http://www.opengis.net/def/crs/EPSG/0/4326', 'epsg')
4326

convert an EPSG compound CRS (WGS84 horizontal + EGM96 vertical)

>>> crsConvert('EPSG:4326+5773', 'proj4')
'+proj=longlat +datum=WGS84 +geoidgrids=egm96_15.gtx +vunits=m +no_defs '
spatialist.auxil.gdal_rasterize(src, dst, options)[source]

a simple wrapper for gdal.Rasterize

Parameters:
spatialist.auxil.gdal_translate(src, dst, options)[source]

a simple wrapper for gdal.Translate

Parameters:
spatialist.auxil.gdalbuildvrt(src, dst, options=None, void=True)[source]

a simple wrapper for gdal.BuildVRT

Parameters:
  • src (str, list, ogr.DataSource or gdal.Dataset) – the input data set(s)
  • dst (str) – the output data set
  • options (dict) – additional parameters passed to gdal.BuildVRT; see gdal.BuildVRTOptions
  • void (bool) – just write the results and don’t return anything? If not, the spatial object is returned
spatialist.auxil.gdalwarp(src, dst, options, pbar=False)[source]

a simple wrapper for gdal.Warp

Parameters:
spatialist.auxil.haversine(lat1, lon1, lat2, lon2)[source]

compute the distance in meters between two points in latlon

Parameters:
  • lat1 (int or float) – the latitude of point 1
  • lon1 (int or float) – the longitude of point 1
  • lat2 (int or float) – the latitude of point 2
  • lon2 (int or float) – the longitude of point 2
Returns:

the distance between point 1 and point 2 in meters

Return type:

float

spatialist.auxil.ogr2ogr(src, dst, options)[source]

a simple wrapper for gdal.VectorTranslate aka ogr2ogr

Parameters:
spatialist.auxil.utm_autodetect(spatial, crsOut)[source]

get the UTM CRS for a spatial object

The bounding box of the object is extracted, reprojected to EPSG:4326 and its center coordinate used for computing the best UTM zone fit.

Parameters:
  • spatial (Raster or Vector) – a spatial object in an arbitrary CRS
  • crsOut (str) – the output CRS type; see function crsConvert() for options
Returns:

the output CRS

Return type:

int or str or osr.SpatialReference

Database Tools

spatialist.sqlite_util.sqlite_setup(driver=':memory:', extensions=None, verbose=False)[source]

Setup a sqlite3 connection and load extensions to it. This function intends to simplify the process of loading extensions to sqlite3, which can be quite difficult depending on the version used. Particularly loading spatialite has caused quite some trouble. In recent distributions of Ubuntu this has become much easier due to a new apt package libsqlite3-mod-spatialite. For use in Windows, spatialist comes with its own spatialite DLL distribution. See here for more details on loading spatialite as an sqlite3 extension.

Parameters:
  • driver (str) – the database file or (by default) an in-memory database
  • extensions (list) – a list of extensions to load
  • verbose (bool) – print loading information?
Returns:

the database connection

Return type:

sqlite3.Connection

Example

>>> from spatialist.sqlite_util import sqlite_setup
>>> conn = sqlite_setup(extensions=['spatialite'])

Ancillary Functions

This script gathers central functions and classes for general applications

spatialist.ancillary.dissolve(inlist)[source]

list and tuple flattening

Parameters:inlist (list) – the list with sub-lists or tuples to be flattened
Returns:the flattened result
Return type:list

Examples

>>> dissolve([[1, 2], [3, 4]])
[1, 2, 3, 4]
>>> dissolve([(1, 2, (3, 4)), [5, (6, 7)]])
[1, 2, 3, 4, 5, 6, 7]
spatialist.ancillary.finder(target, matchlist, foldermode=0, regex=False, recursive=True)[source]

function for finding files/folders in folders and their subdirectories

Parameters:
  • target (str or list of str) – a directory, zip- or tar-archive or a list of them to be searched
  • matchlist (list) – a list of search patterns
  • foldermode (int) –
    • 0: only files
    • 1: files and folders
    • 2: only folders
  • regex (bool) – are the search patterns in matchlist regular expressions or unix shell standard (default)?
  • recursive (bool) – search target recursively into all subdirectories or only in the top level? This is currently only implemented for parameter target being a directory.
Returns:

the absolute names of files/folders matching the patterns

Return type:

list of str

class spatialist.ancillary.HiddenPrints[source]

Bases: object

Suppress console stdout prints, i.e. redirect them to a temporary string object.

Examples

>>> with HiddenPrints():
>>>     print('foobar')
>>> print('foobar')
spatialist.ancillary.multicore(function, cores, multiargs, **singleargs)[source]

wrapper for multicore process execution

Parameters:
  • function – individual function to be applied to each process item
  • cores (int) – the number of subprocesses started/CPUs used; this value is reduced in case the number of subprocesses is smaller
  • multiargs (dict) – a dictionary containing sub-function argument names as keys and lists of arguments to be distributed among the processes as values
  • singleargs – all remaining arguments which are invariant among the subprocesses
Returns:

the return of the function for all subprocesses

Return type:

None or list

Notes

  • all multiargs value lists must be of same length, i.e. all argument keys must be explicitly defined for each subprocess
  • all function arguments passed via singleargs must be provided with the full argument name and its value (i.e. argname=argval); default function args are not accepted
  • if the processes return anything else than None, this function will return a list of results
  • if all processes return None, this function will be of type void

Examples

>>> def add(x, y, z):
>>>     return x + y + z
>>> multicore(add, cores=2, multiargs={'x': [1, 2]}, y=5, z=9)
[15, 16]
>>> multicore(add, cores=2, multiargs={'x': [1, 2], 'y': [5, 6]}, z=9)
[15, 17]
spatialist.ancillary.parse_literal(x)[source]

return the smallest possible data type for a string or list of strings

Parameters:x (str or list) – a string to be parsed
Returns:the parsing result
Return type:int, float or str

Examples

>>> isinstance(parse_literal('1.5'), float)
True
>>> isinstance(parse_literal('1'), int)
True
>>> isinstance(parse_literal('foobar'), str)
True
spatialist.ancillary.run(cmd, outdir=None, logfile=None, inlist=None, void=True, errorpass=False, env=None)[source]
wrapper for subprocess execution including logfile writing and command prompt piping
this is a convenience wrapper around the subprocess module and calls its class Popen internally.
Parameters:
  • cmd (list) – the command arguments
  • outdir (str) – the directory to execute the command in
  • logfile (str) – a file to write stdout to
  • inlist (list) – a list of arguments passed to stdin, i.e. arguments passed to interactive input of the program
  • void (bool) – return stdout and stderr?
  • errorpass (bool) – if False, a subprocess.CalledProcessError is raised if the command fails
  • env (dict or None) – the environment to be passed to the subprocess
Returns:

a tuple of (stdout, stderr) if void is False otherwise None

Return type:

None or Tuple

spatialist.ancillary.which(program, mode=1)[source]
mimics UNIX’s which
can be replaced by shutil.which() starting from Python 3.3
Parameters:
  • program (str) – the program to be found
  • mode (os.F_OK or os.X_OK) – the mode of the found file, i.e. file exists or file is executable; see os.access()
Returns:

the full path and name of the command

Return type:

str or None

spatialist.ancillary.parallel_apply_along_axis(func1d, axis, arr, cores=4, *args, **kwargs)[source]

Like numpy.apply_along_axis() but using multiple threads. Adapted from here.

Parameters:
  • func1d (function) – the function to be applied
  • axis (int) – the axis along which to apply func1d
  • arr (numpy.ndarray) – the input array
  • cores (int) – the number of parallel cores
  • args (any) – Additional arguments to func1d.
  • kwargs (any) – Additional named arguments to func1d.
Returns:

Return type:

numpy.ndarray

ENVI HDR file manipulation

This module offers functionality for editing ENVI header files

class spatialist.envi.HDRobject(data=None)[source]

ENVI HDR info handler

Parameters:data (str, dict or None) – the file or dictionary to get the info from; If None (default), an object with default values for an empty raster file is returned

Examples

>>> from spatialist.envi import HDRobject
>>> with HDRobject('E:/test.hdr') as hdr:
>>>     hdr.band_names = ['one', 'two']
>>>     print(hdr)
>>>     hdr.write()
write(filename='same')[source]

write object to an ENVI header file

spatialist.envi.hdr(data, filename)[source]

write ENVI header files

Parameters:
  • data (str or dict) – the file or dictionary to get the info from
  • filename (str) – the HDR file to write

Data Exploration

Visualization tools using Jupyter notebooks

class spatialist.explorer.RasterViewer(filename, cmap='jet', band_indices=None, band_names=None, pmin=2, pmax=98, zmin=None, zmax=None, ts_convert=None, title=None, datalabel='data', spectrumlabel='time', fontsize=8, custom=None)[source]
Plotting utility for displaying a geocoded image stack file.
On moving the slider, the band at the slider position is read from the file and displayed.
By clicking on the band image display, you can display time series profiles.
The collected profiles can be saved to a csv file.
Parameters:
  • filename (str) – the name of the file to display
  • cmap (str) – the color map name for displaying the image. See matplotlib.colors.Colormap.
  • band_indices (list or None) – a list of indices for renaming the individual band indices in filename; e.g. -70:70, instead of the raw band indices, e.g. 1:140. The number of unique elements must be of same length as the number of bands in filename.
  • band_names (list or None) – alternative names to assign to the individual bands
  • pmin (int) – the minimum percentile for linear histogram stretching
  • pmax (int) – the maximum percentile for linear histogram stretching
  • zmin (int or float or None) – the minimum value of the displayed data range; overrides pmin
  • zmax (int or float or None) – the maximum value of the displayed data range; overrides pmax
  • ts_convert (function or None) – a function to read time stamps from the band names
  • title (str or None) – the plot title to be displayed; per default, if set to None: Figure 1, Figure 2, …
  • datalabel (str) – a label for the units of the displayed data. This also supports LaTeX mathematical notation. See Text rendering With LaTeX.
  • spectrumlabel (str) – a label for the x-axis of the vertical spectra
  • fontsize (int) – the label text font size
  • custom (list or None) –

    Custom functions for plotting figures in additional subplots. Each figure will be updated upon click on the major map display. Each function is required to take at least an argument axis. Furthermore, the following optional arguments are supported:

    • values (list): the time series values collected from the last click
    • timestamps (list): the time stamps as returned by ts_convert
    • band (int): the index of the currently displayed band
    • x (float): the x map coordinate in units of the image CRS
    • y (float): the y map coordinate in units of the image CRS

    Additional subplots are automatically added in a row-major order. The list may contain None elements to leave certain subplots empty for later usage. This might be useful for plots which are not to be updated each time the map display is clicked on.

csv(outname=None)[source]

write the collected samples to a CSV file

Parameters:outname (str) – the name of the file to write; if left at the default None, a graphical file selection dialog is opened
get_current_profile()[source]
Returns:the values of the most recently plotted time series
Return type:list
shp(outname=None)[source]

write the collected samples to a CSV file

Parameters:outname (str) – the name of the file to write; if left at the default None, a graphical file selection dialog is opened

Some general examples

in-memory vector object rasterization

Here we create a new raster data set with the same geo-information and extent as a reference data set and burn the geometries from a shapefile into it.
In this example, the shapefile contains an attribute Site_name and one of the geometries in the shapefile has a value of my_testsite for this attribute.
We use the expressions parameter to subset the shapefile and burn a value of 1 in the raster at all locations where the geometry selection overlaps. Multiple expressions can be defined together with multiple burn values.
Also, burn values can be appended to an already existing raster data set. In this case, the rasterization is performed in-memory to further use it for e.g. plotting. Alternatively, an outname can be defined to directly write the result to disk as a GeoTiff.
See spatialist.raster.rasterize() for further reference.
>>> from spatialist import Vector, Raster
>>> from spatialist.raster import rasterize
>>> import matplotlib.pyplot as plt
>>>
>>> shapefile = 'testsites.shp'
>>> rasterfile = 'extent.tif'
>>>
>>> with Raster(rasterfile) as ras:
>>>     with Vector(shapefile) as vec:
>>>         mask = rasterize(vec, reference=ras, burn_values=1, expressions=["Site_Name='my testsite'"])
>>>         plt.imshow(mask.matrix())
>>>         plt.show()