Skip to content

Python API

xcube_resampling.spatial.resample_in_space(source_ds, target_gm=None, source_gm=None, variables=None, interp_methods=None, agg_methods=None, recover_nans=False, fill_values=None, tile_size=None)

Resample the spatial dimensions of a dataset to match a target grid mapping.

Depending on the regularity and compatibility of the source and target grid mappings, this function will either rectify, reproject, or affine-transform the spatial dimensions of source_ds.

Parameters:

Name Type Description Default
source_ds Dataset

The input xarray.Dataset. Data variables must have dimensions in the following order: optional third dimension followed by the y-dimension (e.g., "y" or "lat") and the x-dimension (e.g., "x" or "lon").

required
target_gm GridMapping | None

The target GridMapping to which the dataset should be resampled. Must be regular. If not provided, a default regular grid is derived from source_gm using to_regular(tile_size).

None
source_gm GridMapping | None

The GridMapping describing the source dataset's spatial layout. If not provided, it is inferred from source_ds using GridMapping.from_dataset(source_ds).

None
variables str | Iterable[str] | None

A single variable name or iterable of variable names to be resampled. If None, all data variables will be processed.

None
interp_methods InterpMethods | None

Optional interpolation method to be used for upsampling spatial data variables. Can be a single interpolation method for all variables or a dictionary mapping variable names or dtypes to interpolation method. Supported methods include:

  • 0 (nearest neighbor)
  • 1 (linear / bilinear)
  • "nearest"
  • "triangular"
  • "bilinear"

The default is 0 for integer arrays, else 1.

None
agg_methods AggMethods | None

Optional aggregation methods for downsampling spatial variables. Can be a single method for all variables or a dictionary mapping variable names or dtypes to methods. Supported methods include: "center", "count", "first", "last", "max", "mean", "median", "mode", "min", "prod", "std", "sum", and "var". Defaults to "center" for integer arrays, else "mean".

None
recover_nans RecoverNans

Optional boolean or mapping to enable NaN recovery during upsampling (only applies when interpolation method is not nearest). Can be a single boolean or a dictionary mapping variable names or dtypes to booleans. Defaults to False.

False
fill_values FillValues | None

Optional fill value(s) for areas outside input coverage. Can be a single value or dictionary by variable or type. If not provided, defaults based on data type are used:

  • float: NaN
  • uint8: 255
  • uint16: 65535
  • other ints: -1
None
tile_size int | tuple[int, int] | None

Optional tile size used when generating a regular grid from an irregular source grid mapping. Only used if target_gm is not provided.

None

Returns:

Type Description
Dataset

A new dataset that has been spatially resampled to match the target grid mapping.

Notes
  • If source_gm is not provided, it is inferred from source_ds.
  • If target_gm is not provided, and the source is irregular, it is derived from source_gm.to_regular(tile_size=tile_size).
  • If both grid mappings are regular and approximately equal, the original dataset is returned unchanged.
  • If the transformation can be represented as an affine mapping, it is applied directly for performance.
  • If the source is irregular, rectification is applied.
  • Otherwise, a reprojection is performed.
  • See the xcube-resampling documentation for more details.
Source code in xcube_resampling/spatial.py
 40
 41
 42
 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
def resample_in_space(
    source_ds: xr.Dataset,
    target_gm: GridMapping | None = None,
    source_gm: GridMapping | None = None,
    variables: str | Iterable[str] | None = None,
    interp_methods: InterpMethods | None = None,
    agg_methods: AggMethods | None = None,
    recover_nans: RecoverNans = False,
    fill_values: FillValues | None = None,
    tile_size: int | tuple[int, int] | None = None,
) -> xr.Dataset:
    """
    Resample the spatial dimensions of a dataset to match a target grid mapping.

    Depending on the regularity and compatibility of the source and target grid
    mappings, this function will either rectify, reproject, or affine-transform the
    spatial dimensions of `source_ds`.

    Args:
        source_ds: The input xarray.Dataset. Data variables must have dimensions
            in the following order: optional third dimension followed by the
            y-dimension (e.g., "y" or "lat") and the x-dimension (e.g., "x" or "lon").
        target_gm: The target GridMapping to which the dataset should be resampled.
            Must be regular. If not provided, a default regular grid is derived
            from `source_gm` using `to_regular(tile_size)`.
        source_gm: The GridMapping describing the source dataset's spatial layout.
            If not provided, it is inferred from `source_ds` using
            `GridMapping.from_dataset(source_ds)`.
        variables: A single variable name or iterable of variable names to be
            resampled. If None, all data variables will be processed.
        interp_methods: Optional interpolation method to be used for upsampling spatial
            data variables. Can be a single interpolation method for all variables or a
            dictionary mapping variable names or dtypes to interpolation method.
            Supported methods include:

            - `0` (nearest neighbor)
            - `1` (linear / bilinear)
            - `"nearest"`
            - `"triangular"`
            - `"bilinear"`

            The default is `0` for integer arrays, else `1`.
        agg_methods: Optional aggregation methods for downsampling spatial variables.
            Can be a single method for all variables or a dictionary mapping variable
            names or dtypes to methods. Supported methods include:
                "center", "count", "first", "last", "max", "mean", "median",
                "mode", "min", "prod", "std", "sum", and "var".
            Defaults to "center" for integer arrays, else "mean".
        recover_nans: Optional boolean or mapping to enable NaN recovery during
            upsampling (only applies when interpolation method is not nearest).
            Can be a single boolean or a dictionary mapping variable names or dtypes
            to booleans. Defaults to False.
        fill_values: Optional fill value(s) for areas outside input coverage.
            Can be a single value or dictionary by variable or type. If not provided,
            defaults based on data type are used:

            - float: NaN
            - uint8: 255
            - uint16: 65535
            - other ints: -1

        tile_size: Optional tile size used when generating a regular grid from
            an irregular source grid mapping. Only used if `target_gm` is not provided.

    Returns:
        A new dataset that has been spatially resampled to match the target grid
            mapping.

    Notes:
        - If `source_gm` is not provided, it is inferred from `source_ds`.
        - If `target_gm` is not provided, and the source is irregular, it is
          derived from `source_gm.to_regular(tile_size=tile_size)`.
        - If both grid mappings are regular and approximately equal, the original
          dataset is returned unchanged.
        - If the transformation can be represented as an affine mapping, it is
          applied directly for performance.
        - If the source is irregular, rectification is applied.
        - Otherwise, a reprojection is performed.
        - See the [xcube-resampling documentation](https://xcube-dev.github.io/xcube-resampling/)
          for more details.
    """
    if source_gm is None:
        source_gm = GridMapping.from_dataset(source_ds)

    if not source_gm.is_regular:
        return rectify_dataset(
            source_ds,
            target_gm=target_gm,
            source_gm=source_gm,
            variables=variables,
            interp_methods=interp_methods,
            agg_methods=agg_methods,
            recover_nans=recover_nans,
            fill_values=fill_values,
            tile_size=tile_size,
        )
    else:
        if target_gm is None:
            LOG.warning(
                "If source grid mapping is regular `target_gm` must be given. "
                "Source dataset is returned."
            )
            return source_ds
        GridMapping.assert_regular(target_gm, name="target_gm")
        if source_gm.is_close(target_gm):
            return source_ds

        if _can_apply_affine_transform(source_gm, target_gm):
            return affine_transform_dataset(
                source_ds,
                target_gm,
                source_gm=source_gm,
                variables=variables,
                interp_methods=interp_methods,
                agg_methods=agg_methods,
                recover_nans=recover_nans,
                fill_values=fill_values,
            )
        else:
            return reproject_dataset(
                source_ds,
                target_gm,
                source_gm=source_gm,
                variables=variables,
                interp_methods=interp_methods,
                agg_methods=agg_methods,
                recover_nans=recover_nans,
                fill_values=fill_values,
            )

xcube_resampling.affine.affine_transform_dataset(source_ds, target_gm, source_gm=None, variables=None, interp_methods=None, agg_methods=None, recover_nans=False, fill_values=None)

Apply an affine transformation to the spatial dimensions of a dataset, transforming it from the source grid mapping to the target grid mapping.

Parameters:

Name Type Description Default
source_ds Dataset

The input dataset to be transformed.

required
target_gm GridMapping

The target grid mapping defining the spatial reference and output geometry.

required
source_gm GridMapping | None

The grid mapping of the input dataset. If None, it is inferred from the dataset.

None
variables str | Iterable[str] | None

Optional variable(s) to transform. If None, all variables are used.

None
interp_methods InterpMethods | None

Optional interpolation method to be used for upsampling spatial data variables. Can be a single interpolation method for all variables or a dictionary mapping variable names or dtypes to interpolation method. Supported methods include:

  • 0 (nearest neighbor)
  • 1 (linear / bilinear)
  • "nearest"
  • "bilinear"

The default is 0 for integer arrays, else 1.

None
agg_methods AggMethods | None

Optional aggregation methods for downsampling spatial variables. Can be a single method for all variables or a dictionary mapping variable names or dtypes to methods. Supported methods include: "center", "count", "first", "last", "max", "mean", "median", "mode", "min", "prod", "std", "sum", and "var". Defaults to "center" for integer arrays, else "mean".

None
recover_nans RecoverNans

Optional boolean or mapping to enable NaN recovery during upsampling (only applies when interpolation method is not nearest). Can be a single boolean or a dictionary mapping variable names or dtypes to booleans. Defaults to False.

False
fill_values FillValues | None

Optional fill value(s) for areas outside the input bounds. Can be a single value or a dictionary mapping variable names or dtypes to fill values. If not provided, defaults are:

  • float: NaN
  • uint8: 255
  • uint16: 65535
  • other integers: -1
None

Returns:

Type Description
Dataset

A new dataset resampled and aligned to the target grid mapping. Data variables without spatial dimensions are copied to the output. Data variables with only one spatial dimension are ignored.

Source code in xcube_resampling/affine.py
 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
