Skip to content

Dataset Class#

  • Detailed class diagram for the Dataset class and related components:
Hold "Ctrl" to enable pan & zoom
classDiagram
    %% configuration class
    class Config {
    }

    %% abstract base class for rasters
    class AbstractDataset {
        +__init__(src, access)
        +__str__()
        +__repr__()
        +access()
        +raster()
        +raster(value)
        +values()
        +rows()
        +columns()
        +shape()
        +geotransform()
        +top_left_corner()
        +epsg()
        +epsg(value)
        +crs()
        +crs(value)
        +cell_size()
        +no_data_value()
        +no_data_value(value)
        +meta_data()
        +meta_data(value)
        +block_size()
        +block_size(value)
        +file_name()
        +driver_type()
        +read_file(path, read_only)
        +read_array(band, window)
        +_read_block(band, window)
        +plot(band, exclude_value, rgb, surface_reflectance, cutoff, overview, overview_index, **kwargs)
    }

    %% concrete raster class
    class Dataset {
        +__init__(src, access)
        +__str__()
        +__repr__()
        +access()
        +raster()
        +raster(value)
        +values()
        +rows()
        +columns()
        +shape()
        +geotransform()
        +epsg()
        +epsg(value)
        +crs()
        +crs(value)
        +cell_size()
        +band_count()
        +band_names()
        +band_names(name_list)
        +band_units()
        +band_units(value)
        +no_data_value()
        +no_data_value(value)
        +meta_data()
        +meta_data(value)
        +block_size()
        +block_size(value)
        +file_name()
        +driver_type()
        +scale()
        +scale(value)
        +offset()
        +offset(value)
        +read_file(path, read_only)
        +create_from_array(arr, top_left_corner, cell_size, epsg)
        +read_array(band, window)
        +_read_block(band, window)
        +plot(band, exclude_value, rgb, surface_reflectance, cutoff, overview, overview_index, **kwargs)
        +to_file(path, driver, band)
        +to_crs(to_epsg, method, maintain_alignment)
        +resample(cell_size, method)
        +align(alignment_src)
        +crop(mask, touch)
        +merge(src, dst, no_data_value, init, n)
        +apply(ufunc)
        +overlay(classes_map, exclude_value)
    }



    %% Driver catalog
    class _utils_Catalog {
    }

    %% NetCDF
    class NetCDF {
    }

    %% error classes
    class _errors_ReadOnlyError
    class _errors_DatasetNotFoundError
    class _errors_NoDataValueError
    class _errors_AlignmentError
    class _errors_DriverNotExistError
    class _errors_FileFormatNotSupportedError
    class _errors_OptionalPackageDoesNotExist
    class _errors_FailedToSaveError
    class _errors_OutOfBoundsError

    %% inheritance relations
    AbstractDataset <|-- Dataset
    Dataset <|-- NetCDF

    %% composition/usage relations
    AbstractDataset ..> _utils_Catalog : "uses Catalog constant"
    AbstractDataset ..> feature_FeatureCollection : "vector ops"
    Dataset ..> feature_FeatureCollection : "vector ops"
    Dataset ..> _errors_ReadOnlyError : "raises"
    Dataset ..> _errors_AlignmentError : "raises"
    Dataset ..> _errors_NoDataValueError : "raises"
    Dataset ..> _errors_FailedToSaveError : "raises"
    Dataset ..> _errors_OutOfBoundsError : "raises"
    NetCDF ..> _errors_OptionalPackageDoesNotExist : "raises"
    Config ..> Dataset : "initialises raster settings"
Hold "Ctrl" to enable pan & zoom
classDiagram

    %% Central dataset class with its main attributes
    class Dataset {
        +raster
        +cell_size
        +values
        +shape
        +rows
        +columns
        +pivot_point
        +geotransform
        +bounds
        +bbox
        +epsg
        +crs
        +lon
        +lat
        +x
        +y
        +band_count
        +band_names
        +variables
        +no_data_value
        +meta_data
        +dtype
        +gdal_dtype
        +numpy_dtype
        +file_name
        +time_stamp
        +driver_type
    }

    %% Group: visualisation functionality
    class Visualization {
        +plot()
        +overview_count()
        +read_overview_array()
        +create_overviews()
        +recreate_overviews()
        +get_overview()
    }
    Dataset --> Visualization : «visualisation»

    %% Group: data access methods
    class AccessData {
        +read_array()
        +get_variables()
        +count_domain_cells()
        +get_band_names()
        +extract()
        +stats()
    }
    Dataset --> AccessData : «data access»

    %% Group: mathematical operations on raster values
    class MathOperations {
        +apply()
        +fill()
        +normalize()
        +cluster()
        +cluster2()
        +get_tile()
        +groupNeighbours()
    }
    Dataset --> MathOperations : «math ops»

    %% Group: spatial operations and reprojection
    class SpatialOperations {
        +to_crs()
        +resample()
        +align()
        +crop()
        +locate_points()
        +overlay()
        +extract()
        +footprint()
    }
    Dataset --> SpatialOperations : «spatial ops»

    %% Group: conversion to other data types
    class Conversion {
        +to_feature_collection()
    }
    Dataset --> Conversion : «conversion»

    %% Group: coordinate system handling
    class OSR {
        +create_sr_from_epsg()
    }
    Dataset --> OSR : «osr»

    %% Group: bounding‐box and bounds calculations
    class BBoxBounds {
        +calculate_bbox()
        +calculate_bounds()
    }
    Dataset --> BBoxBounds : «bbox/bounds»

    %% Group: CRS/EPSG getters
    class CrsEpsg {
        +get_crs()
        +get_epsg()
    }
    Dataset --> CrsEpsg : «crs/epsg»

    %% Group: latitude/longitude getters
    class LatLon {
        +get_lat_lon()
    }
    Dataset --> LatLon : «lat/lon»

    %% Group: band names management
    class BandNames {
        +get_band_names_internal()
        +set_band_names()
    }
    Dataset --> BandNames : «band names»

    %% Group: timestamp handling
    class TimeStamp {
        +get_time_variable()
        +read_variable()
    }
    Dataset --> TimeStamp : «time»

    %% Group: handling of no‐data values
    class NoDataValue {
        +set_no_data_value()
        +set_no_data_value_backend()
        +change_no_data_value_attr()
    }
    Dataset --> NoDataValue : «no data value»

    %% Group: helpers for creating GDAL datasets
    class GdalDataset {
        +create_empty_driver()
        +create_driver_from_scratch()
        +create_mem_gtiff_dataset()
    }
    Dataset --> GdalDataset : «gdal creation»

    %% Group: factory methods for creating Dataset objects
    class CreateObject {
        +from_gdal_dataset()
        +read_file()
        +create_from_array()
        +dataset_like()
    }
    Dataset --> CreateObject : «object factory»

