Skip to content

georeader.rasterize

This module provides functions to convert vector geometries (polygons, lines) into raster format, aligned to an existing GeoTensor or with custom georeferencing.

Overview

Rasterization is essential for:

  • Converting shapefile/GeoJSON boundaries to masks
  • Creating training labels for machine learning from vector annotations
  • Burning attribute values from GeoDataFrames into raster grids

Quick Start

from georeader import rasterize
from georeader.geotensor import GeoTensor
from shapely.geometry import Polygon
import geopandas as gpd

# Create a reference GeoTensor
gt = GeoTensor(np.zeros((100, 100)), transform, crs="EPSG:4326")

# Rasterize a single polygon to match a GeoTensor
polygon = Polygon([(-122.5, 37.5), (-122.0, 37.5), (-122.0, 38.0), (-122.5, 38.0)])
mask = rasterize.rasterize_geometry_like(polygon, gt)

# Rasterize a GeoDataFrame with attribute values
gdf = gpd.read_file("boundaries.geojson")
raster = rasterize.rasterize_geopandas_like(gdf, gt, column="class_id")

Key Functions

Function Description
rasterize_geometry_like Rasterize a single geometry to match a GeoTensor
rasterize_from_geometry Rasterize with custom transform and shape
rasterize_geopandas_like Rasterize GeoDataFrame column to match a GeoTensor
rasterize_from_geopandas Rasterize GeoDataFrame with custom georeferencing

Rasterize Module: Convert vector geometries to raster format.

This module provides functions to burn vector geometries (polygons, lines) into raster grids aligned with existing GeoTensor objects. Essential for creating masks, labels, and region-of-interest maps.

Rasterization Concepts

Converting vector shapes to pixel grids::

β”Œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”
β”‚                    RASTERIZATION PROCESS                                 β”‚
β”œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€
β”‚                                                                          β”‚
β”‚  Vector (Polygon)                    Raster (Grid)                      β”‚
β”‚  ─────────────────                   ─────────────                      β”‚
β”‚                                                                          β”‚
β”‚       ╔═══════════╗                  β”Œβ”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”                  β”‚
β”‚      ╔╝           β•šβ•—                 β”‚β–‘β”‚β–‘β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚β–‘β”‚                  β”‚
β”‚     ╔╝             β•šβ•—                β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                  β”‚
β”‚    ╔╝               β•šβ•—   ═══════►   β”‚β–‘β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚                  β”‚
β”‚    β•‘     Polygon     β•‘   Rasterize  β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                  β”‚
β”‚    β•šβ•—               ╔╝               β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚                  β”‚
β”‚     β•šβ•—             ╔╝                β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                  β”‚
β”‚      β•šβ•—           ╔╝                 β”‚β–‘β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚β–‘β”‚                  β”‚
β”‚       β•šβ•β•β•β•β•β•β•β•β•β•β•β•                  β””β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”˜                  β”‚
β”‚                                                                          β”‚
β”‚  β–‘ = fill value (outside polygon)                                       β”‚
β”‚  β–“ = burn value (inside polygon)                                        β”‚
β””β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”˜

all_touched Option

Controls which pixels are included::

β”Œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”
β”‚              all_touched=False vs all_touched=True                       β”‚
β”œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€
β”‚                                                                          β”‚
β”‚  all_touched=False (default)         all_touched=True                   β”‚
β”‚  ───────────────────────────         ──────────────────                 β”‚
β”‚                                                                          β”‚
β”‚  Only pixels with CENTER inside      ALL pixels that TOUCH the polygon  β”‚
β”‚                                                                          β”‚
β”‚  β”Œβ”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”                       β”Œβ”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”                      β”‚
β”‚  β”‚β–‘β”‚β–‘β”‚β–‘β”‚β–‘β”‚β–‘β”‚β–‘β”‚  β•± polygon edge       β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚                      β”‚
β”‚  β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€ β•±                     β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                      β”‚
β”‚  β”‚β–‘β”‚β–‘β”‚Β·β”‚β–“β”‚β–“β”‚β–‘β”‚β•±                      β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚                      β”‚
β”‚  β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                       β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                      β”‚
β”‚  β”‚β–‘β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚    Β· = center         β”‚β–‘β”‚β–“β”‚β–“β”‚β–“β”‚β–“β”‚β–‘β”‚                      β”‚
β”‚  β””β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”˜    β–“ = included       β””β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”˜                      β”‚
β”‚                                                                          β”‚
β”‚  Best for:                           Best for:                          β”‚
β”‚  β€’ Area calculations                 β€’ Inclusive masks                   β”‚
β”‚  β€’ Avoiding edge pixels              β€’ No gaps at boundaries             β”‚
β”‚  β€’ Conservative estimates            β€’ Complete coverage                 β”‚
β””β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”˜