def affine_transform_dataset(
    source_ds: xr.Dataset,
    target_gm: GridMapping,
    source_gm: GridMapping | None = None,
    variables: str | Iterable[str] | None = None,
    interp_methods: InterpMethods | None = None,
    agg_methods: AggMethods | None = None,
    recover_nans: RecoverNans = False,
    fill_values: FillValues | None = None,
) -> xr.Dataset:
    """
    Apply an affine transformation to the spatial dimensions of a dataset,
    transforming it from the source grid mapping to the target grid mapping.

    Args:
        source_ds: The input dataset to be transformed.
        target_gm: The target grid mapping defining the spatial reference and
            output geometry.
        source_gm: The grid mapping of the input dataset. If None, it is inferred
            from the dataset.
        variables: Optional variable(s) to transform. If None, all variables are used.
        interp_methods: Optional interpolation method to be used for upsampling spatial
            data variables. Can be a single interpolation method for all variables or a
            dictionary mapping variable names or dtypes to interpolation method.
            Supported methods include:

            - `0` (nearest neighbor)
            - `1` (linear / bilinear)
            - `"nearest"`
            - `"bilinear"`

            The default is `0` for integer arrays, else `1`.
        agg_methods: Optional aggregation methods for downsampling spatial variables.
            Can be a single method for all variables or a dictionary mapping variable
            names or dtypes to methods. Supported methods include:
                "center", "count", "first", "last", "max", "mean", "median",
                "mode", "min", "prod", "std", "sum", and "var".
            Defaults to "center" for integer arrays, else "mean".
        recover_nans: Optional boolean or mapping to enable NaN recovery during
            upsampling (only applies when interpolation method is not nearest).
            Can be a single boolean or a dictionary mapping variable names or dtypes
            to booleans. Defaults to False.
        fill_values: Optional fill value(s) for areas outside the input bounds.
            Can be a single value or a dictionary mapping variable names or dtypes
            to fill values. If not provided, defaults are:

            - float: NaN
            - uint8: 255
            - uint16: 65535
            - other integers: -1

    Returns:
        A new dataset resampled and aligned to the target grid mapping.
            Data variables without spatial dimensions are copied to the output.
            Data variables with only one spatial dimension are ignored.
    """
    if source_gm is None:
        source_gm = GridMapping.from_dataset(source_ds)
    source_ds = normalize_grid_mapping(source_ds, source_gm)

    assert _can_apply_affine_transform(source_gm, target_gm), (
        f"Affine transformation cannot be applied to source CRS "
        f"{source_gm.crs.name!r} and target CRS {target_gm.crs.name!r}"
    )

    source_ds = _select_variables(source_ds, variables)

    target_ds = resample_dataset(
        source_ds,
        target_gm.ij_transform_to(source_gm),
        (source_gm.xy_dim_names[1], source_gm.xy_dim_names[0]),
        target_gm.size,
        target_gm.tile_size,
        interp_methods,
        agg_methods,
        recover_nans,
        fill_values,
    )

    # assign coordinates from target grid-mapping
    x_name, y_name = target_gm.xy_var_names
    target_ds = target_ds.assign_coords(
        {x_name: target_gm.x_coords, y_name: target_gm.y_coords}
    )

    return target_ds

xcube_resampling.reproject.reproject_dataset(source_ds, target_gm, source_gm=None, variables=None, interp_methods=None, agg_methods=None, recover_nans=False, fill_values=None)

Reproject a dataset from one coordinate reference system (CRS) to another.

This function transforms a dataset’s 2D spatial variables to match a new CRS and grid layout defined by target_gm. It handles interpolation, optional downsampling, and fill values for areas outside the input bounds.

Parameters:

Name Type Description Default
source_ds Dataset

The input dataset to be reprojected.

required
target_gm GridMapping

The target grid mapping that defines the spatial reference and grid layout in the target CRS.

required
source_gm GridMapping | None

Optional source grid mapping of the input dataset. If not provided, it is inferred from the dataset.

None
variables str | Iterable[str] | None

Optional variable name or list of variable names to reproject. If None, all suitable variables are processed.

None
interp_methods InterpMethods | None

Optional interpolation method to be used for upsampling spatial data variables. Can be a single interpolation method for all variables or a dictionary mapping variable names or dtypes to interpolation method. Supported methods include:

  • 0 (nearest neighbor)
  • 1 (linear / bilinear)
  • "nearest"
  • "triangular"
  • "bilinear"

The default is 0 for integer arrays, else 1.

None
agg_methods AggMethods | None

Optional aggregation methods for downsampling spatial variables. Can be a single method for all variables or a dictionary mapping variable names or dtypes to methods. Supported methods include: "center", "count", "first", "last", "max", "mean", "median", "mode", "min", "prod", "std", "sum", and "var". Defaults to "center" for integer arrays, else "mean".

None
recover_nans RecoverNans

Optional boolean or mapping to enable NaN recovery during upsampling (only applies when interpolation method is not nearest). Can be a single boolean or a dictionary mapping variable names or dtypes to booleans. Defaults to False.

False
fill_values FillValues | None

Optional fill value(s) to assign outside source coverage. Can be a single value or dictionary by variable or type. If not set, defaults are:

  • float: NaN
  • uint8: 255
  • uint16: 65535
  • other integers: -1
None

Returns:

Type Description
Dataset

A new dataset with variables reprojected to the target CRS and grid. Variables without 2D spatial dimensions are copied as-is. 1D spatial coordinate variables are ignored in the output.

Source code in xcube_resampling/reproject.py
 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
def reproject_dataset(
    source_ds: xr.Dataset,
    target_gm: GridMapping,
    source_gm: GridMapping | None = None,
    variables: str | Iterable[str] | None = None,
    interp_methods: InterpMethods | None = None,
    agg_methods: AggMethods | None = None,
    recover_nans: RecoverNans = False,
    fill_values: FillValues | None = None,
) -> xr.Dataset:
    """
    Reproject a dataset from one coordinate reference system (CRS) to another.

    This function transforms a dataset’s 2D spatial variables to match a new CRS and
    grid layout defined by `target_gm`. It handles interpolation, optional
    downsampling, and fill values for areas outside the input bounds.

    Args:
        source_ds: The input dataset to be reprojected.
        target_gm: The target grid mapping that defines the spatial reference and
            grid layout in the target CRS.
        source_gm: Optional source grid mapping of the input dataset. If not
            provided, it is inferred from the dataset.
        variables: Optional variable name or list of variable names to reproject.
            If None, all suitable variables are processed.
        interp_methods: Optional interpolation method to be used for upsampling spatial
            data variables. Can be a single interpolation method for all variables or a
            dictionary mapping variable names or dtypes to interpolation method.
            Supported methods include:

            - `0` (nearest neighbor)
            - `1` (linear / bilinear)
            - `"nearest"`
            - `"triangular"`
            - `"bilinear"`

            The default is `0` for integer arrays, else `1`.
        agg_methods: Optional aggregation methods for downsampling spatial variables.
            Can be a single method for all variables or a dictionary mapping variable
            names or dtypes to methods. Supported methods include:
                "center", "count", "first", "last", "max", "mean", "median",
                "mode", "min", "prod", "std", "sum", and "var".
            Defaults to "center" for integer arrays, else "mean".
        recover_nans: Optional boolean or mapping to enable NaN recovery during
            upsampling (only applies when interpolation method is not nearest).
            Can be a single boolean or a dictionary mapping variable names or dtypes
            to booleans. Defaults to False.
        fill_values: Optional fill value(s) to assign outside source coverage.
            Can be a single value or dictionary by variable or type. If not set,
            defaults are:

            - float: NaN
            - uint8: 255
            - uint16: 65535
            - other integers: -1

    Returns:
        A new dataset with variables reprojected to the target CRS and
            grid. Variables without 2D spatial dimensions are copied as-is.
            1D spatial coordinate variables are ignored in the output.
    """

    if source_gm is None:
        source_gm = GridMapping.from_dataset(source_ds)
    if source_gm.is_j_axis_up:
        v_var = source_gm.xy_var_names[1]
        source_ds = source_ds.isel({v_var: slice(None, None, -1)})
        source_gm = GridMapping.from_dataset(source_ds)

    source_ds = normalize_grid_mapping(source_ds, source_gm)

    source_ds = _select_variables(source_ds, variables)

    transformer = pyproj.Transformer.from_crs(
        target_gm.crs, source_gm.crs, always_xy=True
    )

    # If source has higher resolution than target, downscale first, then reproject
    source_ds, source_gm = _downscale_source_dataset(
        source_ds,
        source_gm,
        target_gm,
        transformer,
        interp_methods,
        agg_methods,
        recover_nans,
    )

    # For each bounding box in the target grid mapping:
    # - determine the indices of the bbox in the source dataset
    # - extract the corresponding coordinates for each bbox in the source dataset
    # - compute the pad_width to handle areas requested by target_gm that exceed the
    #   bounds of source_gm.
    scr_ij_bboxes, x_coords, y_coords, pad_width = _get_scr_bboxes_indices(
        transformer, source_gm, target_gm
    )

    # transform grid points from target grid mapping to source grid mapping
    source_xx, source_yy = _transform_gridpoints(transformer, target_gm)

    # reproject dataset
    x_name, y_name = source_gm.xy_var_names
    coords = source_ds.coords.to_dataset()
    coords = coords.drop_vars((x_name, y_name))
    x_name, y_name = target_gm.xy_var_names
    coords[x_name] = target_gm.x_coords
    coords[y_name] = target_gm.y_coords
    coords["spatial_ref"] = xr.DataArray(0, attrs=target_gm.crs.to_cf())
    target_ds = xr.Dataset(coords=coords, attrs=source_ds.attrs)

    yx_dims = (source_gm.xy_dim_names[1], source_gm.xy_dim_names[0])
    for var_name, data_array in source_ds.items():
        if data_array.dims[-2:] == yx_dims:
            assert len(data_array.dims) in (
                2,
                3,
            ), f"Data variable {var_name} has {len(data_array.dims)} dimensions."

            target_ds[var_name] = _reproject_data_array(
                data_array,
                var_name,
                source_gm,
                target_gm,
                source_xx,
                source_yy,
                x_coords,
                y_coords,
                scr_ij_bboxes,
                pad_width,
                interp_methods,
                fill_values,
            )
        elif yx_dims[0] not in data_array.dims and yx_dims[1] not in data_array.dims:
            target_ds[var_name] = data_array

    return target_ds

xcube_resampling.rectify.rectify_dataset(source_ds, target_gm=None, source_gm=None, variables=None, interp_methods=None, agg_methods=None, recover_nans=False, fill_values=None, tile_size=None)

Rectify a dataset with non-regular grid to a regular grid defined by a target grid mapping.

This function transforms spatial coordinates to a regular grid while preserving data values. It optionally downsamples high-resolution inputs prior to rectifying.

Parameters:

Name Type Description Default
source_ds Dataset