pyramids.dataset.Dataset #

Bases: BandMetadata, IO, Spatial, Analysis, Vectorize, Cell, AbstractDataset

Single-band or multi-band raster dataset (GeoTIFF, etc.).

Wraps a GDAL dataset with spatial operations (crop, reproject, align, mosaic), band-level I/O, and no-data handling. For NetCDF files use the :class:~pyramids.netcdf.NetCDF subclass; for temporal stacks of rasters use :class:~pyramids.dataset.DatasetCollection.

Source code in src/pyramids/dataset/dataset.py
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
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
251
252
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
302
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
381
382
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
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
class Dataset(  # type: ignore[misc]
    BandMetadata,
    IO,
    Spatial,
    Analysis,
    Vectorize,
    Cell,
    AbstractDataset,
):
    """Single-band or multi-band raster dataset (GeoTIFF, etc.).

    Wraps a GDAL dataset with spatial operations (crop, reproject, align,
    mosaic), band-level I/O, and no-data handling.  For NetCDF files use
    the :class:`~pyramids.netcdf.NetCDF` subclass; for temporal stacks of
    rasters use :class:`~pyramids.dataset.DatasetCollection`.
    """

    def __init__(self, src: gdal.Dataset, access: str = "read_only"):
        """__init__."""
        self.logger = logging.getLogger(__name__)
        super().__init__(src, access=access)

        self._no_data_value = [
            src.GetRasterBand(i).GetNoDataValue() for i in range(1, self.band_count + 1)
        ]
        self._band_names = self._get_band_names()
        self._band_units = [
            src.GetRasterBand(i).GetUnitType() for i in range(1, self.band_count + 1)
        ]

    def _update_inplace(self, src: gdal.Dataset, access: str | None = None) -> None:
        """Swap internal state from a new GDAL dataset.

        Creates a fresh Dataset and copies its internal data
        into this instance, similar to pandas' _update_inplace.
        """
        new = Dataset(src, access=access or self._access)
        self.__dict__.update(new.__dict__)

    def __str__(self) -> str:
        """__str__."""
        message = f"""
            Top Left Corner: {self.top_left_corner}
            Cell size: {self.cell_size}
            Dimension: {self.rows} * {self.columns}
            EPSG: {self.epsg}
            Number of Bands: {self.band_count}
            Band names: {self.band_names}
            Band colors: {self.band_color}
            Band units: {self.band_units}
            Scale: {self.scale}
            Offset: {self.offset}
            Mask: {self._no_data_value[0]}
            Data type: {self.dtype[0]}
            File: {self.file_name}
        """
        return message

    def __repr__(self) -> str:
        """__repr__."""
        return str(gdal.Info(self.raster))

    @property
    def access(self) -> str:
        """
        Access mode.

        Returns:
            str:
                The access mode of the dataset (read_only/write).
        """
        return str(super().access)

    @property
    def raster(self) -> gdal.Dataset:
        """Base GDAL Dataset (read-only)."""
        return super().raster

    @property
    def rows(self) -> int:
        """Number of rows in the raster array."""
        return int(self._rows)

    @property
    def columns(self) -> int:
        """Number of columns in the raster array."""
        return int(self._columns)

    @property
    def shape(self) -> tuple[int, int, int]:
        """Shape (bands, rows, columns)."""
        return self.band_count, self.rows, self.columns

    @property
    def geotransform(self) -> tuple[float, float, float, float, float, float]:
        """WKT projection.

        (top left corner X/lon coordinate, cell_size, 0, top left corner y/lat coordinate, 0, -cell_size).

        See Also:
            - Dataset.top_left_corner: Coordinate of the top left corner of the dataset.
            - Dataset.epsg: EPSG number of the dataset coordinate reference system.
        """
        gt: tuple[float, float, float, float, float, float] = self._geotransform
        return gt

    @property
    def epsg(self) -> int:
        """EPSG number."""
        crs = self.raster.GetProjection()
        return FeatureCollection.get_epsg_from_prj(crs)

    @epsg.setter
    def epsg(self, value: int):
        """EPSG number."""
        sr = Dataset._create_sr_from_epsg(value)
        self.raster.SetProjection(sr.ExportToWkt())
        self._update_inplace(self._raster)

    @property
    def crs(self) -> str:
        """Coordinate reference system.

        Returns:
            str:
                the coordinate reference system of the dataset.

        See Also:
            Dataset.set_crs : Set the Coordinate Reference System (CRS).
            Dataset.to_crs : Reproject the dataset to any projection.
            Dataset.epsg : epsg number of the dataset coordinate reference system.
        """
        return self._get_crs()

    @crs.setter
    def crs(self, value: str):
        """Coordinate reference system.

        Args:
            value (str):
                WellKnownText (WKT) string.

        See Also:
            - Dataset.set_crs: Set the Coordinate Reference System (CRS).
            - Dataset.to_crs: Reproject the dataset to any projection.
            - Dataset.epsg: EPSG number of the dataset coordinate reference system.
        """
        self.set_crs(value)

    @property
    def cell_size(self) -> float:
        """Cell size."""
        return float(self._cell_size)

    @property
    def band_count(self) -> int:
        """Number of bands in the raster."""
        return int(self._band_count)

    @property
    def band_names(self) -> list[str]:
        """Band names."""
        return self._get_band_names()

    @band_names.setter
    def band_names(self, name_list: list):
        """Band names."""
        self._set_band_names(name_list)

    @property
    def band_units(self) -> list[str]:
        """Band units."""
        return self._band_units

    @band_units.setter
    def band_units(self, value: list[str]):
        """Band units setter."""
        self._band_units = value
        for i, val in enumerate(value):
            self._iloc(i).SetUnitType(val)

    @property
    def no_data_value(self) -> list:
        """No data value that marks the cells out of the domain."""
        return list(self._no_data_value)

    @no_data_value.setter
    def no_data_value(self, value: list | Number):
        """no_data_value.

        No data value that marks the cells out of the domain

        Notes:
            - The setter does not change the values of the cells to the new no_data_value, it only changes the
            `no_data_value` attribute.
            - Use this method to change the `no_data_value` attribute to match the value that is stored in the cells.
            - To change the values of the cells, to the new no_data_value, use the `change_no_data_value` method.

        See Also:
            - Dataset.change_no_data_value: Change the No Data Value.
        """
        if isinstance(value, list):
            for i, val in enumerate(value):
                self._change_no_data_value_attr(i, val)
        else:
            self._change_no_data_value_attr(0, value)

    @property
    def meta_data(self):
        """Meta-data."""
        return super().meta_data

    @meta_data.setter
    def meta_data(self, value: dict[str, str]):
        """Meta-data."""
        for key, val in value.items():
            self._raster.SetMetadataItem(key, val)

    @property
    def block_size(self) -> list[tuple[int, int]]:
        """Block Size.

        The block size is the size of the block that the raster is divided into, the block size is used to
        read and write the raster data in blocks.

        See Also:
            - Dataset.get_block_arrangement: Get block arrangement to read the dataset in chunks.
            - Dataset.get_tile: Get tiles.
            - Dataset.read_array: Read the data stored in the dataset bands.
        """
        return self._block_size

    @block_size.setter
    def block_size(self, value: list[tuple[int, int]]):
        """Block Size.

        Args:
            value (List[Tuple[int, int]]):
                block size for each band in the raster(512, 512).
        """
        if len(value[0]) != 2:
            raise ValueError("block size should be a tuple of 2 integers")

        self._block_size = value

    @property
    def file_name(self):
        """File name."""
        return super().file_name

    @property
    def driver_type(self):
        """Driver Type."""
        return super().driver_type

    @property
    def scale(self) -> list[float]:
        """Scale.

        The value of the scale is used to convert the pixel values to the real-world values.
        """
        scale_list = []
        for i in range(self.band_count):
            band_scale = self._iloc(i).GetScale()
            scale_list.append(band_scale if band_scale is not None else 1.0)
        return scale_list

    @scale.setter
    def scale(self, value: list[float]):
        """Scale."""
        for i, val in enumerate(value):
            self._iloc(i).SetScale(val)

    @property
    def offset(self):
        """Offset.

        The value of the offset is used to convert the pixel values to the real-world values.
        """
        offset_list = []
        for i in range(self.band_count):
            band_offset = self._iloc(i).GetOffset()
            offset_list.append(band_offset if band_offset is not None else 0)
        return offset_list

    @offset.setter
    def offset(self, value: list[float]):
        """Offset."""
        for i, val in enumerate(value):
            self._iloc(i).SetOffset(val)

    @property
    def top_left_corner(self):
        """Top left corner coordinates.

        See Also:
            - Dataset.geotransform: Dataset geotransform.
        """
        return super().top_left_corner

    @property
    def bounds(self) -> GeoDataFrame:
        """Bounds - the bbox as a geodataframe with a polygon geometry.

        See Also:
            - Dataset.bbox: Dataset bounding box.
        """
        return self._calculate_bounds()

    @property
    def bbox(self) -> list:
        """Bound box [xmin, ymin, xmax, ymax].

        See Also:
            - Dataset.bounds: Dataset bounding polygon.
        """
        return self._calculate_bbox()

    @property
    def lon(self) -> np.ndarray:
        """Longitude coordinates.

        See Also:
            - Dataset.x: Dataset x coordinates.
            - Dataset.lat: Dataset latitude.
        """
        x_coords = self.get_x_lon_dimension_array(
            self.top_left_corner[0], self.cell_size, self.columns
        )
        return x_coords

    @property
    def lat(self) -> np.ndarray:
        """Latitude-coordinate.

        See Also:
            - Dataset.x: Dataset x coordinates.
            - Dataset.y: Dataset y coordinates.
            - Dataset.lon: Dataset longitude.
        """
        y_coords = self.get_y_lat_dimension_array(
            self.top_left_corner[1], self.cell_size, self.rows
        )
        return y_coords

    @property
    def x(self) -> np.ndarray:
        """X-coordinate/Longitude.

        See Also:
            - Dataset.lat: Dataset latitude.
            - Dataset.y: Dataset y coordinates.
            - Dataset.lon: Dataset longitude.
        """
        # X_coordinate = upper-left corner x + index * cell size + cell-size/2
        return self.lon

    @property
    def y(self) -> np.ndarray:
        """Y-coordinate/Latitude.

        See Also:
            - Dataset.x: Dataset y coordinates.
            - Dataset.lat: Dataset latitude.
            - Dataset.lon: Dataset longitude.
        """
        # Y_coordinate = upper-left corner y - index * cell size - cell-size/2
        return self.lat

    @property
    def gdal_dtype(self):
        """Data Type."""
        return [
            self.raster.GetRasterBand(i).DataType for i in range(1, self.band_count + 1)
        ]

    @property
    def numpy_dtype(self) -> list[type]:
        """List of the numpy data Type of each band, the data type is a numpy function."""
        return [
            DTYPE_CONVERSION_DF.loc[DTYPE_CONVERSION_DF["gdal"] == i, "numpy"].values[0]
            for i in self.gdal_dtype
        ]

    @property
    def dtype(self) -> list[str]:
        """List of the data Type of each band as strings."""
        return [
            DTYPE_CONVERSION_DF.loc[DTYPE_CONVERSION_DF["gdal"] == i, "name"].values[0]
            for i in self.gdal_dtype
        ]

    @classmethod
    def read_file(
        cls,
        path: str | Path,
        read_only=True,
        file_i: int = 0,
    ) -> Dataset:
        """read_file.

        Args:
            path (str):
                Path of file to open.
            read_only (bool):
                File mode, set to False, to open in "update" mode.
            file_i (int):
                Index to the file inside the compressed file you want to read, if the compressed file
                has only one file. Default is 0.

        Returns:
            Dataset:
                Opened dataset instance.

        See Also:
            - Dataset.read_array: Read the values stored in a dataset band.
        """
        src = _io.read_file(path, read_only=read_only, file_i=file_i)
        return cls(src, access="read_only" if read_only else "write")

    def copy(self, path: str | Path | None = None) -> Dataset:
        """Deep copy.

        Args:
            path (str, optional):
                Destination path to save the copied dataset. If None is passed, the copied dataset
                will be created in memory.
        """
        if path is None:
            path = ""
            driver = "MEM"
        else:
            driver = "GTiff"

        src = gdal.GetDriverByName(driver).CreateCopy(str(path), self._raster)

        return Dataset(src, access="write")

    def close(self) -> None:
        """Close the dataset.

        Safe to call multiple times — subsequent calls after the first are no-ops.
        """
        if self._raster is not None:
            self._raster.FlushCache()
            self._raster = None

    @staticmethod
    def _create_dataset(
        cols: int,
        rows: int,
        bands: int,
        dtype: int,
        driver: str = "MEM",
        path: str | Path | None = None,
    ) -> gdal.Dataset:
        """Create a GDAL driver.

            creates a driver and save it to disk and in memory if the path is not given.

        Args:
            cols (int):
                Number of columns.
            rows (int):
                Number of rows.
            bands (int):
                Number of bands.
            dtype:
                GDAL data type.
            driver (str):
                Driver type ["GTiff", "MEM"].
            path (str):
                Path to save the GTiff driver.

        Returns:
            gdal driver
        """
        if path:
            driver = "GTiff" if driver == "MEM" else driver
            if not isinstance(path, (str, Path)):
                raise TypeError(
                    f"The path input should be string or Path type, given: {type(path)}"
                )
            path = Path(path)
            if driver == "GTiff" and path.suffix != ".tif":
                raise TypeError(
                    "The path to save the created raster should end with .tif"
                )
            # LZW is a lossless compression method achieve the highest compression but with a lot
            # of computations.
            src = gdal.GetDriverByName(driver).Create(
                str(path), cols, rows, bands, dtype, ["COMPRESS=LZW"]
            )
        else:
            # for memory drivers
            driver = "MEM"
            src = gdal.GetDriverByName(driver).Create("", cols, rows, bands, dtype)
        return src

    @classmethod
    def _build_dataset(
        cls,
        cols: int,
        rows: int,
        bands: int,
        dtype: int,
        geo: tuple,
        crs: str,
        no_data_value,
        driver: str = "MEM",
        path: str | Path | None = None,
        access: str = "write",
    ) -> Dataset:
        """Create a GDAL dataset, set its spatial metadata, and wrap it as a Dataset.

        Consolidates the repeated pattern of _create_dataset + SetGeoTransform +
        SetProjection + wrap + _set_no_data_value into a single helper.

        Args:
            cols (int): Number of columns.
            rows (int): Number of rows.
            bands (int): Number of bands.
            dtype (int): GDAL data type.
            geo (tuple): Geotransform tuple.
            crs (str): Projection as WKT string.
            no_data_value: No-data value. Scalar (broadcast to all bands) or list (one per band).
            driver (str): Driver type. Default is "MEM".
            path (str | Path | None): Path for disk-based drivers.
            access (str): Access mode for the Dataset wrapper. Default is "write".
                Note: MEM driver datasets can be written to regardless of access mode since
                the access flag is enforced at the pyramids level, not by GDAL.

        Returns:
            Dataset: A fully configured Dataset object.
        """
        dst = cls._create_dataset(cols, rows, bands, dtype, driver=driver, path=path)
        dst.SetGeoTransform(geo)
        dst.SetProjection(crs)
        dst_obj = cls(dst, access=access)
        dst_obj._set_no_data_value(no_data_value=no_data_value)
        return dst_obj

    @classmethod
    def create(
        cls,
        cell_size: int | float,
        rows: int,
        columns: int,
        dtype: str,
        bands: int,
        top_left_corner: tuple,
        epsg: int,
        no_data_value: Any | None = None,
        path: str | Path | None = None,
    ) -> Dataset:
        """Create a new dataset and fill it with the no_data_value.

        The new dataset will have an array filled with the no_data_value.

        Args:
            cell_size (int|float):
                Cell size.
            rows (int):
                Number of rows.
            columns (int):
                Number of columns.
            dtype (str):
                Data type.
            bands (int|None):
                Number of bands to create in the output raster.
            top_left_corner (Tuple):
                Coordinates of the top left corner point.
            epsg (int):
                EPSG number to identify the projection of the coordinates in the created raster.
            no_data_value (float|None):
                No data value.
            path (str, optional):
                Path on disk; if None, the dataset is created in memory. Default is None.

        Returns:
            Dataset: A new dataset
        """
        # Create the driver.
        gdal_dtype = numpy_to_gdal_dtype(dtype)
        dst = Dataset._create_dataset(columns, rows, bands, gdal_dtype, path=path)
        sr = Dataset._create_sr_from_epsg(epsg)
        geotransform = (
            top_left_corner[0],
            cell_size,
            0,
            top_left_corner[1],
            0,
            -1 * cell_size,
        )
        dst.SetGeoTransform(geotransform)
        # Set the projection.
        dst.SetProjection(sr.ExportToWkt())

        dst = cls(dst, access="write")
        if no_data_value is not None:
            dst._set_no_data_value(no_data_value=no_data_value)

        return dst

    @classmethod
    def create_from_array(  # type: ignore[override]
        cls,
        arr: np.ndarray,
        top_left_corner: tuple[float, float] | None = None,
        cell_size: int | float | None = None,
        geo: tuple[float, float, float, float, float, float] | None = None,
        epsg: str | int = 4326,
        no_data_value: Any | list = DEFAULT_NO_DATA_VALUE,
        driver_type: str = "MEM",
        path: str | Path | None = None,
    ) -> Dataset:
        """Create a new dataset from an array.

        Args:
            arr (np.ndarray):
                Numpy array.
            top_left_corner (Tuple[float, float], optional):
                The coordinates of the top left corner of the dataset.
            cell_size (int|float, optional):
                Cell size in the same units of the coordinate reference system defined by the `epsg`
                parameter.
            geo (Tuple[float, float, float, float, float, float], optional):
                Geotransform tuple (minimum lon/x, pixel-size, rotation, maximum lat/y, rotation,
                pixel-size).
            epsg (int):
                Integer reference number to the projection (https://epsg.io/).
            no_data_value (Any, optional):
                No data value to mask the cells out of the domain. The default is -9999.
            driver_type (str, optional):
                Driver type ["GTiff", "MEM", "netcdf"]. Default is "MEM".
            path (str, optional):
                Path to save the driver.

        Returns:
            Dataset:
                Dataset object will be returned.
        """
        if geo is None:
            if top_left_corner is None or cell_size is None:
                raise ValueError(
                    "Either top_left_corner and cell_size or geo should be provided."
                )
            geo = (
                top_left_corner[0],
                cell_size,
                0,
                top_left_corner[1],
                0,
                -1 * cell_size,
            )

        if arr.ndim == 2:
            bands = 1
            rows = int(arr.shape[0])
            cols = int(arr.shape[1])
        else:
            bands = arr.shape[0]
            rows = int(arr.shape[1])
            cols = int(arr.shape[2])

        dst_obj = cls._create_gtiff_from_array(
            arr,
            cols,
            rows,
            bands,
            geo,
            epsg,
            no_data_value,
            driver_type=driver_type,
            path=path,
        )

        return dst_obj

    @staticmethod
    def _create_gtiff_from_array(
        arr: np.ndarray,
        cols: int,
        rows: int,
        bands: int = 1,
        geo: tuple[float, float, float, float, float, float] | None = None,
        epsg: str | int | None = None,
        no_data_value: Any | list = DEFAULT_NO_DATA_VALUE,
        driver_type: str = "MEM",
        path: str | Path | None = None,
    ) -> Dataset:
        dtype = numpy_to_gdal_dtype(arr)
        dst_ds = Dataset._create_dataset(
            cols, rows, bands, dtype, driver=driver_type, path=path
        )

        if epsg is None:
            raise ValueError("epsg must be provided")

        srse = Dataset._create_sr_from_epsg(epsg=int(epsg))
        dst_ds.SetProjection(srse.ExportToWkt())
        dst_ds.SetGeoTransform(geo)

        dst_obj = Dataset(dst_ds, access="write")
        dst_obj._set_no_data_value(no_data_value=no_data_value)
        dst_obj._raster.FlushCache()

        if bands == 1:
            dst_obj.raster.GetRasterBand(1).WriteArray(arr)
        else:
            for i in range(bands):
                dst_obj.raster.GetRasterBand(i + 1).WriteArray(arr[i, :, :])

        dst_obj._raster.FlushCache()
        return dst_obj

    @classmethod
    def dataset_like(
        cls,
        src: Dataset,
        array: np.ndarray,
        path: str | Path | None = None,
    ) -> Dataset:
        """Create a new dataset like another dataset.

        dataset_like method creates a Dataset from an array like another source dataset. The new dataset
        will have the same `projection`, `coordinates` or the `top left corner` of the original dataset,
        `cell size`, `no_data_velue`, and number of `rows` and `columns`.
        the array and the source dataset should have the same number of columns and rows

        Args:
            src (Dataset):
                source raster to get the spatial information
            array (ndarray):
                data to store in the new dataset.
            path (str, optional):
                path to save the new dataset, if not given, the method will return in-memory dataset.

        Returns:
            Dataset:
                if the `path` is given, the method will save the new raster to the given path, else the
                method will return an in-memory dataset.
        """
        if not isinstance(array, np.ndarray):
            raise TypeError("array should be of type numpy array")

        if array.ndim == 2:
            bands = 1
        else:
            bands = array.shape[0]

        dtype = numpy_to_gdal_dtype(array)

        dst_obj = cls._build_dataset(
            src.columns, src.rows, bands, dtype, src.geotransform, src.crs,
            src.no_data_value[0], path=path,
        )

        if bands == 1:
            dst_obj.raster.GetRasterBand(1).WriteArray(array)
        else:
            for band_i in range(bands):
                dst_obj.raster.GetRasterBand(band_i + 1).WriteArray(array[band_i, :, :])

        if path is not None:
            dst_obj.raster.FlushCache()

        return dst_obj

    @staticmethod
    def _create_sr_from_epsg(epsg: int) -> SpatialReference:
        """Create a spatial reference object from EPSG number.

        https://gdal.org/tutorials/osr_api_tut.html

        Args:
            epsg (int):
                EPSG number.

        Returns:
            SpatialReference:
                SpatialReference object.
        """
        sr = osr.SpatialReference()
        sr.ImportFromEPSG(int(epsg))
        return sr