GeoDataFrame Rasterization

Burn multiple geometries with attribute values::

β”Œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”
β”‚              GEODATAFRAME RASTERIZATION                                  β”‚
β”œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€
β”‚                                                                          β”‚
β”‚  GeoDataFrame:                       Output Raster:                     β”‚
β”‚  β”Œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”            β”Œβ”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”¬β”€β”                  β”‚
β”‚  β”‚ geometry β”‚ class_id  β”‚            β”‚0β”‚0β”‚0β”‚0β”‚0β”‚0β”‚0β”‚0β”‚                  β”‚
β”‚  β”œβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”Όβ”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€            β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                  β”‚
β”‚  β”‚ Poly A   β”‚    1      β”‚  ══════►   β”‚0β”‚1β”‚1β”‚1β”‚0β”‚2β”‚2β”‚0β”‚                  β”‚
β”‚  β”‚ Poly B   β”‚    2      β”‚            β”œβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”Όβ”€β”€                  β”‚
β”‚  β”‚ Poly C   β”‚    3      β”‚            β”‚0β”‚1β”‚1β”‚0β”‚0β”‚0β”‚3β”‚0β”‚                  β”‚
β”‚  β””β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”΄β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”˜            β””β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”΄β”€β”˜                  β”‚
β”‚                                                                          β”‚
β”‚  Usage:                                                                  β”‚
β”‚    rasterize_geodataframe(gdf, data_like, attribute="class_id")         β”‚
β””β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”€β”˜

Module Functions Overview

Single Geometry
  • :func:rasterize_geometry_like: Rasterize one geometry to match GeoData
  • :func:rasterize_geometry: Rasterize with explicit transform/shape
GeoDataFrame
  • :func:rasterize_geodataframe: Burn multiple geometries with attributes

Quick Start

Create a mask from a polygon::

from georeader import rasterize
from shapely.geometry import box

# Area of interest polygon
aoi = box(-122.5, 37.0, -122.0, 37.5)

# Create mask aligned with existing raster
mask = rasterize.rasterize_geometry_like(
    aoi,
    data_like=my_raster,
    crs_geometry="EPSG:4326",
    value=1,
    fill=0
)

Rasterize a GeoDataFrame with class labels::

import geopandas as gpd

# GeoDataFrame with land cover polygons
gdf = gpd.read_file("landcover.geojson")

# Burn class_id values into raster
labels = rasterize.rasterize_geodataframe(
    gdf,
    data_like=my_raster,
    attribute="class_id",
    all_touched=True
)

See Also

georeader.vectorize : Inverse operation (raster β†’ vector) georeader.read : Reading raster data rasterio.features.rasterize : Underlying implementation

References

  • Rasterio rasterize: https://rasterio.readthedocs.io/en/latest/api/rasterio.features.html
  • GDAL rasterize: https://gdal.org/programs/gdal_rasterize.html

rasterize_from_geometry(geometry, bounds=None, transform=None, resolution=None, window_out=None, value=1, dtype=np.uint8, crs_geom_bounds=None, fill=0, all_touched=False, return_only_data=False)

Rasterise the provided geometry over the bounds with the specified resolution, transform, shape and crs.

Parameters:

Name Type Description Default
geometry Union[Polygon, MultiPolygon, LineString]

geometry to rasterise (with crs crs_geom_bounds)

required
bounds Optional[Tuple[float, float, float, float]]

bounds where the polygons will be rasterised. (with crs crs_geom_bounds)

None
transform Optional[Affine]

if transform is provided it will use this instead of resolution (with crs crs_geom_bounds)

None
resolution Optional[Union[float, Tuple[float, float]]]