The source dataset with 2D spatial coordinate variables.

required
target_gm GridMapping | None

Optional target grid mapping defining the output regular grid. If not provided, one is derived from the source grid mapping.

None
source_gm GridMapping | None

Optional grid mapping of the source dataset. If not given, it is inferred from the dataset.

None
variables str | Iterable[str] | None

Optional variable(s) to rectify. If None, all eligible variables are processed.

None
interp_methods InterpMethods | None

Optional interpolation method to be used for upsampling spatial data variables. Can be a single interpolation method for all variables or a dictionary mapping variable names or dtypes to interpolation method. Supported methods include:

  • 0 (nearest neighbor)
  • 1 (linear / bilinear)
  • "nearest"
  • "triangular"
  • "bilinear"

The default is 0 for integer arrays, else 1.

None
agg_methods AggMethods | None

Optional aggregation methods for downsampling spatial variables. Can be a single method for all variables or a dictionary mapping variable names or dtypes to methods. Supported methods include: "center", "count", "first", "last", "max", "mean", "median", "mode", "min", "prod", "std", "sum", and "var". Defaults to "center" for integer arrays, else "mean".

None
recover_nans RecoverNans

Optional boolean or mapping to enable NaN recovery during upsampling (only applies when interpolation method is not nearest). Can be a single boolean or a dictionary mapping variable names or dtypes to booleans. Defaults to False.

False
fill_values FillValues | None

Optional fill value(s) for areas outside input coverage. Can be a single value or dictionary by variable or type. If not provided, defaults based on data type are used:

  • float: NaN
  • uint8: 255
  • uint16: 65535
  • other ints: -1
None
tile_size int | tuple[int, int] | None

Optional tile size for inferring a regular grid, if target_gm is not provided.

None

Returns:

Type Description
Dataset

A new dataset with spatial variables rectified to a regular grid. Variables not having 2D spatial dimensions are copied as-is. 1D spatial coordinate variables are ignored in the output.

Source code in xcube_resampling/rectify.py
 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
def rectify_dataset(
    source_ds: xr.Dataset,
    target_gm: GridMapping | None = None,
    source_gm: GridMapping | None = None,
    variables: str | Iterable[str] | None = None,
    interp_methods: InterpMethods | None = None,
    agg_methods: AggMethods | None = None,
    recover_nans: RecoverNans = False,
    fill_values: FillValues | None = None,
    tile_size: int | tuple[int, int] | None = None,
) -> xr.Dataset:
    """
    Rectify a dataset with non-regular grid to a regular grid defined by a target
    grid mapping.

    This function transforms spatial coordinates to a regular grid while preserving
    data values. It optionally downsamples high-resolution inputs prior to rectifying.

    Args:
        source_ds: The source dataset with 2D spatial coordinate variables.
        target_gm: Optional target grid mapping defining the output regular grid.
            If not provided, one is derived from the source grid mapping.
        source_gm: Optional grid mapping of the source dataset. If not given, it is
            inferred from the dataset.
        variables: Optional variable(s) to rectify. If None, all eligible variables
            are processed.
        interp_methods: Optional interpolation method to be used for upsampling spatial
            data variables. Can be a single interpolation method for all variables or a
            dictionary mapping variable names or dtypes to interpolation method.
            Supported methods include:

            - `0` (nearest neighbor)
            - `1` (linear / bilinear)
            - `"nearest"`
            - `"triangular"`
            - `"bilinear"`

            The default is `0` for integer arrays, else `1`.
        agg_methods: Optional aggregation methods for downsampling spatial variables.
            Can be a single method for all variables or a dictionary mapping variable
            names or dtypes to methods. Supported methods include:
                "center", "count", "first", "last", "max", "mean", "median",
                "mode", "min", "prod", "std", "sum", and "var".
            Defaults to "center" for integer arrays, else "mean".
        recover_nans: Optional boolean or mapping to enable NaN recovery during
            upsampling (only applies when interpolation method is not nearest).
            Can be a single boolean or a dictionary mapping variable names or dtypes
            to booleans. Defaults to False.
        fill_values: Optional fill value(s) for areas outside input coverage.
            Can be a single value or dictionary by variable or type. If not provided,
            defaults based on data type are used:

            - float: NaN
            - uint8: 255
            - uint16: 65535
            - other ints: -1

        tile_size: Optional tile size for inferring a regular grid, if `target_gm` is
            not provided.

    Returns:
        A new dataset with spatial variables rectified to a regular grid.
            Variables not having 2D spatial dimensions are copied as-is. 1D spatial
            coordinate variables are ignored in the output.
    """
    if source_gm is None:
        source_gm = GridMapping.from_dataset(source_ds)
    source_ds = normalize_grid_mapping(source_ds, source_gm)

    if target_gm is None:
        target_gm = source_gm.to_regular(tile_size=tile_size)

    # transform 2d spatial coordinate of source dataset to target CRS
    if not _is_equal_crs(source_gm, target_gm):
        source_ds = _transform_coords(source_ds, source_gm, target_gm)
        source_gm = GridMapping.from_dataset(source_ds)

    source_ds = _select_variables(source_ds, variables)

    # ToDo: clip dataset

    # If source has higher resolution than target, downscale first, then rectify
    source_ds, source_gm = _downscale_source_dataset(
        source_ds,
        source_gm,
        target_gm,
        interp_methods,
        agg_methods,
        recover_nans,
    )

    # calculate source indices in target grid-mapping
    target_source_ij = _compute_target_source_ij(source_gm, target_gm, UV_DELTA)

    # rectify dataset
    x_name, y_name = source_gm.xy_var_names
    coords = source_ds.coords.to_dataset()
    coords = coords.drop_vars((x_name, y_name))
    x_name, y_name = target_gm.xy_var_names
    target_coords = target_gm.to_coords()
    coords[x_name] = target_coords[x_name]
    coords[y_name] = target_coords[y_name]
    coords["spatial_ref"] = xr.DataArray(0, attrs=target_gm.crs.to_cf())
    target_ds = xr.Dataset(coords=coords, attrs=source_ds.attrs)

    yx_dims = (source_gm.xy_dim_names[1], source_gm.xy_dim_names[0])
    for var_name, data_array in source_ds.data_vars.items():
        if data_array.dims[-2:] == yx_dims:
            assert len(data_array.dims) in (
                2,
                3,
            ), f"Data variable {var_name} has {len(data_array.dims)} dimensions."

            target_ds[var_name] = _rectify_data_array(
                data_array,
                var_name,
                target_gm,
                target_source_ij,
                interp_methods,
                fill_values,
            )

        elif yx_dims[0] not in data_array.dims and yx_dims[1] not in data_array.dims:
            target_ds[var_name] = data_array

    return target_ds

xcube_resampling.gridmapping.GridMapping

Bases: ABC

An abstract base class for grid mappings that define an image grid and a transformation from image pixel coordinates to spatial Earth coordinates defined in a well-known coordinate reference system (CRS).

This class cannot be instantiated directly. Use one of its factory methods to create instances:

  • :meth:regular
  • :meth:from_dataset
  • :meth:from_coords

Some instance methods can be used to derive new instances:

  • :meth:derive
  • :meth:scale
  • :meth:transform
  • :meth:to_regular

This class is thread-safe.