access property #

Access mode.

Returns:

Name Type Description
str str

The access mode of the dataset (read_only/write).

raster property #

Base GDAL Dataset (read-only).

rows property #

Number of rows in the raster array.

columns property #

Number of columns in the raster array.

shape property #

Shape (bands, rows, columns).

geotransform property #

WKT projection.

(top left corner X/lon coordinate, cell_size, 0, top left corner y/lat coordinate, 0, -cell_size).

See Also
  • Dataset.top_left_corner: Coordinate of the top left corner of the dataset.
  • Dataset.epsg: EPSG number of the dataset coordinate reference system.

crs property writable #

Coordinate reference system.

Returns:

Name Type Description
str str

the coordinate reference system of the dataset.

See Also

Dataset.set_crs : Set the Coordinate Reference System (CRS). Dataset.to_crs : Reproject the dataset to any projection. Dataset.epsg : epsg number of the dataset coordinate reference system.

band_count property #

Number of bands in the raster.

no_data_value property writable #

No data value that marks the cells out of the domain.

block_size property writable #

Block Size.

The block size is the size of the block that the raster is divided into, the block size is used to read and write the raster data in blocks.

See Also
  • Dataset.get_block_arrangement: Get block arrangement to read the dataset in chunks.
  • Dataset.get_tile: Get tiles.
  • Dataset.read_array: Read the data stored in the dataset bands.