spatial resolution of the rasterised array. It won't be used if transform is provided (with crs crs_geom_bounds)

None
window_out Optional[Window]

Window out in crs_geom_bounds. If not provided it is computed from the bounds.

None
value Number

column to take the values for rasterisation.

1
dtype Any

dtype of the rasterise raster.

uint8
crs_geom_bounds Optional[Any]

CRS of geometry and bounds

None
fill Union[int, float]

fill option for rasterio.features.rasterize. Value for pixels not covered by the geometries.

0
all_touched bool

all_touched option for rasterio.features.rasterize. If True, all pixels touched by geometries will be burned in. If false, only pixels whose center is within the polygon or that are selected by Bresenham's line algorithm will be burned in.

False
return_only_data bool

if True returns only the np.ndarray without georref info.

False

Returns:

Type Description
Union[GeoTensor, ndarray]

GeoTensor or np.ndarray with shape (H, W) with the rasterised polygon

Source code in georeader/rasterize.py
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
def rasterize_from_geometry(geometry:Union[Polygon, MultiPolygon, LineString],
                            bounds:Optional[Tuple[float, float, float, float]]=None,
                            transform:Optional[rasterio.Affine]=None,
                            resolution:Optional[Union[float, Tuple[float, float]]]=None,
                            window_out:Optional[rasterio.windows.Window]=None,
                            value:Number=1,
                            dtype:Any=np.uint8,
                            crs_geom_bounds:Optional[Any]=None, fill:Union[int, float]=0, all_touched:bool=False,
                            return_only_data:bool=False)-> Union[GeoTensor, np.ndarray]:
    """
    Rasterise the provided geometry over the bounds with the specified resolution, transform, shape and crs.

    Args:
        geometry: geometry to rasterise (with crs `crs_geom_bounds`)
        bounds: bounds where the polygons will be rasterised. (with crs `crs_geom_bounds`)
        transform: if transform is provided it will use this instead of `resolution` (with crs `crs_geom_bounds`)
        resolution: spatial resolution of the rasterised array. It won't be used if transform is provided (with crs `crs_geom_bounds`)
        window_out: Window out in `crs_geom_bounds`. If not provided it is computed from the bounds.
        value: column to take the values for rasterisation.
        dtype: dtype of the rasterise raster.
        crs_geom_bounds: CRS of geometry and bounds
        fill: fill option for `rasterio.features.rasterize`. Value for pixels not covered by the geometries.
        all_touched: all_touched option for `rasterio.features.rasterize`. If True, all pixels touched 
            by geometries will be burned in.  If false, only pixels whose center is within the polygon or that
            are selected by Bresenham's line algorithm will be burned in.
        return_only_data: if `True` returns only the np.ndarray without georref info.

    Returns:
        `GeoTensor` or `np.ndarray` with shape `(H, W)` with the rasterised polygon
    """

    transform = window_utils.figure_out_transform(transform=transform, bounds=bounds,
                                                  resolution_dst=resolution)
    if window_out is None:
        window_out = rasterio.windows.from_bounds(*bounds,
                                                  transform=transform).round_lengths(op="ceil",
                                                                                     pixel_precision=PIXEL_PRECISION)


    chip_label = rasterio.features.rasterize(shapes=[(geometry, value)],
                                             out_shape=(window_out.height, window_out.width),
                                             transform=transform,
                                             dtype=dtype,
                                             fill=fill,
                                             all_touched=all_touched)
    if return_only_data:
        return chip_label

    return GeoTensor(chip_label, transform=transform, crs=crs_geom_bounds, fill_value_default=fill)

rasterize_from_geopandas(dataframe, column, bounds=None, transform=None, window_out=None, resolution=None, crs_out=None, fill=0, all_touched=False, return_only_data=False)

Rasterise the provided geodataframe over the bounds with the specified resolution.

Parameters:

Name Type Description Default
dataframe GeoDataFrame

GeoDataFrame with columns geometry and column. The 'geometry' column is expected to have shapely geometries.

required
bounds Optional[Tuple[float, float, float, float]]

bounds where the polygons will be rasterised with CRS crs_out.

None
transform Optional[Affine]

if transform is provided if will use this for the resolution.

None
resolution Optional[Union[float, Tuple[float, float]]]

spatial resolution of the rasterised array