Source code in xcube_resampling/gridmapping/base.py
 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
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
class GridMapping(abc.ABC):
    """An abstract base class for grid mappings that define an image grid and
    a transformation from image pixel coordinates to spatial Earth coordinates
    defined in a well-known coordinate reference system (CRS).

    This class cannot be instantiated directly. Use one of its factory methods
    to create instances:

    * :meth:`regular`
    * :meth:`from_dataset`
    * :meth:`from_coords`

    Some instance methods can be used to derive new instances:

    * :meth:`derive`
    * :meth:`scale`
    * :meth:`transform`
    * :meth:`to_regular`

    This class is thread-safe.

    """

    def __init__(
        self,
        /,
        size: int | tuple[int, int],
        tile_size: int | tuple[int, int] | None,
        xy_bbox: tuple[FloatInt, FloatInt, FloatInt, FloatInt],
        xy_res: FloatInt | tuple[FloatInt, FloatInt],
        crs: pyproj.crs.CRS,
        xy_var_names: tuple[str, str],
        xy_dim_names: tuple[str, str],
        is_regular: bool | None = None,
        is_lon_360: bool | None = None,
        is_j_axis_up: bool | None = None,
        x_coords: xr.DataArray | None = None,
        y_coords: xr.DataArray | None = None,
    ):
        width, height = _normalize_int_pair(size, name="size")
        assert_true(width > 1 and height > 1, "invalid size")

        tile_width, tile_height = _normalize_int_pair(
            tile_size, default=(width, height)
        )
        assert_true(tile_width > 1 and tile_height > 1, "invalid tile_size")

        assert_given(xy_bbox, name="xy_bbox")
        assert_given(xy_res, name="xy_res")
        _assert_valid_xy_names(xy_var_names, name="xy_var_names")
        _assert_valid_xy_names(xy_dim_names, name="xy_dim_names")
        assert_instance(crs, pyproj.crs.CRS, name="crs")

        if x_coords is not None:
            assert_instance(x_coords, xr.DataArray, name="x_coords")
            assert_true(
                x_coords.ndim in (1, 2),
                message=f"x_coords.ndim must be 1 or 2, was {x_coords.ndim}",
            )
        if y_coords is not None:
            assert_instance(y_coords, xr.DataArray, name="y_coords")
            assert_true(
                y_coords.ndim in (1, 2),
                message=f"y_coords.ndim must be 1 or 2, was {y_coords.ndim}",
            )

        x_min, y_min, x_max, y_max = xy_bbox
        x_res, y_res = _normalize_number_pair(xy_res, name="xy_res")
        assert_true(x_res > 0 and y_res > 0, "invalid xy_res")

        self._lock = threading.RLock()

        self._size = width, height
        self._tile_size = tile_width, tile_height
        self._xy_bbox = x_min, y_min, x_max, y_max
        self._xy_res = x_res, y_res
        self._crs = crs
        self._xy_var_names = xy_var_names
        self._xy_dim_names = xy_dim_names
        self._is_regular = is_regular
        self._is_lon_360 = is_lon_360
        self._is_j_axis_up = is_j_axis_up
        self._x_coords = x_coords
        self._y_coords = y_coords
        self._xy_coords = None

    def derive(
        self,
        /,
        xy_var_names: tuple[str, str] = None,
        xy_dim_names: tuple[str, str] = None,
        tile_size: int | tuple[int, int] = None,
        is_j_axis_up: bool = None,
    ) -> "GridMapping":
        """Derive a new grid mapping from this one with some properties changed.

        Args:
            xy_var_names: The new x-, and y-variable names.
            xy_dim_names: The new x-, and y-dimension names.
            tile_size: The new tile size
            is_j_axis_up: Whether j-axis points up.

        Returns:
            A new, derived grid mapping.
        """
        other = copy.copy(self)
        if xy_var_names is not None:
            _assert_valid_xy_names(xy_var_names, name="xy_var_names")
            other._xy_var_names = xy_var_names
        if xy_dim_names is not None:
            _assert_valid_xy_names(xy_dim_names, name="xy_dim_names")
            other._xy_dim_names = xy_dim_names
        if tile_size is not None:
            tile_width, tile_height = _normalize_int_pair(tile_size, name="tile_size")
            assert_true(tile_width > 1 and tile_height > 1, "invalid tile_size")
            tile_size = tile_width, tile_height
            if other.tile_size != tile_size:
                other._tile_size = tile_width, tile_height
                with self._lock:
                    # if other._xy_coords has not been initialized before, we will do it
                    # in the next line. Otherwise, the following lines raise an error
                    if other._xy_coords is None:
                        _ = other.xy_coords
                    other._xy_coords = other._xy_coords.chunk(
                        {
                            dim: size
                            for (dim, size) in zip(
                                other._xy_coords.dims, other.xy_coords_chunks
                            )
                        }
                    )
        if is_j_axis_up is not None and is_j_axis_up != other._is_j_axis_up:
            other._is_j_axis_up = is_j_axis_up
            if other._y_coords is not None:
                other._y_coords = other._y_coords[::-1]
            if other._xy_coords is not None:
                other._xy_coords = other._xy_coords[:, ::-1, :]
                other._xy_coords = other._xy_coords.chunk(
                    {
                        dim: size
                        for (dim, size) in zip(
                            other._xy_coords.dims, other.xy_coords_chunks
                        )
                    }
                )

        return other

    def scale(
        self,
        xy_scale: FloatInt | tuple[FloatInt, FloatInt],
        tile_size: int | tuple[int, int] | None = None,
    ) -> "GridMapping":
        """Derive a scaled version of this regular grid mapping.

        Scaling factors larger than one correspond to up-scaling
        (pixels sizes decrease, image size increases).

        Scaling factors lower than one correspond to down-scaling.
        (pixels sizes increase, image size decreases).

        Args:
            xy_scale: The x-, and y-scaling factors. May be a single
                number or tuple.
            tile_size: The new tile size

        Returns:
            A new, scaled grid mapping.
        """
        self._assert_regular()
        x_scale, y_scale = _normalize_number_pair(xy_scale)
        new_xy_res, new_size = scale_xy_res_and_size(
            self.xy_res, self.size, (x_scale, y_scale)
        )
        if tile_size is not None:
            tile_width, tile_height = _normalize_int_pair(tile_size, name="tile_size")
        else:
            tile_width, tile_height = self.tile_size
        tile_width = min(new_size[0], tile_width)
        tile_height = min(new_size[1], tile_height)
        return self.regular(
            new_size,
            (self.x_min, self.y_min),
            new_xy_res,
            self.crs,
            tile_size=(tile_width, tile_height),
            is_j_axis_up=self.is_j_axis_up,
        ).derive(xy_dim_names=self.xy_dim_names, xy_var_names=self.xy_var_names)

    @property
    def size(self) -> tuple[int, int]:
        """Image size (width, height) in pixels."""
        return self._size

    @property
    def width(self) -> int:
        """Image width in pixels."""
        return self.size[0]

    @property
    def height(self) -> int:
        """Image height in pixels."""
        return self.size[1]

    @property
    def tile_size(self) -> tuple[int, int]:
        """Image tile size (width, height) in pixels."""
        return self._tile_size

    @property
    def is_tiled(self) -> bool:
        """Whether the image is tiled."""
        return self.size != self.tile_size

    @property
    def tile_width(self) -> int:
        """Image tile width in pixels."""
        return self.tile_size[0]

    @property
    def tile_height(self) -> int:
        """Image tile height in pixels."""
        return self.tile_size[1]

    @property
    def x_coords(self):
        """The 1D or 2D x-coordinate array of
        shape (width,) or (height, width).
        """
        return self._get_computed_attribute("_x_coords", self._new_x_coords)

    @abc.abstractmethod
    def _new_x_coords(self) -> xr.DataArray:
        """Create new 1D or 2D x-coordinate array of
        shape (width,) or (height, width).
        """

    @property
    def y_coords(self):
        """The 1D or 2D y-coordinate array of
        shape (width,) or (height, width).
        """
        return self._get_computed_attribute("_y_coords", self._new_y_coords)

    @abc.abstractmethod
    def _new_y_coords(self) -> xr.DataArray:
        """Create new 1D or 2D y-coordinate array of
        shape (width,) or (height, width).
        """

    @property
    def xy_coords(self) -> xr.DataArray:
        """The x,y coordinates as data array of shape (2, height, width).
        Coordinates are given in units of the CRS.
        """
        xy_coords = self._get_computed_attribute("_xy_coords", self._new_xy_coords)
        _assert_valid_xy_coords(xy_coords)
        return xy_coords

    @property
    def xy_coords_chunks(self) -> tuple[int, int, int]:
        """Get the chunks for the *xy_coords* array."""
        return 2, self.tile_height, self.tile_width

    @abc.abstractmethod
    def _new_xy_coords(self) -> xr.DataArray:
        """Create new coordinate array of shape (2, height, width)."""

    def _get_computed_attribute(self, name: str, computer: Callable[[], Any]) -> Any:
        """Get the value for a computed attribute.
        Utility to be used by this and derived classes.
        """
        value = getattr(self, name)
        if value is not None:
            return value
        with self._lock:
            # Double null check
            value = getattr(self, name)
            if value is not None:
                return value
            value = computer()
            setattr(self, name, value)
            return value

    @property
    def xy_var_names(self) -> tuple[str, str]:
        """The variable names of the x,y coordinates as
        tuple (x_var_name, y_var_name).
        """
        return self._xy_var_names

    @property
    def xy_dim_names(self) -> tuple[str, str]:
        """The dimension names of the x,y coordinates as
        tuple (x_dim_name, y_dim_name).
        """
        return self._xy_dim_names

    @property
    def xy_bbox(self) -> tuple[float, float, float, float]:
        """The image's bounding box in CRS coordinates."""
        return self._xy_bbox

    @property
    def x_min(self) -> FloatInt:
        """Minimum x-coordinate in CRS units."""
        return self._xy_bbox[0]

    @property
    def y_min(self) -> FloatInt:
        """Minimum y-coordinate in CRS units."""
        return self._xy_bbox[1]

    @property
    def x_max(self) -> FloatInt:
        """Maximum x-coordinate in CRS units."""
        return self._xy_bbox[2]

    @property
    def y_max(self) -> FloatInt:
        """Maximum y-coordinate in CRS units."""
        return self._xy_bbox[3]

    @property
    def xy_res(self) -> tuple[FloatInt, FloatInt]:
        """Pixel size in x and y direction."""
        return self._xy_res

    @property
    def x_res(self) -> FloatInt:
        """Pixel size in CRS units per pixel in x-direction."""
        return self._xy_res[0]

    @property
    def y_res(self) -> FloatInt:
        """Pixel size in CRS units per pixel in y-direction."""
        return self._xy_res[1]

    @property
    def crs(self) -> pyproj.crs.CRS:
        """The coordinate reference system."""
        return self._crs

    @property
    def spatial_unit_name(self) -> str:
        return self._crs.axis_info[0].unit_name

    @property
    def is_lon_360(self) -> bool | None:
        """Check whether *x_max* is greater than 180 degrees.
        Effectively tests whether the range *x_min*, *x_max* crosses
        the anti-meridian at 180 degrees.
        Works only for geographical coordinate reference systems.
        """
        return self._is_lon_360

    @property
    def is_regular(self) -> bool | None:
        """Do the x,y coordinates for a regular grid?
        A regular grid has a constant delta in both
        x- and y-directions of the x- and y-coordinates.

        Returns: None, if this property cannot be determined,
            True or False otherwise.
        """
        return self._is_regular

    @property
    def is_j_axis_up(self) -> bool | None:
        """Does the positive image j-axis point up?
        By default, the positive image j-axis points down.

        Returns: None, if this property cannot be determined,
            True or False otherwise.
        """
        return self._is_j_axis_up

    @property
    def ij_to_xy_transform(self) -> AffineTransformMatrix:
        """The affine transformation matrix from image to CRS coordinates.
        Defined only for grid mappings with rectified x,y coordinates.
        """
        self._assert_regular()
        if self.is_j_axis_up:
            return (
                (self.x_res, 0.0, self.x_min),
                (0.0, self.y_res, self.y_min),
            )
        else:
            return (
                (self.x_res, 0.0, self.x_min),
                (0.0, -self.y_res, self.y_max),
            )

    @property
    def xy_to_ij_transform(self) -> AffineTransformMatrix:
        """The affine transformation matrix from CRS to image coordinates.
        Defined only for grid mappings with rectified x,y coordinates.
        """
        self._assert_regular()
        return _from_affine(~_to_affine(self.ij_to_xy_transform))

    def ij_transform_to(self, other: "GridMapping") -> AffineTransformMatrix:
        """Get the affine transformation matrix that transforms
        image coordinates of *other* into image coordinates
        of this image geometry.

        Defined only for grid mappings with rectified x,y coordinates.

        Args:
            other: The other image geometry

        Returns:
            Affine transformation matrix
        """
        self._assert_regular()
        self.assert_regular(other, name="other")
        a = _to_affine(self.ij_to_xy_transform)
        b = _to_affine(other.xy_to_ij_transform)
        return _from_affine(b * a)

    def ij_transform_from(self, other: "GridMapping") -> AffineTransformMatrix:
        """Get the affine transformation matrix that transforms
        image coordinates of this image geometry to image
        coordinates of *other*.

        Defined only for grid mappings with rectified x,y coordinates.

        Args:
            other: The other image geometry

        Returns:
            Affine transformation matrix
        """
        self._assert_regular()
        self.assert_regular(other, name="other")
        a = _to_affine(self.ij_transform_to(other))
        return _from_affine(~a)

    @property
    def ij_bbox(self) -> tuple[int, int, int, int]:
        """The image's bounding box in pixel coordinates."""
        return 0, 0, self.width, self.height

    @property
    def ij_bboxes(self) -> np.ndarray:
        """The image tiles' bounding boxes in image pixel coordinates."""
        chunk_sizes = get_chunk_sizes(
            (self.height, self.width), (self.tile_height, self.tile_width)
        )
        _, _, block_slices = get_block_iterators(chunk_sizes)
        block_slices = tuple(block_slices)
        n = len(block_slices)
        ij_bboxes = np.ndarray((n, 4), dtype=np.int64)
        for i in range(n):
            y_slice, x_slice = block_slices[i]
            ij_bboxes[i, 0] = x_slice.start
            ij_bboxes[i, 1] = y_slice.start
            ij_bboxes[i, 2] = x_slice.stop
            ij_bboxes[i, 3] = y_slice.stop
        return ij_bboxes

    @property
    def xy_bboxes(self) -> np.ndarray:
        """The image tiles' bounding boxes in CRS coordinates."""
        if self.is_j_axis_up:
            xy_offset = np.array([self.x_min, self.y_min, self.x_min, self.y_min])
            xy_scale = np.array([self.x_res, self.y_res, self.x_res, self.y_res])
            xy_bboxes = xy_offset + xy_scale * self.ij_bboxes
        else:
            xy_offset = np.array([self.x_min, self.y_max, self.x_min, self.y_max])
            xy_scale = np.array([self.x_res, -self.y_res, self.x_res, -self.y_res])
            xy_bboxes = xy_offset + xy_scale * self.ij_bboxes
            xy_bboxes[:, [1, 3]] = xy_bboxes[:, [3, 1]]
        return xy_bboxes

    def ij_bbox_from_xy_bbox(
        self,
        xy_bbox: tuple[float, float, float, float],
        xy_border: float = 0.0,
        ij_border: int = 0,
    ) -> tuple[int, int, int, int]:
        """Compute bounding box in i,j pixel coordinates given a
        bounding box *xy_bbox* in x,y coordinates.

        Args:
            xy_bbox: Box (x_min, y_min, x_max, y_max) given in the same
                CS as x and y.
            xy_border: If non-zero, grows the bounding box *xy_bbox*
                before using it for comparisons. Defaults to 0.
            ij_border: If non-zero, grows the returned i,j bounding box
                and clips it to size. Defaults to 0.

        Returns:
            Bounding box in (i_min, j_min, i_max, j_max) in pixel
            coordinates. Returns ``(-1, -1, -1, -1)`` if *xy_bbox* isn't
            intersecting any of the x,y coordinates.
        """
        xy_bboxes = np.array([xy_bbox], dtype=np.float64)
        ij_bboxes = np.full_like(xy_bboxes, -1, dtype=np.int64)
        self.ij_bboxes_from_xy_bboxes(
            xy_bboxes, xy_border=xy_border, ij_border=ij_border, ij_bboxes=ij_bboxes
        )
        # noinspection PyTypeChecker
        return tuple(map(int, ij_bboxes[0]))

    def ij_bboxes_from_xy_bboxes(
        self,
        xy_bboxes: np.ndarray,
        xy_border: float = 0.0,
        ij_border: int = 0,
        ij_bboxes: np.ndarray = None,
    ) -> np.ndarray:
        """Compute bounding boxes in pixel coordinates given bounding boxes
        *xy_bboxes* [[x_min, y_min, x_max, y_max], ...] in x,y coordinates.

        The returned array in i,j pixel coordinates
        has the same shape as *xy_bboxes*. The value ranges in the
        returned array [[i_min, j_min, i_max, j_max], ..]] are:

        * i_min from 0 to width-1, i_max from 1 to width;
        * j_min from 0 to height-1, j_max from 1 to height;

        so the i,j pixel coordinates can be used as array index slices.

        Args:
            xy_bboxes: Numpy array of x,y bounding boxes [[x_min, y_min,
                x_max, y_max], ...] given in the same CS as x and y.
            xy_border: If non-zero, grows the bounding box *xy_bbox*
                before using it for comparisons. Defaults to 0.
            ij_border: If non-zero, grows the returned i,j bounding box
                and clips it to size. Defaults to 0.
            ij_bboxes: Numpy array of pixel i,j bounding boxes [[x_min,
                y_min, x_max, y_max], ...]. If given, must have same
                shape as *xy_bboxes*.

        Returns:
            Bounding boxes in [[i_min, j_min, i_max, j_max], ..]] in
            pixel coordinates.
        """
        if ij_bboxes is None:
            ij_bboxes = np.full_like(xy_bboxes, -1, dtype=np.int64)
        else:
            ij_bboxes[:, :] = -1
        xy_coords = self.xy_coords
        return self._compute_ij_bboxes_dask(
            xy_coords[0], xy_coords[1], xy_bboxes, xy_border, ij_border, ij_bboxes
        )

    def _compute_ij_bboxes_dask(
        self,
        x_coords: xr.DataArray,
        y_coords: xr.DataArray,
        xy_bboxes: np.ndarray,
        xy_border: float,
        ij_border: int,
        ij_bboxes: np.ndarray,
    ):
        from .bboxes import compute_ij_bboxes

        da.map_blocks(
            compute_ij_bboxes,
            x_coords.values,
            y_coords.values,
            xy_bboxes,
            xy_border,
            ij_border,
            ij_bboxes,
            dtype=ij_bboxes.dtype,
        ).compute()
        return ij_bboxes

    def to_coords(
        self,
        xy_var_names: tuple[str, str] = None,
        xy_dim_names: tuple[str, str] = None,
        exclude_bounds: bool = False,
        reuse_coords: bool = False,
    ) -> Mapping[str, xr.DataArray]:
        """Get CF-compliant axis coordinate variables and cell boundary
        coordinate variables.

        Defined only for grid mappings with regular x,y coordinates.

        Args:
            xy_var_names: Optional coordinate variable names
                (x_var_name, y_var_name).
            xy_dim_names: Optional coordinate dimensions names
                (x_dim_name, y_dim_name).
            exclude_bounds: If True, do not create bounds coordinates.
                Defaults to False.
            reuse_coords: Whether to either reuse target coordinate
                arrays from target_gm or to compute new ones.

        Returns:
            dictionary with coordinate variables
        """
        self._assert_regular()
        from .coords import grid_mapping_to_coords

        return grid_mapping_to_coords(
            self,
            xy_var_names=xy_var_names,
            xy_dim_names=xy_dim_names,
            exclude_bounds=exclude_bounds,
            reuse_coords=reuse_coords,
        )

    def transform(
        self,
        crs: str | pyproj.crs.CRS,
        *,
        xy_res: FloatInt | tuple[FloatInt, FloatInt] = None,
        tile_size: int | tuple[int, int] = None,
        xy_var_names: tuple[str, str] = None,
        tolerance: float = DEFAULT_TOLERANCE,
    ) -> "GridMapping":
        """Transform this grid mapping so it uses the given
        spatial coordinate reference system into another *crs*.

        Args:
            crs: The new spatial coordinate reference system.
            xy_res: Optional resolution in x- and y-directions.
                If given, speeds up the method by avoiding time-consuming
                spatial resolution estimation.
            tile_size: Optional new tile size.
            xy_var_names: Optional new coordinate names.
            tolerance: Absolute tolerance used when comparing
                coordinates with each other. Must be in the units of the
                *crs* and must be greater zero.

        Returns:
            A new grid mapping that uses *crs*.
        """
        from .transform import transform_grid_mapping

        return transform_grid_mapping(
            self,
            crs,
            xy_res=xy_res,
            tile_size=tile_size,
            xy_var_names=xy_var_names,
            tolerance=tolerance,
        )

    @classmethod
    def regular(
        cls,
        size: int | tuple[int, int],
        xy_min: tuple[float, float],
        xy_res: float | tuple[float, float],
        crs: str | pyproj.crs.CRS,
        *,
        tile_size: int | tuple[int, int] = None,
        is_j_axis_up: bool = False,
    ) -> "GridMapping":
        """Create a new regular grid mapping.

        Args:
            size: Size in pixels.
            xy_min: Minimum x- and y-coordinates.
            xy_res: Resolution in x- and y-directions.
            crs: Spatial coordinate reference system.
            tile_size: Optional tile size.
            is_j_axis_up: Whether positive j-axis points up. Defaults to
                false.

        Returns:
            A new regular grid mapping.
        """
        from .regular import new_regular_grid_mapping

        return new_regular_grid_mapping(
            size=size,
            xy_min=xy_min,
            xy_res=xy_res,
            crs=crs,
            tile_size=tile_size,
            is_j_axis_up=is_j_axis_up,
        )

    def to_regular(
        self, tile_size: int | tuple[int, int] | None = None, is_j_axis_up: bool = False
    ) -> "GridMapping":
        """Transform this grid mapping into one that is regular.

        Args:
            tile_size: Optional tile size.
            is_j_axis_up: Whether positive j-axis points up. Defaults to
                false.

        Returns:
            A new regular grid mapping or this grid mapping, if it is
            already regular.
        """
        from .regular import to_regular_grid_mapping

        return to_regular_grid_mapping(
            self, tile_size=tile_size, is_j_axis_up=is_j_axis_up
        )

    @classmethod
    def from_dataset(
        cls,
        dataset: xr.Dataset,
        *,
        crs: str | pyproj.crs.CRS | None = None,
        tile_size: int | tuple[int, int] | None = None,
        prefer_is_regular: bool = True,
        prefer_crs: str | pyproj.crs.CRS | None = None,
        emit_warnings: bool = False,
        tolerance: float = DEFAULT_TOLERANCE,
    ) -> "GridMapping":
        """Create a grid mapping for the given *dataset*.

        Args:
            dataset: The dataset.
            crs: Optional spatial coordinate reference system.
            tile_size: Optional tile size
            prefer_is_regular: Whether to prefer a regular grid mapping
                if multiple found. Default is True.
            prefer_crs: The preferred CRS of a grid mapping if multiple
                found.
            emit_warnings: Whether to emit warning for non-CF compliant
                datasets.
            tolerance: Absolute tolerance used when comparing
                coordinates with each other. Must be in the units of the
                *crs* and must be greater zero.

        Returns:
            a new grid mapping instance.
        """
        from .dataset import new_grid_mapping_from_dataset

        return new_grid_mapping_from_dataset(
            dataset=dataset,
            crs=crs,
            tile_size=tile_size,
            prefer_is_regular=prefer_is_regular,
            prefer_crs=prefer_crs,
            emit_warnings=emit_warnings,
            tolerance=tolerance,
        )

    @classmethod
    def from_coords(
        cls,
        x_coords: xr.DataArray,
        y_coords: xr.DataArray,
        crs: str | pyproj.crs.CRS,
        *,
        tile_size: int | tuple[int, int] | None = None,
        tolerance: float = DEFAULT_TOLERANCE,
    ) -> "GridMapping":
        """Create a grid mapping from given x- and y-coordinates
        *x_coords*, *y_coords* and spatial coordinate reference
        system *crs*.

        Args:
            x_coords: The x-coordinates.
            y_coords: The y-coordinates.
            crs: The spatial coordinate reference system.
            tile_size: Optional tile size.
            tolerance: Absolute tolerance used when comparing
                coordinates with each other. Must be in the units of the
                *crs* and must be greater zero.

        Returns:
            A new grid mapping.
        """
        from .coords import new_grid_mapping_from_coords

        return new_grid_mapping_from_coords(
            x_coords=x_coords,
            y_coords=y_coords,
            crs=crs,
            tile_size=tile_size,
            tolerance=tolerance,
        )

    def is_close(
        self, other: "GridMapping", tolerance: float = DEFAULT_TOLERANCE
    ) -> bool:
        """Tests whether this grid mapping is close to *other*.

        Args:
            other: The other grid mapping.
            tolerance: Absolute tolerance used when comparing
                coordinates with each other. Must be in the units of the
                *crs* and must be greater zero.

        Returns:
            True, if so, False otherwise.
        """
        if self is other:
            return True
        if (
            self.is_j_axis_up == other.is_j_axis_up
            and self.is_lon_360 == other.is_lon_360
            and self.is_regular == other.is_regular
            and self.size == other.size
            and self.tile_size == other.tile_size
            and self.crs == other.crs
        ):
            sxr, syr = self.xy_res
            oxr, oyr = other.xy_res
            if math.isclose(sxr, oxr, abs_tol=tolerance) and math.isclose(
                syr, oyr, abs_tol=tolerance
            ):
                sx1, sy1, sx2, sy2 = self.xy_bbox
                ox1, oy1, ox2, oy2 = other.xy_bbox
                return (
                    math.isclose(sx1, ox1, abs_tol=tolerance)
                    and math.isclose(sy1, oy1, abs_tol=tolerance)
                    and math.isclose(sx2, ox2, abs_tol=tolerance)
                    and math.isclose(sy2, oy2, abs_tol=tolerance)
                )
        return False

    @classmethod
    def assert_regular(cls, value: Any, name: str = None):
        assert_instance(value, GridMapping, name=name)
        if not value.is_regular:
            raise ValueError(f"{name or 'value'} must be a regular grid mapping")

    def _assert_regular(self):
        if not self.is_regular:
            raise NotImplementedError(
                "Operation not implemented for non-regular grid mappings"
            )

    def _repr_markdown_(self) -> str:
        """Generate an IPython Notebook Markdown representation."""
        is_regular = self.is_regular if self.is_regular is not None else "_unknown_"
        is_j_axis_up = (
            self.is_j_axis_up if self.is_j_axis_up is not None else "_unknown_"
        )
        is_lon_360 = self.is_lon_360 if self.is_lon_360 is not None else "_unknown_"
        xy_res = repr(self.xy_res) + ("" if self.is_regular else "  _estimated_")
        return "\n".join(
            [
                f"class: **{self.__class__.__name__}**",
                f"* is_regular: {is_regular}",
                f"* is_j_axis_up: {is_j_axis_up}",
                f"* is_lon_360: {is_lon_360}",
                f"* crs: {self.crs}",
                f"* xy_res: {xy_res}",
                f"* xy_bbox: {self.xy_bbox}",
                f"* ij_bbox: {self.ij_bbox}",
                f"* xy_dim_names: {self.xy_dim_names}",
                f"* xy_var_names: {self.xy_var_names}",
                f"* size: {self.size}",
                f"* tile_size: {self.tile_size}",
            ]
        )