scale property writable #

Scale.

The value of the scale is used to convert the pixel values to the real-world values.

offset property writable #

Offset.

The value of the offset is used to convert the pixel values to the real-world values.

top_left_corner property #

Top left corner coordinates.

See Also
  • Dataset.geotransform: Dataset geotransform.

bounds property #

Bounds - the bbox as a geodataframe with a polygon geometry.

See Also
  • Dataset.bbox: Dataset bounding box.

bbox property #

Bound box [xmin, ymin, xmax, ymax].

See Also
  • Dataset.bounds: Dataset bounding polygon.

lon property #

Longitude coordinates.

See Also
  • Dataset.x: Dataset x coordinates.
  • Dataset.lat: Dataset latitude.

lat property #

Latitude-coordinate.

See Also
  • Dataset.x: Dataset x coordinates.
  • Dataset.y: Dataset y coordinates.
  • Dataset.lon: Dataset longitude.

x property #

X-coordinate/Longitude.

See Also
  • Dataset.lat: Dataset latitude.
  • Dataset.y: Dataset y coordinates.
  • Dataset.lon: Dataset longitude.

y property #

Y-coordinate/Latitude.

See Also
  • Dataset.x: Dataset y coordinates.
  • Dataset.lat: Dataset latitude.
  • Dataset.lon: Dataset longitude.