None
window_out Optional[Window]

Window out in crs_geom_bounds. If not provided it is computed from the bounds.

None
column str

column to take the values for rasterisation.

required
crs_out Optional[Any]

defaults to dataframe.crs. This function will transform the geometries from dataframe.crs to this crs before rasterisation. bounds are in this crs.

None
fill Union[int, float]

fill option for rasterio.features.rasterize. Value for pixels not covered by the geometries.

0
all_touched bool

all_touched option for rasterio.features.rasterize. If True, all pixels touched by geometries will be burned in. If false, only pixels whose center is within the polygon or that are selected by Bresenham's line algorithm will be burned in.

False
return_only_data bool

if True returns only the np.ndarray.

False

Returns:

Type Description
Union[GeoTensor, ndarray]

GeoTensor or np.ndarray with shape (H, W) with the rasterised polygons of the dataframe

Source code in georeader/rasterize.py
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
def rasterize_from_geopandas(dataframe:gpd.GeoDataFrame,
                             column:str,
                             bounds:Optional[Tuple[float, float, float, float]]=None,
                             transform:Optional[rasterio.Affine]=None,
                             window_out:Optional[rasterio.windows.Window]=None,
                             resolution:Optional[Union[float, Tuple[float, float]]]=None,
                             crs_out:Optional[Any]=None, fill:Union[int, float]=0, all_touched:bool=False,
                             return_only_data:bool=False) -> Union[GeoTensor, np.ndarray]:
    """
    Rasterise the provided geodataframe over the bounds with the specified resolution.

    Args:
        dataframe: `GeoDataFrame` with columns `geometry` and `column`. 
            The 'geometry' column is expected to have shapely geometries.
        bounds: bounds where the polygons will be rasterised with CRS `crs_out`.
        transform: if transform is provided if will use this for the resolution.
        resolution: spatial resolution of the rasterised array
        window_out: Window out in `crs_geom_bounds`. If not provided it is computed from the bounds.
        column: column to take the values for rasterisation.
        crs_out: defaults to dataframe.crs. This function will transform the geometries from dataframe.crs to this crs
            before rasterisation. `bounds` are in this crs.
        fill: fill option for `rasterio.features.rasterize`. Value for pixels not covered by the geometries.
        all_touched: all_touched option for `rasterio.features.rasterize`. If True, all pixels touched 
            by geometries will be burned in.  If false, only pixels whose center is within the polygon or that
            are selected by Bresenham's line algorithm will be burned in.
        return_only_data: if `True` returns only the `np.ndarray`.

    Returns:
        `GeoTensor` or `np.ndarray` with shape `(H, W)` with the rasterised polygons of the dataframe
    """

    if crs_out is None:
        crs_out = str(dataframe.crs).lower()
    else:
        data_crs = str(dataframe.crs).lower()
        crs_out = str(crs_out).lower().replace("+init=","")
        if data_crs != crs_out:
            dataframe = dataframe.to_crs(crs=crs_out)

    transform = window_utils.figure_out_transform(transform=transform, bounds=bounds,
                                                  resolution_dst=resolution)
    if window_out is None:
        window_out = rasterio.windows.from_bounds(*bounds,
                                                  transform=transform).round_lengths(op="ceil",
                                                                                     pixel_precision=PIXEL_PRECISION)

    dtype = dataframe[column].dtype
    chip_label = rasterio.features.rasterize(shapes=zip(dataframe.geometry, dataframe[column]),
                                             out_shape=(window_out.height, window_out.width),
                                             transform=transform,
                                             dtype=dtype,
                                             fill=fill,
                                             all_touched=all_touched)
    if return_only_data:
        return chip_label

    return GeoTensor(chip_label, transform=transform, crs=crs_out, fill_value_default=fill)

rasterize_geometry_like(geometry, data_like, value=1, dtype=np.uint8, crs_geometry=None, fill=0, all_touched=False, return_only_data=False)

Rasterize a geometry to match an existing GeoData object's grid.

Creates a raster mask from a vector geometry, aligned to the same extent, resolution, and CRS as the reference data_like object. This is the recommended function when you have an existing raster and want to create a corresponding mask or label layer.

The function automatically reprojects the geometry if its CRS differs from the target raster.

Parameters:

Name Type Description Default
geometry Union[Polygon, MultiPolygon, LineString]

Shapely geometry to rasterize. Polygons burn filled areas, LineStrings burn line pixels.

required
data_like GeoData

Reference raster defining the output grid. The result will have matching shape[-2:], transform, and CRS.

required
value Number

Pixel value to burn inside the geometry. Default 1. Use different values for multi-class rasterization.

1
dtype Any

Output array data type. Default np.uint8. Use np.float32 for continuous values, np.int32 for large class IDs.

uint8
crs_geometry Optional[Any]

CRS of the input geometry. If provided and different from data_like.crs, geometry is reprojected automatically. Accepts EPSG codes, WKT strings, or pyproj CRS objects.

None
fill Union[int, float]

Background value for pixels outside geometry. Default 0.

0
all_touched bool

Pixel inclusion rule. Default False.

  • False: Only pixels whose center falls inside the geometry
  • True: All pixels that touch the geometry boundary
False
return_only_data bool

If True, return raw numpy array instead of GeoTensor. Default False.

False

Returns:

Type Description
Union[GeoTensor, ndarray]

Union[GeoTensor, np.ndarray]: Rasterized geometry as 2D array (H, W). GeoTensor includes georeferencing; np.ndarray if return_only_data=True.

Examples:

Create a mask from an area of interest polygon:

>>> from georeader import rasterize
>>> from shapely.geometry import box
>>>
>>> # Define AOI in WGS84
>>> aoi = box(-122.5, 37.0, -122.0, 37.5)
>>>
>>> # Create mask matching existing raster
>>> mask = rasterize.rasterize_geometry_like(
...     aoi,
...     data_like=my_raster,
...     crs_geometry="EPSG:4326",
...     value=1,
...     fill=0
... )
>>> mask.shape  # Matches my_raster spatial dims
(1000, 1000)

Rasterize with all_touched for inclusive mask:

>>> mask_inclusive = rasterize.rasterize_geometry_like(
...     aoi,
...     data_like=my_raster,
...     crs_geometry="EPSG:4326",
...     all_touched=True  # Include all edge pixels
... )

Rasterize a line feature:

>>> from shapely.geometry import LineString
>>> road = LineString([(-122.4, 37.2), (-122.3, 37.3), (-122.2, 37.25)])
>>> road_mask = rasterize.rasterize_geometry_like(
...     road,
...     data_like=my_raster,
...     crs_geometry="EPSG:4326",
...     value=255,
...     dtype=np.uint8
... )
See Also

rasterize_geopandas_like: Rasterize multiple geometries with attributes. rasterize_from_geometry: Rasterize with explicit transform/bounds. georeader.vectorize: Inverse operation (raster β†’ vector).