size property

Image size (width, height) in pixels.

width property

Image width in pixels.

height property

Image height in pixels.

tile_size property

Image tile size (width, height) in pixels.

is_tiled property

Whether the image is tiled.

tile_width property

Image tile width in pixels.

tile_height property

Image tile height in pixels.

x_coords property

The 1D or 2D x-coordinate array of shape (width,) or (height, width).

y_coords property

The 1D or 2D y-coordinate array of shape (width,) or (height, width).

xy_coords property

The x,y coordinates as data array of shape (2, height, width). Coordinates are given in units of the CRS.

xy_coords_chunks property

Get the chunks for the xy_coords array.

xy_var_names property

The variable names of the x,y coordinates as tuple (x_var_name, y_var_name).

xy_dim_names property

The dimension names of the x,y coordinates as tuple (x_dim_name, y_dim_name).

xy_bbox property

The image's bounding box in CRS coordinates.

x_min property

Minimum x-coordinate in CRS units.

y_min property

Minimum y-coordinate in CRS units.

x_max property

Maximum x-coordinate in CRS units.

y_max property

Maximum y-coordinate in CRS units.

xy_res property

Pixel size in x and y direction.

x_res property

Pixel size in CRS units per pixel in x-direction.

y_res property

Pixel size in CRS units per pixel in y-direction.