numpy_dtype property #

List of the numpy data Type of each band, the data type is a numpy function.

dtype property #

List of the data Type of each band as strings.

__init__(src, access='read_only') #

init.

Source code in src/pyramids/dataset/dataset.py
def __init__(self, src: gdal.Dataset, access: str = "read_only"):
    """__init__."""
    self.logger = logging.getLogger(__name__)
    super().__init__(src, access=access)

    self._no_data_value = [
        src.GetRasterBand(i).GetNoDataValue() for i in range(1, self.band_count + 1)
    ]
    self._band_names = self._get_band_names()
    self._band_units = [
        src.GetRasterBand(i).GetUnitType() for i in range(1, self.band_count + 1)
    ]

__str__() #

str.

Source code in src/pyramids/dataset/dataset.py
def __str__(self) -> str:
    """__str__."""
    message = f"""
        Top Left Corner: {self.top_left_corner}
        Cell size: {self.cell_size}
        Dimension: {self.rows} * {self.columns}
        EPSG: {self.epsg}
        Number of Bands: {self.band_count}
        Band names: {self.band_names}
        Band colors: {self.band_color}
        Band units: {self.band_units}
        Scale: {self.scale}
        Offset: {self.offset}
        Mask: {self._no_data_value[0]}
        Data type: {self.dtype[0]}
        File: {self.file_name}
    """
    return message

