xrspatial.rasterize.rasterize#

xrspatial.rasterize.rasterize(geometries: Any, width: int | None = None, height: int | None = None, bounds: Tuple[float, float, float, float] | None = None, column: str | None = None, columns: Sequence[str] | None = None, fill: float = nan, dtype: dtype | None = None, all_touched: bool = False, use_cuda: bool = False, name: str = 'rasterize', resolution: float | Tuple[float, float] | None = None, like: DataArray | None = None, merge: str | Callable = 'last', chunks: int | Tuple[int, int] | None = None, max_pixels: int = 1000000000) DataArray[source]#

Rasterize vector geometries into a 2D DataArray.

Converts geometries from a GeoDataFrame or a list of (geometry, value) pairs into a regularly gridded raster. No GDAL dependency.

Supported geometry types:

  • Polygon / MultiPolygon – scanline fill

  • LineString / MultiLineString – Bresenham line rasterization

  • Point / MultiPoint – direct pixel burn

  • GeometryCollection – recursively unpacked

Parameters:
  • geometries (GeoDataFrame or iterable of (geometry, value)) – Input vector data. If a GeoDataFrame, the column parameter selects which attribute to burn into the raster (defaults to the first numeric column). If an iterable, each element must be a (shapely.geometry, numeric_value) pair.

  • width (int, optional) – Number of columns in the output raster. Required unless resolution or like is given. Must be passed together with height: a partial override (only one of the two) is rejected with ValueError rather than silently filling the missing dimension from like or resolution.

  • height (int, optional) – Number of rows in the output raster. Required unless resolution or like is given. Must be passed together with width (see above).

  • bounds (tuple of (xmin, ymin, xmax, ymax), optional) – Geographic extent of the output raster. Inferred from the geometries (or like) if omitted.

  • column (str, optional) – Name of the GeoDataFrame column whose values are burned into the raster. Ignored when geometries is a list of pairs. Mutually exclusive with columns.

  • columns (list of str, optional) – Names of multiple GeoDataFrame columns to pass as a properties array to the merge function. Mutually exclusive with column. When given, the merge function receives a 1D float64 array of length len(columns) as its props argument.

  • fill (float, default np.nan) – Value for pixels not covered by any geometry. When dtype resolves to an integer type (either via the dtype argument or via like carrying an integer dtype), the default NaN fill is rejected with ValueError because NaN has no integer representation and the cast would silently land on a platform-specific sentinel with no _FillValue attr to mark it; pass an explicit integer sentinel (e.g. fill=0 or fill=-9999) or use a floating dtype.

  • dtype (numpy dtype, optional) – Data type of the output array. Defaults to np.float64, or to the dtype of like if provided.

  • all_touched (bool, default False) – If True, every pixel a polygon boundary passes through is burned in addition to the normal center-fill, using a supercover (Amanatides & Woo) grid traversal. Output matches rasterio.features.rasterize(..., all_touched=True) pixel-for-pixel up to rasterization tie-breaking on shared edges. If False, only pixels whose centers fall inside a polygon are burned.

  • use_cuda (bool, default False) – If True, use the CuPy/CUDA backend.

  • name (str, default 'rasterize') – Name for the output DataArray.

  • resolution (float or (x_res, y_res), optional) – Pixel size. Must be finite and > 0. When given with bounds, computes width and height automatically. A single float uses the same resolution for both axes. The requested cell size is honored exactly. When the bounds don’t divide evenly by the resolution, the output grid extends past the original xmax / ymin (anchored at xmin / ymax) rather than shrinking the cell to fit. Matches the rasterio.transform.from_origin convention.

  • like (xr.DataArray, optional) – Template raster. Width, height, bounds, and dtype are copied from this array. Bounds and dtype can be overridden one at a time; width and height must be overridden together (passing only one raises ValueError). Must have uniformly spaced x and y dim coords – the rasterizer only writes to a regular grid, so a non-uniform like is rejected with ValueError rather than silently producing pixel labels that don’t match the data. Both descending (top-down, ymax first) and ascending (bottom-up, ymin first) y coords are supported – the burned rows are flipped to match so result.sel(y=...) lines up with the geometry in world coordinates either way, and result.y always equals like.y exactly.

  • merge (str or callable, default 'last') –

    How to combine values when geometries overlap.

    Built-in modes (pass as string):

    • 'last' – last geometry in input order wins

    • 'first' – first geometry wins

    • 'max' / 'min' – keep the larger / smaller value

    • 'sum' – add values together

    • 'count' – count overlapping geometries

    Custom merge function (pass a callable):

    For CPU backends, pass a @ngjit-decorated function. For GPU backends (use_cuda=True), pass a @numba.cuda.jit(device=True) function. Signature:

    merge_fn(pixel, props, is_first) -> float64
    
    • pixel: current pixel value (fill value on first write)

    • props: 1D float64 array of property values for the geometry

    • is_first: 1 on first write to this pixel, 0 otherwise

  • chunks (int or (int, int), optional) – If given, use the dask backend and split the output raster into tiles of this size (row_chunk, col_chunk). Both axes must be > 0. A single int uses the same chunk size for both axes. Combined with use_cuda to select dask+numpy vs dask+cupy.

  • max_pixels (int, default 1_000_000_000) – Safety cap on the resolved output size (width * height). The function raises ValueError before any host or device allocation if the cap is exceeded. Raise this explicitly when rasterizing a legitimately large grid.

Returns:

2D raster with dims ('y', 'x').

Return type:

xr.DataArray

Examples

>>> from shapely.geometry import box
>>> result = rasterize([(box(0, 0, 5, 5), 1.0)],
...                    width=10, height=10)

>>> # Using resolution instead of width/height:
>>> result = rasterize(gdf, resolution=0.5,
...                    bounds=(0, 0, 10, 10), column='value')

>>> # Match an existing raster grid:
>>> zones = rasterize(gdf, like=elevation, column='zone')

>>> # Sum overlapping values:
>>> density = rasterize(gdf, width=100, height=100,
...                     column='pop', merge='sum', fill=0)