crs property

The coordinate reference system.

is_lon_360 property

Check whether x_max is greater than 180 degrees. Effectively tests whether the range x_min, x_max crosses the anti-meridian at 180 degrees. Works only for geographical coordinate reference systems.

is_regular property

Do the x,y coordinates for a regular grid? A regular grid has a constant delta in both x- and y-directions of the x- and y-coordinates.

None, if this property cannot be determined,

Type Description
bool | None

True or False otherwise.

is_j_axis_up property

Does the positive image j-axis point up? By default, the positive image j-axis points down.

None, if this property cannot be determined,

Type Description
bool | None

True or False otherwise.

ij_to_xy_transform property

The affine transformation matrix from image to CRS coordinates. Defined only for grid mappings with rectified x,y coordinates.

xy_to_ij_transform property

The affine transformation matrix from CRS to image coordinates. Defined only for grid mappings with rectified x,y coordinates.

ij_bbox property

The image's bounding box in pixel coordinates.

ij_bboxes property

The image tiles' bounding boxes in image pixel coordinates.

xy_bboxes property

The image tiles' bounding boxes in CRS coordinates.

derive(xy_var_names=None, xy_dim_names=None, tile_size=None, is_j_axis_up=None)

Derive a new grid mapping from this one with some properties changed.

Parameters:

Name Type Description Default
xy_var_names tuple[str, str]

The new x-, and y-variable names.

None
xy_dim_names tuple[str, str]

The new x-, and y-dimension names.

None
tile_size int | tuple[int, int]

The new tile size

None
is_j_axis_up bool

Whether j-axis points up.

None

Returns:

Type Description
GridMapping

A new, derived grid mapping.

Source code in xcube_resampling/gridmapping/base.py
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
def derive(
    self,
    /,
    xy_var_names: tuple[str, str] = None,
    xy_dim_names: tuple[str, str] = None,
    tile_size: int | tuple[int, int] = None,
    is_j_axis_up: bool = None,
) -> "GridMapping":
    """Derive a new grid mapping from this one with some properties changed.

    Args:
        xy_var_names: The new x-, and y-variable names.
        xy_dim_names: The new x-, and y-dimension names.
        tile_size: The new tile size
        is_j_axis_up: Whether j-axis points up.

    Returns:
        A new, derived grid mapping.
    """
    other = copy.copy(self)
    if xy_var_names is not None:
        _assert_valid_xy_names(xy_var_names, name="xy_var_names")
        other._xy_var_names = xy_var_names
    if xy_dim_names is not None:
        _assert_valid_xy_names(xy_dim_names, name="xy_dim_names")
        other._xy_dim_names = xy_dim_names
    if tile_size is not None:
        tile_width, tile_height = _normalize_int_pair(tile_size, name="tile_size")
        assert_true(tile_width > 1 and tile_height > 1, "invalid tile_size")
        tile_size = tile_width, tile_height
        if other.tile_size != tile_size:
            other._tile_size = tile_width, tile_height
            with self._lock:
                # if other._xy_coords has not been initialized before, we will do it
                # in the next line. Otherwise, the following lines raise an error
                if other._xy_coords is None:
                    _ = other.xy_coords
                other._xy_coords = other._xy_coords.chunk(
                    {
                        dim: size
                        for (dim, size) in zip(
                            other._xy_coords.dims, other.xy_coords_chunks
                        )
                    }
                )
    if is_j_axis_up is not None and is_j_axis_up != other._is_j_axis_up:
        other._is_j_axis_up = is_j_axis_up
        if other._y_coords is not None:
            other._y_coords = other._y_coords[::-1]
        if other._xy_coords is not None:
            other._xy_coords = other._xy_coords[:, ::-1, :]
            other._xy_coords = other._xy_coords.chunk(
                {
                    dim: size
                    for (dim, size) in zip(
                        other._xy_coords.dims, other.xy_coords_chunks
                    )
                }
            )

    return other

scale(xy_scale, tile_size=None)

Derive a scaled version of this regular grid mapping.

Scaling factors larger than one correspond to up-scaling (pixels sizes decrease, image size increases).

Scaling factors lower than one correspond to down-scaling. (pixels sizes increase, image size decreases).

Parameters:

Name Type Description Default
xy_scale FloatInt | tuple[FloatInt, FloatInt]

The x-, and y-scaling factors. May be a single number or tuple.

required
tile_size int | tuple[int, int] | None

The new tile size

None

Returns:

Type Description
GridMapping

A new, scaled grid mapping.