__repr__() #

repr.

Source code in src/pyramids/dataset/dataset.py
def __repr__(self) -> str:
    """__repr__."""
    return str(gdal.Info(self.raster))

read_file(path, read_only=True, file_i=0) classmethod #

read_file.

Parameters:

Name Type Description Default
path str

Path of file to open.

required
read_only bool

File mode, set to False, to open in "update" mode.

True
file_i int

Index to the file inside the compressed file you want to read, if the compressed file has only one file. Default is 0.

0

Returns:

Name Type Description
Dataset Dataset

Opened dataset instance.

See Also
  • Dataset.read_array: Read the values stored in a dataset band.
Source code in src/pyramids/dataset/dataset.py
@classmethod
def read_file(
    cls,
    path: str | Path,
    read_only=True,
    file_i: int = 0,
) -> Dataset:
    """read_file.

    Args:
        path (str):
            Path of file to open.
        read_only (bool):
            File mode, set to False, to open in "update" mode.
        file_i (int):
            Index to the file inside the compressed file you want to read, if the compressed file
            has only one file. Default is 0.

    Returns:
        Dataset:
            Opened dataset instance.

    See Also:
        - Dataset.read_array: Read the values stored in a dataset band.
    """
    src = _io.read_file(path, read_only=read_only, file_i=file_i)
    return cls(src, access="read_only" if read_only else "write")