Source code in georeader/rasterize.py
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
def rasterize_geometry_like(geometry:Union[Polygon, MultiPolygon, LineString], 
                            data_like: GeoData, value:Number=1,
                            dtype:Any=np.uint8,
                            crs_geometry:Optional[Any]=None, fill:Union[int, float]=0, all_touched:bool=False,
                            return_only_data:bool=False)-> Union[GeoTensor, np.ndarray]:
    """
    Rasterize a geometry to match an existing GeoData object's grid.

    Creates a raster mask from a vector geometry, aligned to the same extent,
    resolution, and CRS as the reference `data_like` object. This is the
    recommended function when you have an existing raster and want to create
    a corresponding mask or label layer.

    The function automatically reprojects the geometry if its CRS differs from
    the target raster.

    Args:
        geometry (Union[Polygon, MultiPolygon, LineString]): Shapely geometry to
            rasterize. Polygons burn filled areas, LineStrings burn line pixels.
        data_like (GeoData): Reference raster defining the output grid. The result
            will have matching shape[-2:], transform, and CRS.
        value (Number): Pixel value to burn inside the geometry. Default 1.
            Use different values for multi-class rasterization.
        dtype (Any): Output array data type. Default np.uint8.
            Use np.float32 for continuous values, np.int32 for large class IDs.
        crs_geometry (Optional[Any]): CRS of the input geometry. If provided and
            different from data_like.crs, geometry is reprojected automatically.
            Accepts EPSG codes, WKT strings, or pyproj CRS objects.
        fill (Union[int, float]): Background value for pixels outside geometry.
            Default 0.
        all_touched (bool): Pixel inclusion rule. Default False.

            - False: Only pixels whose center falls inside the geometry
            - True: All pixels that touch the geometry boundary

        return_only_data (bool): If True, return raw numpy array instead of
            GeoTensor. Default False.

    Returns:
        Union[GeoTensor, np.ndarray]: Rasterized geometry as 2D array (H, W).
            GeoTensor includes georeferencing; np.ndarray if return_only_data=True.

    Examples:
        Create a mask from an area of interest polygon:

        >>> from georeader import rasterize
        >>> from shapely.geometry import box
        >>>
        >>> # Define AOI in WGS84
        >>> aoi = box(-122.5, 37.0, -122.0, 37.5)
        >>>
        >>> # Create mask matching existing raster
        >>> mask = rasterize.rasterize_geometry_like(
        ...     aoi,
        ...     data_like=my_raster,
        ...     crs_geometry="EPSG:4326",
        ...     value=1,
        ...     fill=0
        ... )
        >>> mask.shape  # Matches my_raster spatial dims
        (1000, 1000)

        Rasterize with all_touched for inclusive mask:

        >>> mask_inclusive = rasterize.rasterize_geometry_like(
        ...     aoi,
        ...     data_like=my_raster,
        ...     crs_geometry="EPSG:4326",
        ...     all_touched=True  # Include all edge pixels
        ... )

        Rasterize a line feature:

        >>> from shapely.geometry import LineString
        >>> road = LineString([(-122.4, 37.2), (-122.3, 37.3), (-122.2, 37.25)])
        >>> road_mask = rasterize.rasterize_geometry_like(
        ...     road,
        ...     data_like=my_raster,
        ...     crs_geometry="EPSG:4326",
        ...     value=255,
        ...     dtype=np.uint8
        ... )

    See Also:
        rasterize_geopandas_like: Rasterize multiple geometries with attributes.
        rasterize_from_geometry: Rasterize with explicit transform/bounds.
        georeader.vectorize: Inverse operation (raster β†’ vector).
    """
    shape_out = data_like.shape
    if crs_geometry and not window_utils.compare_crs(data_like.crs, crs_geometry):
        geometry = window_utils.polygon_to_crs(geometry, crs_geometry, data_like.crs)

    return rasterize_from_geometry(geometry, crs_geom_bounds=data_like.crs,
                                   transform=data_like.transform,
                                   window_out=rasterio.windows.Window(0, 0, width=shape_out[-1], height=shape_out[-2]),
                                   return_only_data=return_only_data,dtype=dtype, value=value,
                                   fill=fill, all_touched=all_touched)

rasterize_geopandas_like(dataframe, data_like, column, fill=0, all_touched=False, return_only_data=False)

Rasterize a GeoDataFrame to match an existing GeoData object's grid.

Burns attribute values from the specified column into a raster aligned with the reference data_like object. Ideal for creating labeled training data from vector annotations or converting land cover polygons to raster format.

The GeoDataFrame is automatically reprojected if its CRS differs from the target raster.

Parameters:

Name Type Description Default
dataframe GeoDataFrame

GeoDataFrame with geometry column and value column. Must have a valid CRS set.

required
data_like GeoData

Reference raster defining output grid (extent, resolution, CRS). Output will have matching shape[-2:] and transform.

required
column str

Column name containing values to burn. Values are cast to the output dtype. Example: 'class_id', 'land_cover', 'priority'.

required
fill Union[int, float]

Background value for pixels not covered by any geometry. Default 0.

0
all_touched bool

Pixel inclusion rule. Default False.

  • False: Burn only pixels whose center falls inside a geometry
  • True: Burn all pixels that touch any geometry boundary
False
return_only_data bool

Return raw numpy array instead of GeoTensor. Default False.

False

Returns:

Type Description
Union[GeoTensor, ndarray]

Union[GeoTensor, np.ndarray]: Rasterized geometries as 2D array (H, W). Pixel values come from the specified column. Overlapping geometries use the last geometry's value (order matters).

Examples:

Rasterize land cover polygons:

>>> import geopandas as gpd
>>> from georeader import rasterize
>>>
>>> # Load land cover polygons with class_id column
>>> gdf = gpd.read_file("landcover.geojson")
>>> print(gdf[['class_id', 'geometry']].head())
   class_id                                           geometry