Source code in xcube_resampling/gridmapping/base.py
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
def scale(
    self,
    xy_scale: FloatInt | tuple[FloatInt, FloatInt],
    tile_size: int | tuple[int, int] | None = None,
) -> "GridMapping":
    """Derive a scaled version of this regular grid mapping.

    Scaling factors larger than one correspond to up-scaling
    (pixels sizes decrease, image size increases).

    Scaling factors lower than one correspond to down-scaling.
    (pixels sizes increase, image size decreases).

    Args:
        xy_scale: The x-, and y-scaling factors. May be a single
            number or tuple.
        tile_size: The new tile size

    Returns:
        A new, scaled grid mapping.
    """
    self._assert_regular()
    x_scale, y_scale = _normalize_number_pair(xy_scale)
    new_xy_res, new_size = scale_xy_res_and_size(
        self.xy_res, self.size, (x_scale, y_scale)
    )
    if tile_size is not None:
        tile_width, tile_height = _normalize_int_pair(tile_size, name="tile_size")
    else:
        tile_width, tile_height = self.tile_size
    tile_width = min(new_size[0], tile_width)
    tile_height = min(new_size[1], tile_height)
    return self.regular(
        new_size,
        (self.x_min, self.y_min),
        new_xy_res,
        self.crs,
        tile_size=(tile_width, tile_height),
        is_j_axis_up=self.is_j_axis_up,
    ).derive(xy_dim_names=self.xy_dim_names, xy_var_names=self.xy_var_names)

ij_transform_to(other)

Get the affine transformation matrix that transforms image coordinates of other into image coordinates of this image geometry.

Defined only for grid mappings with rectified x,y coordinates.

Parameters:

Name Type Description Default
other GridMapping

The other image geometry

required

Returns:

Type Description
AffineTransformMatrix

Affine transformation matrix

Source code in xcube_resampling/gridmapping/base.py
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
def ij_transform_to(self, other: "GridMapping") -> AffineTransformMatrix:
    """Get the affine transformation matrix that transforms
    image coordinates of *other* into image coordinates
    of this image geometry.

    Defined only for grid mappings with rectified x,y coordinates.

    Args:
        other: The other image geometry

    Returns:
        Affine transformation matrix
    """
    self._assert_regular()
    self.assert_regular(other, name="other")
    a = _to_affine(self.ij_to_xy_transform)
    b = _to_affine(other.xy_to_ij_transform)
    return _from_affine(b * a)

ij_transform_from(other)

Get the affine transformation matrix that transforms image coordinates of this image geometry to image coordinates of other.

Defined only for grid mappings with rectified x,y coordinates.

Parameters:

Name Type Description Default
other GridMapping

The other image geometry

required

Returns:

Type Description
AffineTransformMatrix

Affine transformation matrix

Source code in xcube_resampling/gridmapping/base.py
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
def ij_transform_from(self, other: "GridMapping") -> AffineTransformMatrix:
    """Get the affine transformation matrix that transforms
    image coordinates of this image geometry to image
    coordinates of *other*.

    Defined only for grid mappings with rectified x,y coordinates.

    Args:
        other: The other image geometry

    Returns:
        Affine transformation matrix
    """
    self._assert_regular()
    self.assert_regular(other, name="other")
    a = _to_affine(self.ij_transform_to(other))
    return _from_affine(~a)

ij_bbox_from_xy_bbox(xy_bbox, xy_border=0.0, ij_border=0)

Compute bounding box in i,j pixel coordinates given a bounding box xy_bbox in x,y coordinates.

Parameters:

Name Type Description Default
xy_bbox tuple[float, float, float, float]

Box (x_min, y_min, x_max, y_max) given in the same CS as x and y.

required
xy_border float

If non-zero, grows the bounding box xy_bbox before using it for comparisons. Defaults to 0.

0.0
ij_border int

If non-zero, grows the returned i,j bounding box and clips it to size. Defaults to 0.

0

Returns:

Type Description
int

Bounding box in (i_min, j_min, i_max, j_max) in pixel

int

coordinates. Returns (-1, -1, -1, -1) if xy_bbox isn't

int

intersecting any of the x,y coordinates.

Source code in xcube_resampling/gridmapping/base.py
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
def ij_bbox_from_xy_bbox(
    self,
    xy_bbox: tuple[float, float, float, float],
    xy_border: float = 0.0,
    ij_border: int = 0,
) -> tuple[int, int, int, int]:
    """Compute bounding box in i,j pixel coordinates given a
    bounding box *xy_bbox* in x,y coordinates.

    Args:
        xy_bbox: Box (x_min, y_min, x_max, y_max) given in the same
            CS as x and y.
        xy_border: If non-zero, grows the bounding box *xy_bbox*
            before using it for comparisons. Defaults to 0.
        ij_border: If non-zero, grows the returned i,j bounding box
            and clips it to size. Defaults to 0.

    Returns:
        Bounding box in (i_min, j_min, i_max, j_max) in pixel
        coordinates. Returns ``(-1, -1, -1, -1)`` if *xy_bbox* isn't
        intersecting any of the x,y coordinates.
    """
    xy_bboxes = np.array([xy_bbox], dtype=np.float64)
    ij_bboxes = np.full_like(xy_bboxes, -1, dtype=np.int64)
    self.ij_bboxes_from_xy_bboxes(
        xy_bboxes, xy_border=xy_border, ij_border=ij_border, ij_bboxes=ij_bboxes
    )
    # noinspection PyTypeChecker
    return tuple(map(int, ij_bboxes[0]))

ij_bboxes_from_xy_bboxes(xy_bboxes, xy_border=0.0, ij_border=0, ij_bboxes=None)

Compute bounding boxes in pixel coordinates given bounding boxes xy_bboxes [[x_min, y_min, x_max, y_max], ...] in x,y coordinates.

The returned array in i,j pixel coordinates has the same shape as xy_bboxes. The value ranges in the returned array [[i_min, j_min, i_max, j_max], ..]] are:

  • i_min from 0 to width-1, i_max from 1 to width;
  • j_min from 0 to height-1, j_max from 1 to height;

so the i,j pixel coordinates can be used as array index slices.

Parameters:

Name Type Description Default
xy_bboxes ndarray

Numpy array of x,y bounding boxes [[x_min, y_min, x_max, y_max], ...] given in the same CS as x and y.

required
xy_border float

If non-zero, grows the bounding box xy_bbox before using it for comparisons. Defaults to 0.

0.0
ij_border int

If non-zero, grows the returned i,j bounding box and clips it to size. Defaults to 0.

0
ij_bboxes ndarray

Numpy array of pixel i,j bounding boxes [[x_min, y_min, x_max, y_max], ...]. If given, must have same shape as xy_bboxes.

None

Returns:

Type Description
ndarray

Bounding boxes in [[i_min, j_min, i_max, j_max], ..]] in

ndarray

pixel coordinates.

Source code in xcube_resampling/gridmapping/base.py
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
def ij_bboxes_from_xy_bboxes(
    self,
    xy_bboxes: np.ndarray,
    xy_border: float = 0.0,
    ij_border: int = 0,
    ij_bboxes: np.ndarray = None,
) -> np.ndarray:
    """Compute bounding boxes in pixel coordinates given bounding boxes
    *xy_bboxes* [[x_min, y_min, x_max, y_max], ...] in x,y coordinates.

    The returned array in i,j pixel coordinates
    has the same shape as *xy_bboxes*. The value ranges in the
    returned array [[i_min, j_min, i_max, j_max], ..]] are:

    * i_min from 0 to width-1, i_max from 1 to width;
    * j_min from 0 to height-1, j_max from 1 to height;

    so the i,j pixel coordinates can be used as array index slices.

    Args:
        xy_bboxes: Numpy array of x,y bounding boxes [[x_min, y_min,
            x_max, y_max], ...] given in the same CS as x and y.
        xy_border: If non-zero, grows the bounding box *xy_bbox*
            before using it for comparisons. Defaults to 0.
        ij_border: If non-zero, grows the returned i,j bounding box
            and clips it to size. Defaults to 0.
        ij_bboxes: Numpy array of pixel i,j bounding boxes [[x_min,
            y_min, x_max, y_max], ...]. If given, must have same
            shape as *xy_bboxes*.

    Returns:
        Bounding boxes in [[i_min, j_min, i_max, j_max], ..]] in
        pixel coordinates.
    """
    if ij_bboxes is None:
        ij_bboxes = np.full_like(xy_bboxes, -1, dtype=np.int64)
    else:
        ij_bboxes[:, :] = -1
    xy_coords = self.xy_coords
    return self._compute_ij_bboxes_dask(
        xy_coords[0], xy_coords[1], xy_bboxes, xy_border, ij_border, ij_bboxes
    )

to_coords(xy_var_names=None, xy_dim_names=None, exclude_bounds=False, reuse_coords=False)

Get CF-compliant axis coordinate variables and cell boundary coordinate variables.

Defined only for grid mappings with regular x,y coordinates.

Parameters:

Name Type Description Default
xy_var_names tuple[str, str]

Optional coordinate variable names (x_var_name, y_var_name).

None
xy_dim_names tuple[str, str]

Optional coordinate dimensions names (x_dim_name, y_dim_name).

None
exclude_bounds bool

If True, do not create bounds coordinates. Defaults to False.

False
reuse_coords bool

Whether to either reuse target coordinate arrays from target_gm or to compute new ones.

False

Returns:

Type Description
Mapping[str, DataArray]

dictionary with coordinate variables

Source code in xcube_resampling/gridmapping/base.py
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
def to_coords(
    self,
    xy_var_names: tuple[str, str] = None,
    xy_dim_names: tuple[str, str] = None,
    exclude_bounds: bool = False,
    reuse_coords: bool = False,
) -> Mapping[str, xr.DataArray]:
    """Get CF-compliant axis coordinate variables and cell boundary
    coordinate variables.

    Defined only for grid mappings with regular x,y coordinates.

    Args:
        xy_var_names: Optional coordinate variable names
            (x_var_name, y_var_name).
        xy_dim_names: Optional coordinate dimensions names
            (x_dim_name, y_dim_name).
        exclude_bounds: If True, do not create bounds coordinates.
            Defaults to False.
        reuse_coords: Whether to either reuse target coordinate
            arrays from target_gm or to compute new ones.

    Returns:
        dictionary with coordinate variables
    """
    self._assert_regular()
    from .coords import grid_mapping_to_coords

    return grid_mapping_to_coords(
        self,
        xy_var_names=xy_var_names,
        xy_dim_names=xy_dim_names,
        exclude_bounds=exclude_bounds,
        reuse_coords=reuse_coords,
    )

transform(crs, *, xy_res=None, tile_size=None, xy_var_names=None, tolerance=DEFAULT_TOLERANCE)

Transform this grid mapping so it uses the given spatial coordinate reference system into another crs.

Parameters:

Name Type Description Default
crs str | CRS

The new spatial coordinate reference system.