copy(path=None) #

Deep copy.

Parameters:

Name Type Description Default
path str

Destination path to save the copied dataset. If None is passed, the copied dataset will be created in memory.

None
Source code in src/pyramids/dataset/dataset.py
def copy(self, path: str | Path | None = None) -> Dataset:
    """Deep copy.

    Args:
        path (str, optional):
            Destination path to save the copied dataset. If None is passed, the copied dataset
            will be created in memory.
    """
    if path is None:
        path = ""
        driver = "MEM"
    else:
        driver = "GTiff"

    src = gdal.GetDriverByName(driver).CreateCopy(str(path), self._raster)

    return Dataset(src, access="write")

close() #

Close the dataset.

Safe to call multiple times — subsequent calls after the first are no-ops.

Source code in src/pyramids/dataset/dataset.py
def close(self) -> None:
    """Close the dataset.

    Safe to call multiple times — subsequent calls after the first are no-ops.
    """
    if self._raster is not None:
        self._raster.FlushCache()
        self._raster = None

create(cell_size, rows, columns, dtype, bands, top_left_corner, epsg, no_data_value=None, path=None) classmethod #

Create a new dataset and fill it with the no_data_value.

The new dataset will have an array filled with the no_data_value.

Parameters:

Name Type Description Default
cell_size int | float

Cell size.

required
rows int

Number of rows.

required
columns int

Number of columns.

required
dtype str

Data type.

required
bands int | None

Number of bands to create in the output raster.

required
top_left_corner Tuple

Coordinates of the top left corner point.

required
epsg int

EPSG number to identify the projection of the coordinates in the created raster.

required
no_data_value float | None

No data value.

None
path str

Path on disk; if None, the dataset is created in memory. Default is None.

None

Returns:

Name Type Description
Dataset Dataset

A new dataset

Source code in src/pyramids/dataset/dataset.py
@classmethod
def create(
    cls,
    cell_size: int | float,
    rows: int,
    columns: int,
    dtype: str,
    bands: int,
    top_left_corner: tuple,
    epsg: int,
    no_data_value: Any | None = None,
    path: str | Path | None = None,
) -> Dataset:
    """Create a new dataset and fill it with the no_data_value.

    The new dataset will have an array filled with the no_data_value.

    Args:
        cell_size (int|float):
            Cell size.
        rows (int):
            Number of rows.
        columns (int):
            Number of columns.
        dtype (str):
            Data type.
        bands (int|None):
            Number of bands to create in the output raster.
        top_left_corner (Tuple):
            Coordinates of the top left corner point.
        epsg (int):
            EPSG number to identify the projection of the coordinates in the created raster.
        no_data_value (float|None):
            No data value.
        path (str, optional):
            Path on disk; if None, the dataset is created in memory. Default is None.

    Returns:
        Dataset: A new dataset
    """
    # Create the driver.
    gdal_dtype = numpy_to_gdal_dtype(dtype)
    dst = Dataset._create_dataset(columns, rows, bands, gdal_dtype, path=path)
    sr = Dataset._create_sr_from_epsg(epsg)
    geotransform = (
        top_left_corner[0],
        cell_size,
        0,
        top_left_corner[1],
        0,
        -1 * cell_size,
    )
    dst.SetGeoTransform(geotransform)
    # Set the projection.
    dst.SetProjection(sr.ExportToWkt())

    dst = cls(dst, access="write")
    if no_data_value is not None:
        dst._set_no_data_value(no_data_value=no_data_value)

    return dst