0         1  POLYGON ((-122.5 37.0, -122.4 37.0, ...))
1         2  POLYGON ((-122.3 37.1, -122.2 37.1, ...))
>>>
>>> # Create label raster matching satellite image
>>> labels = rasterize.rasterize_geopandas_like(
...     gdf,
...     data_like=satellite_image,
...     column='class_id'
... )
>>> np.unique(labels.values)
array([0, 1, 2], dtype=uint8)  # 0=background, 1,2=classes

All-touched for inclusive boundaries:

>>> labels_inclusive = rasterize.rasterize_geopandas_like(
...     gdf,
...     data_like=satellite_image,
...     column='class_id',
...     all_touched=True  # Better for thin/small features
... )
See Also

rasterize_geometry_like: Rasterize single geometry with constant value. rasterize_from_geopandas: Rasterize with explicit transform/bounds.

Source code in georeader/rasterize.py
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
def rasterize_geopandas_like(dataframe:gpd.GeoDataFrame,data_like: GeoData, column:str,
                             fill:Union[int, float]=0, all_touched:bool=False,
                             return_only_data:bool=False)-> Union[GeoTensor, np.ndarray]:
    """
    Rasterize a GeoDataFrame to match an existing GeoData object's grid.

    Burns attribute values from the specified column into a raster aligned with
    the reference data_like object. Ideal for creating labeled training data
    from vector annotations or converting land cover polygons to raster format.

    The GeoDataFrame is automatically reprojected if its CRS differs from
    the target raster.

    Args:
        dataframe (gpd.GeoDataFrame): GeoDataFrame with geometry column and
            value column. Must have a valid CRS set.
        data_like (GeoData): Reference raster defining output grid (extent,
            resolution, CRS). Output will have matching shape[-2:] and transform.
        column (str): Column name containing values to burn. Values are cast
            to the output dtype. Example: 'class_id', 'land_cover', 'priority'.
        fill (Union[int, float]): Background value for pixels not covered by
            any geometry. Default 0.
        all_touched (bool): Pixel inclusion rule. Default False.

            - False: Burn only pixels whose center falls inside a geometry
            - True: Burn all pixels that touch any geometry boundary

        return_only_data (bool): Return raw numpy array instead of GeoTensor.
            Default False.

    Returns:
        Union[GeoTensor, np.ndarray]: Rasterized geometries as 2D array (H, W).
            Pixel values come from the specified column. Overlapping geometries
            use the last geometry's value (order matters).

    Examples:
        Rasterize land cover polygons:

        >>> import geopandas as gpd
        >>> from georeader import rasterize
        >>>
        >>> # Load land cover polygons with class_id column
        >>> gdf = gpd.read_file("landcover.geojson")
        >>> print(gdf[['class_id', 'geometry']].head())
           class_id                                           geometry
        0         1  POLYGON ((-122.5 37.0, -122.4 37.0, ...))
        1         2  POLYGON ((-122.3 37.1, -122.2 37.1, ...))
        >>>
        >>> # Create label raster matching satellite image
        >>> labels = rasterize.rasterize_geopandas_like(
        ...     gdf,
        ...     data_like=satellite_image,
        ...     column='class_id'
        ... )
        >>> np.unique(labels.values)
        array([0, 1, 2], dtype=uint8)  # 0=background, 1,2=classes

        All-touched for inclusive boundaries:

        >>> labels_inclusive = rasterize.rasterize_geopandas_like(
        ...     gdf,
        ...     data_like=satellite_image,
        ...     column='class_id',
        ...     all_touched=True  # Better for thin/small features
        ... )

    See Also:
        rasterize_geometry_like: Rasterize single geometry with constant value.
        rasterize_from_geopandas: Rasterize with explicit transform/bounds.
    """

    shape_out = data_like.shape
    return rasterize_from_geopandas(dataframe, column=column,
                                    crs_out=data_like.crs,
                                    transform=data_like.transform,
                                    window_out=rasterio.windows.Window(0, 0, width=shape_out[-1], height=shape_out[-2]),
                                    return_only_data=return_only_data,
                                    fill=fill, all_touched=all_touched)