required
xy_res FloatInt | tuple[FloatInt, FloatInt]

Optional resolution in x- and y-directions. If given, speeds up the method by avoiding time-consuming spatial resolution estimation.

None
tile_size int | tuple[int, int]

Optional new tile size.

None
xy_var_names tuple[str, str]

Optional new coordinate names.

None
tolerance float

Absolute tolerance used when comparing coordinates with each other. Must be in the units of the crs and must be greater zero.

DEFAULT_TOLERANCE

Returns:

Type Description
GridMapping

A new grid mapping that uses crs.

Source code in xcube_resampling/gridmapping/base.py
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
def transform(
    self,
    crs: str | pyproj.crs.CRS,
    *,
    xy_res: FloatInt | tuple[FloatInt, FloatInt] = None,
    tile_size: int | tuple[int, int] = None,
    xy_var_names: tuple[str, str] = None,
    tolerance: float = DEFAULT_TOLERANCE,
) -> "GridMapping":
    """Transform this grid mapping so it uses the given
    spatial coordinate reference system into another *crs*.

    Args:
        crs: The new spatial coordinate reference system.
        xy_res: Optional resolution in x- and y-directions.
            If given, speeds up the method by avoiding time-consuming
            spatial resolution estimation.
        tile_size: Optional new tile size.
        xy_var_names: Optional new coordinate names.
        tolerance: Absolute tolerance used when comparing
            coordinates with each other. Must be in the units of the
            *crs* and must be greater zero.

    Returns:
        A new grid mapping that uses *crs*.
    """
    from .transform import transform_grid_mapping

    return transform_grid_mapping(
        self,
        crs,
        xy_res=xy_res,
        tile_size=tile_size,
        xy_var_names=xy_var_names,
        tolerance=tolerance,
    )

regular(size, xy_min, xy_res, crs, *, tile_size=None, is_j_axis_up=False) classmethod

Create a new regular grid mapping.

Parameters:

Name Type Description Default
size int | tuple[int, int]

Size in pixels.

required
xy_min tuple[float, float]

Minimum x- and y-coordinates.

required
xy_res float | tuple[float, float]

Resolution in x- and y-directions.

required
crs str | CRS

Spatial coordinate reference system.

required
tile_size int | tuple[int, int]

Optional tile size.

None
is_j_axis_up bool

Whether positive j-axis points up. Defaults to false.

False

Returns:

Type Description
GridMapping

A new regular grid mapping.

Source code in xcube_resampling/gridmapping/base.py
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
@classmethod
def regular(
    cls,
    size: int | tuple[int, int],
    xy_min: tuple[float, float],
    xy_res: float | tuple[float, float],
    crs: str | pyproj.crs.CRS,
    *,
    tile_size: int | tuple[int, int] = None,
    is_j_axis_up: bool = False,
) -> "GridMapping":
    """Create a new regular grid mapping.

    Args:
        size: Size in pixels.
        xy_min: Minimum x- and y-coordinates.
        xy_res: Resolution in x- and y-directions.
        crs: Spatial coordinate reference system.
        tile_size: Optional tile size.
        is_j_axis_up: Whether positive j-axis points up. Defaults to
            false.

    Returns:
        A new regular grid mapping.
    """
    from .regular import new_regular_grid_mapping

    return new_regular_grid_mapping(
        size=size,
        xy_min=xy_min,
        xy_res=xy_res,
        crs=crs,
        tile_size=tile_size,
        is_j_axis_up=is_j_axis_up,
    )

to_regular(tile_size=None, is_j_axis_up=False)

Transform this grid mapping into one that is regular.

Parameters:

Name Type Description Default
tile_size int | tuple[int, int] | None

Optional tile size.

None
is_j_axis_up bool

Whether positive j-axis points up. Defaults to false.

False

Returns:

Type Description
GridMapping

A new regular grid mapping or this grid mapping, if it is

GridMapping

already regular.

Source code in xcube_resampling/gridmapping/base.py
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
def to_regular(
    self, tile_size: int | tuple[int, int] | None = None, is_j_axis_up: bool = False
) -> "GridMapping":
    """Transform this grid mapping into one that is regular.

    Args:
        tile_size: Optional tile size.
        is_j_axis_up: Whether positive j-axis points up. Defaults to
            false.

    Returns:
        A new regular grid mapping or this grid mapping, if it is
        already regular.
    """
    from .regular import to_regular_grid_mapping

    return to_regular_grid_mapping(
        self, tile_size=tile_size, is_j_axis_up=is_j_axis_up
    )

from_dataset(dataset, *, crs=None, tile_size=None, prefer_is_regular=True, prefer_crs=None, emit_warnings=False, tolerance=DEFAULT_TOLERANCE) classmethod

Create a grid mapping for the given dataset.

Parameters:

Name Type Description Default
dataset Dataset

The dataset.

required
crs str | CRS | None

Optional spatial coordinate reference system.

None
tile_size int | tuple[int, int] | None

Optional tile size

None
prefer_is_regular bool

Whether to prefer a regular grid mapping if multiple found. Default is True.

True
prefer_crs str | CRS | None

The preferred CRS of a grid mapping if multiple found.

None
emit_warnings bool

Whether to emit warning for non-CF compliant datasets.

False
tolerance float

Absolute tolerance used when comparing coordinates with each other. Must be in the units of the crs and must be greater zero.

DEFAULT_TOLERANCE

Returns:

Type Description
GridMapping

a new grid mapping instance.

Source code in xcube_resampling/gridmapping/base.py
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
@classmethod
def from_dataset(
    cls,
    dataset: xr.Dataset,
    *,
    crs: str | pyproj.crs.CRS | None = None,
    tile_size: int | tuple[int, int] | None = None,
    prefer_is_regular: bool = True,
    prefer_crs: str | pyproj.crs.CRS | None = None,
    emit_warnings: bool = False,
    tolerance: float = DEFAULT_TOLERANCE,
) -> "GridMapping":
    """Create a grid mapping for the given *dataset*.

    Args:
        dataset: The dataset.
        crs: Optional spatial coordinate reference system.
        tile_size: Optional tile size
        prefer_is_regular: Whether to prefer a regular grid mapping
            if multiple found. Default is True.
        prefer_crs: The preferred CRS of a grid mapping if multiple
            found.
        emit_warnings: Whether to emit warning for non-CF compliant
            datasets.
        tolerance: Absolute tolerance used when comparing
            coordinates with each other. Must be in the units of the
            *crs* and must be greater zero.

    Returns:
        a new grid mapping instance.
    """
    from .dataset import new_grid_mapping_from_dataset

    return new_grid_mapping_from_dataset(
        dataset=dataset,
        crs=crs,
        tile_size=tile_size,
        prefer_is_regular=prefer_is_regular,
        prefer_crs=prefer_crs,
        emit_warnings=emit_warnings,
        tolerance=tolerance,
    )

from_coords(x_coords, y_coords, crs, *, tile_size=None, tolerance=DEFAULT_TOLERANCE) classmethod

Create a grid mapping from given x- and y-coordinates x_coords, y_coords and spatial coordinate reference system crs.

Parameters:

Name Type Description Default
x_coords DataArray

The x-coordinates.

required
y_coords DataArray

The y-coordinates.

required
crs str | CRS

The spatial coordinate reference system.

required
tile_size int | tuple[int, int] | None

Optional tile size.

None
tolerance float

Absolute tolerance used when comparing coordinates with each other. Must be in the units of the crs and must be greater zero.

DEFAULT_TOLERANCE

Returns:

Type Description
GridMapping

A new grid mapping.

Source code in xcube_resampling/gridmapping/base.py
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
829
830
831
832
833
834
835
836
837
@classmethod
def from_coords(
    cls,
    x_coords: xr.DataArray,
    y_coords: xr.DataArray,
    crs: str | pyproj.crs.CRS,
    *,
    tile_size: int | tuple[int, int] | None = None,
    tolerance: float = DEFAULT_TOLERANCE,
) -> "GridMapping":
    """Create a grid mapping from given x- and y-coordinates
    *x_coords*, *y_coords* and spatial coordinate reference
    system *crs*.

    Args:
        x_coords: The x-coordinates.
        y_coords: The y-coordinates.
        crs: The spatial coordinate reference system.
        tile_size: Optional tile size.
        tolerance: Absolute tolerance used when comparing
            coordinates with each other. Must be in the units of the
            *crs* and must be greater zero.

    Returns:
        A new grid mapping.
    """
    from .coords import new_grid_mapping_from_coords

    return new_grid_mapping_from_coords(
        x_coords=x_coords,
        y_coords=y_coords,
        crs=crs,
        tile_size=tile_size,
        tolerance=tolerance,
    )

is_close(other, tolerance=DEFAULT_TOLERANCE)

Tests whether this grid mapping is close to other.

Parameters:

Name Type Description Default
other GridMapping

The other grid mapping.

required
tolerance float

Absolute tolerance used when comparing coordinates with each other. Must be in the units of the crs and must be greater zero.

DEFAULT_TOLERANCE

Returns:

Type Description
bool

True, if so, False otherwise.

Source code in xcube_resampling/gridmapping/base.py
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
def is_close(
    self, other: "GridMapping", tolerance: float = DEFAULT_TOLERANCE
) -> bool:
    """Tests whether this grid mapping is close to *other*.

    Args:
        other: The other grid mapping.
        tolerance: Absolute tolerance used when comparing
            coordinates with each other. Must be in the units of the
            *crs* and must be greater zero.

    Returns:
        True, if so, False otherwise.
    """
    if self is other:
        return True
    if (
        self.is_j_axis_up == other.is_j_axis_up
        and self.is_lon_360 == other.is_lon_360
        and self.is_regular == other.is_regular
        and self.size == other.size
        and self.tile_size == other.tile_size
        and self.crs == other.crs
    ):
        sxr, syr = self.xy_res
        oxr, oyr = other.xy_res
        if math.isclose(sxr, oxr, abs_tol=tolerance) and math.isclose(
            syr, oyr, abs_tol=tolerance
        ):
            sx1, sy1, sx2, sy2 = self.xy_bbox
            ox1, oy1, ox2, oy2 = other.xy_bbox
            return (
                math.isclose(sx1, ox1, abs_tol=tolerance)
                and math.isclose(sy1, oy1, abs_tol=tolerance)
                and math.isclose(sx2, ox2, abs_tol=tolerance)
                and math.isclose(sy2, oy2, abs_tol=tolerance)
            )
    return False