create_from_array(arr, top_left_corner=None, cell_size=None, geo=None, epsg=4326, no_data_value=DEFAULT_NO_DATA_VALUE, driver_type='MEM', path=None) classmethod #

Create a new dataset from an array.

Parameters:

Name Type Description Default
arr ndarray

Numpy array.

required
top_left_corner Tuple[float, float]

The coordinates of the top left corner of the dataset.

None
cell_size int | float

Cell size in the same units of the coordinate reference system defined by the epsg parameter.

None
geo Tuple[float, float, float, float, float, float]

Geotransform tuple (minimum lon/x, pixel-size, rotation, maximum lat/y, rotation, pixel-size).

None
epsg int

Integer reference number to the projection (https://epsg.io/).

4326
no_data_value Any

No data value to mask the cells out of the domain. The default is -9999.

DEFAULT_NO_DATA_VALUE
driver_type str

Driver type ["GTiff", "MEM", "netcdf"]. Default is "MEM".

'MEM'
path str

Path to save the driver.

None

Returns:

Name Type Description
Dataset Dataset

Dataset object will be returned.

Source code in src/pyramids/dataset/dataset.py
@classmethod
def create_from_array(  # type: ignore[override]
    cls,
    arr: np.ndarray,
    top_left_corner: tuple[float, float] | None = None,
    cell_size: int | float | None = None,
    geo: tuple[float, float, float, float, float, float] | None = None,
    epsg: str | int = 4326,
    no_data_value: Any | list = DEFAULT_NO_DATA_VALUE,
    driver_type: str = "MEM",
    path: str | Path | None = None,
) -> Dataset:
    """Create a new dataset from an array.

    Args:
        arr (np.ndarray):
            Numpy array.
        top_left_corner (Tuple[float, float], optional):
            The coordinates of the top left corner of the dataset.
        cell_size (int|float, optional):
            Cell size in the same units of the coordinate reference system defined by the `epsg`
            parameter.
        geo (Tuple[float, float, float, float, float, float], optional):
            Geotransform tuple (minimum lon/x, pixel-size, rotation, maximum lat/y, rotation,
            pixel-size).
        epsg (int):
            Integer reference number to the projection (https://epsg.io/).
        no_data_value (Any, optional):
            No data value to mask the cells out of the domain. The default is -9999.
        driver_type (str, optional):
            Driver type ["GTiff", "MEM", "netcdf"]. Default is "MEM".
        path (str, optional):
            Path to save the driver.

    Returns:
        Dataset:
            Dataset object will be returned.
    """
    if geo is None:
        if top_left_corner is None or cell_size is None:
            raise ValueError(
                "Either top_left_corner and cell_size or geo should be provided."
            )
        geo = (
            top_left_corner[0],
            cell_size,
            0,
            top_left_corner[1],
            0,
            -1 * cell_size,
        )

    if arr.ndim == 2:
        bands = 1
        rows = int(arr.shape[0])
        cols = int(arr.shape[1])
    else:
        bands = arr.shape[0]
        rows = int(arr.shape[1])
        cols = int(arr.shape[2])

    dst_obj = cls._create_gtiff_from_array(
        arr,
        cols,
        rows,
        bands,
        geo,
        epsg,
        no_data_value,
        driver_type=driver_type,
        path=path,
    )

    return dst_obj

dataset_like(src, array, path=None) classmethod #

Create a new dataset like another dataset.

dataset_like method creates a Dataset from an array like another source dataset. The new dataset will have the same projection, coordinates or the top left corner of the original dataset, cell size, no_data_velue, and number of rows and columns. the array and the source dataset should have the same number of columns and rows

Parameters:

Name Type Description Default
src Dataset

source raster to get the spatial information

required
array ndarray

data to store in the new dataset.

required
path str

path to save the new dataset, if not given, the method will return in-memory dataset.

None

Returns:

Name Type Description
Dataset Dataset

if the path is given, the method will save the new raster to the given path, else the method will return an in-memory dataset.

Source code in src/pyramids/dataset/dataset.py
@classmethod
def dataset_like(
    cls,
    src: Dataset,
    array: np.ndarray,
    path: str | Path | None = None,
) -> Dataset:
    """Create a new dataset like another dataset.

    dataset_like method creates a Dataset from an array like another source dataset. The new dataset
    will have the same `projection`, `coordinates` or the `top left corner` of the original dataset,
    `cell size`, `no_data_velue`, and number of `rows` and `columns`.
    the array and the source dataset should have the same number of columns and rows

    Args:
        src (Dataset):
            source raster to get the spatial information
        array (ndarray):
            data to store in the new dataset.
        path (str, optional):
            path to save the new dataset, if not given, the method will return in-memory dataset.

    Returns:
        Dataset:
            if the `path` is given, the method will save the new raster to the given path, else the
            method will return an in-memory dataset.
    """
    if not isinstance(array, np.ndarray):
        raise TypeError("array should be of type numpy array")

    if array.ndim == 2:
        bands = 1
    else:
        bands = array.shape[0]

    dtype = numpy_to_gdal_dtype(array)

    dst_obj = cls._build_dataset(
        src.columns, src.rows, bands, dtype, src.geotransform, src.crs,
        src.no_data_value[0], path=path,
    )

    if bands == 1:
        dst_obj.raster.GetRasterBand(1).WriteArray(array)
    else:
        for band_i in range(bands):
            dst_obj.raster.GetRasterBand(band_i + 1).WriteArray(array[band_i, :, :])

    if path is not None:
        dst_obj.raster.FlushCache()

    return dst_obj