xarray_regrid.methods.conservative

Conservative regridding implementation.

Attributes

sparse

Functions

conservative_regrid(…)

Refine a dataset using conservative regridding.

conservative_regrid_dataset(→ xarray.Dataset)

Dataset implementation of the conservative regridding method.

apply_weights(→ xarray.DataArray)

Apply the weights over all regridding dimensions simultaneously with xr.dot.

get_valid_threshold(→ float)

Invert the nan_threshold and coerce it to just above zero and below

get_weights(→ numpy.ndarray)

Determine the weights to map from the old coordinates to the new coordinates.

apply_spherical_correction(→ xarray.DataArray)

Apply a sperical earth correction on the prepared dot product weights.

lat_weight(→ numpy.ndarray)

Return the weight of gridcells based on their latitude.

format_weights(→ xarray.DataArray)

Format the raw weights array such that:

Module Contents

xarray_regrid.methods.conservative.sparse = None[source]
xarray_regrid.methods.conservative.conservative_regrid(data: xarray.DataArray, target_ds: xarray.Dataset, latitude_coord: str | None, skipna: bool = True, nan_threshold: float = 1.0, output_chunks: dict[collections.abc.Hashable, int] | None = None) xarray.DataArray[source]
xarray_regrid.methods.conservative.conservative_regrid(data: xarray.Dataset, target_ds: xarray.Dataset, latitude_coord: str | None, skipna: bool = True, nan_threshold: float = 1.0, output_chunks: dict[collections.abc.Hashable, int] | None = None) xarray.Dataset

Refine a dataset using conservative regridding.

The method implementation is based on a post by Stephan Hoyer; “For the case of interpolation between rectilinear grids (even on the sphere), you can factorize regridding along each axis. This is less general but makes the entire calculation much simpler, because its feasible to store interpolation weights as dense matrices and to use dense matrix multiplication.” https://discourse.pangeo.io/t/conservative-region-aggregation-with-xarray-geopandas-and-sparse/2715/3

Parameters:
  • data – Input dataset.

  • target_ds – Dataset which coordinates the input dataset should be regrid to.

  • latitude_coord – Name of the latitude coordinate. If not provided, attempt to infer it from the first coordinate equaling the string ‘lat’ or ‘latitude’.

  • skipna – If True, enable handling for NaN values. This adds some overhead, so should be disabled for optimal performance on data without NaNs.

  • nan_threshold – Threshold value that will retain any output points containing at least this many non-null input points. The default value is 1.0, which will keep output points containing any non-null inputs. The threshold is applied sequentially to each dimension, and may produce different results than a threshold applied concurrently to all regridding dimensions.

  • output_chunks – Optional dictionary of explicit chunk sizes for the output data. If not provided, the output will be chunked the same as the input data.

Returns:

Regridded input dataset

xarray_regrid.methods.conservative.conservative_regrid_dataset(data: xarray.Dataset, coords: dict[collections.abc.Hashable, xarray.DataArray], latitude_coord: collections.abc.Hashable, skipna: bool, nan_threshold: float, output_chunks: dict[collections.abc.Hashable, int] | None = None) xarray.Dataset[source]

Dataset implementation of the conservative regridding method.

xarray_regrid.methods.conservative.apply_weights(da: xarray.DataArray, weights: dict[collections.abc.Hashable, xarray.DataArray], skipna: bool, nan_threshold: float) xarray.DataArray[source]

Apply the weights over all regridding dimensions simultaneously with xr.dot.

xarray_regrid.methods.conservative.get_valid_threshold(nan_threshold: float) float[source]

Invert the nan_threshold and coerce it to just above zero and below one to handle numerical precision limitations in the weight sum.

xarray_regrid.methods.conservative.get_weights(source_coords: numpy.ndarray, target_coords: numpy.ndarray) numpy.ndarray[source]

Determine the weights to map from the old coordinates to the new coordinates.

Parameters:
  • source_coords – Source coordinates (center points)

  • coordinates (target_coords Target)

Returns:

Weights, which can be used with a dot product to apply the conservative regrid.

xarray_regrid.methods.conservative.apply_spherical_correction(dot_array: xarray.DataArray, latitude_coord: collections.abc.Hashable) xarray.DataArray[source]

Apply a sperical earth correction on the prepared dot product weights.

xarray_regrid.methods.conservative.lat_weight(latitude: numpy.ndarray, latitude_res: float) numpy.ndarray[source]

Return the weight of gridcells based on their latitude.

Parameters:
  • latitude – (Center) latitude values of the gridcells, in degrees.

  • latitude_res – Resolution/width of the grid cells, in degrees.

Returns:

Weights, same shape as latitude input.

xarray_regrid.methods.conservative.format_weights(weights: xarray.DataArray, coord: collections.abc.Hashable, input_dtype: numpy.dtype, input_chunks: tuple[int, Ellipsis] | None, output_chunks: tuple[int, Ellipsis] | int | None) xarray.DataArray[source]

Format the raw weights array such that:

1. Weights match the dtype of the input data 1. Weights are chunked 1:1 with the source data 2. Weights are chunked as requested in the target grid. If no chunks are

provided, the same chunksize as the source grid will be used. See: https://github.com/dask/dask/issues/2225

  1. Weights are converted to a sparse representation (on a per chunk basis)

    if the sparse package is available.