Remapping#

Remapping (also known as regridding) is the process of transferring data defined on one spatial discretization to another. Whether you’re aligning model output to a different grid or comparing datasets on distinct grids, remapping ensures that values are accurately assigned or interpolated between coordinate systems.

For a comprehensive overview of common remapping techniques, see the Climate Data Guide: Regridding Overview.

UXarray currently supports three remapping techniques:

  • Nearest Neighbor

  • Inverse Distance Weighted

  • Bilinear

import os
import urllib.request
import warnings

import cartopy.crs as ccrs
import cmocean
import holoviews as hv

import uxarray as ux

warnings.filterwarnings("ignore")

hv.extension("matplotlib")

common_kwargs = {"cmap": cmocean.cm.deep, "features": ["coastline"]}

Data#

In this notebook, we are using two datasets with different resolutions (480km and 120km) from the MPAS Ocean Model. We will be remapping the bottomDepth variable, which measures the ocean depth.

data_var = "bottomDepth"

grid_filename_480 = "oQU480.grid.nc"
data_filename_480 = "oQU480.data.nc"

grid_filename_120 = "oQU120.grid.nc"
data_filename_120 = "oQU120.data.nc"

filenames = [grid_filename_480, data_filename_480, grid_filename_120, data_filename_120]

for filename in filenames:
    if not os.path.isfile(filename):
        # downloads the files from Cookbook repo, if they haven't been downloaded locally yet
        url = f"https://github.com/ProjectPythia/unstructured-grid-viz-cookbook/raw/main/meshfiles/{filename}"
        _, headers = urllib.request.urlretrieve(url, filename=filename)


file_path_dict = {
    "480km": [grid_filename_480, data_filename_480],
    "120km": [grid_filename_120, data_filename_120],
}
uxds_480 = ux.open_dataset(*file_path_dict["480km"])
uxds_120 = ux.open_dataset(*file_path_dict["120km"])
(
    uxds_480["bottomDepth"].plot(title="Bottom Depth (480km)", **common_kwargs)
    + uxds_120["bottomDepth"].plot(title="Bottom Depth (120km)", **common_kwargs)
).cols(1).opts(fig_size=200)

We can view the supported remapping methods by accessing the .remap accessor that is part of a UxDataArray or UxDataset

uxds_120.remap
<UxDataset.remap>
Supported methods:
  • nearest_neighbor(destination_grid, remap_to='faces')
  • inverse_distance_weighted(destination_grid, remap_to='faces', power=2, k=8)

Nearest Neighbor Remapping#

Nearest-neighbor remapping assigns each point on the destination grid the value of the closest point on the source grid. Under the hood, UXarray leverages a scipy.spatial.KDTree to compute distances efficiently and identify the nearest source location for each destination point.

Use the remap.nearest_neighbor() accessor, which accepts:

  • destination_grid
    The Grid instance you want to interpolate your data onto.

  • remap_to
    The grid element where values should be placed, one of faces, edges, or nodes.

Warning

Nearest-neighbor remapping is fast and simple, but it does not conserve integrated quantities and can introduce discontinuities where grid spacing changes abruptly.

Upsampling#

In this example, we remap data from a coarse 480km resolution grid to a finer 120km grid, which is an example of upsampling.

upsampling = uxds_480["bottomDepth"].remap.nearest_neighbor(
    destination_grid=uxds_120.uxgrid, remap_to="faces"
)
(
    uxds_480["bottomDepth"].plot(title="Bottom Depth (480km)", **common_kwargs)
    + upsampling.plot(title="Remapped Bottom Depth (480km to 120km)", **common_kwargs)
    + uxds_480["bottomDepth"].plot(
        title="Zoomed (480km)", xlim=(-10, 10), ylim=(-5, 5), **common_kwargs
    )
    + upsampling.plot(
        title="Zoomed Remap (480km to 120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(2).opts(fig_size=200)

After applying nearest-neighbor remapping, the data is upsampled from the coarse 480 km grid to the finer 120 km grid. Each destination face simply inherits the value of its closest source face, so no new information is created—existing values are redistributed across the denser mesh. Consequently, you’ll see the same source value repeated on multiple adjacent faces in the 120 km grid.

Downsampling#

In this example, we remap data from a finer 120km resolution grid to a coarser 480km grid, which is an example of downsampling.

downsampling = uxds_120["bottomDepth"].remap.nearest_neighbor(
    destination_grid=uxds_480.uxgrid, remap_to="face centers"
)
(
    uxds_120["bottomDepth"].plot(title="Bottom Depth (120km)", **common_kwargs)
    + downsampling.plot(title="Remapped Bottom Depth (120km to 480km)", **common_kwargs)
    + uxds_120["bottomDepth"].plot(
        title="Zoomed (120km)", xlim=(-10, 10), ylim=(-5, 5), **common_kwargs
    )
    + downsampling.plot(
        title="Zoomed Remap (120km to 480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(2).opts(fig_size=200)

After downsampling, the coarser 480 km grid no longer retains the fine-scale details of the original 120 km dataset. Because the destination grid has far fewer points than the source, information is inevitably lost. Subtle features present at higher resolution will be smoothed out or omitted entirely.

Inverse Distance Weighted Remapping#

Inverse-distance weighted (IDW) remapping computes each destination value as a weighted average of nearby source points, with closer points contributing more strongly. This approach yields a smoother, more continuous field than nearest-neighbor interpolation and can help mitigate isolated outliers.

Use the remap.inverse_distance_weighted() accessor, which accepts the same parameters as nearest-neighbor plus:

  • power
    The exponent governing how rapidly a source point’s influence decays with distance. Larger values localize the interpolation by down-weighting distant points.

  • k
    The number of nearest source points to include in the weighted average.

Note

IDW remapping produces smoother transitions at the cost of additional computation and can blur sharp gradients. It does not conserve integrated quantities by default.

Using the same upsampling and downsampling examples as before, you can compare how IDW preserves continuity relative to nearest-neighbor results.

Upsampling#

Here, we remap data from the coarse 480 km grid onto the finer 120 km grid using IDW interpolation. The resulting field is much smoother than with nearest-neighbor, and the 120 km output more closely resembles the original high-resolution dataset.

upsampling_idw = uxds_480["bottomDepth"].remap.inverse_distance_weighted(
    destination_grid=uxds_120.uxgrid, remap_to="faces"
)
(
    uxds_480["bottomDepth"].plot(title="Bottom Depth (480km)", **common_kwargs)
    + upsampling_idw.plot(
        title="Remapped Bottom Depth (480km to 120km)", **common_kwargs
    )
    + uxds_480["bottomDepth"].plot(
        title="Zoomed (480km)", xlim=(-10, 10), ylim=(-5, 5), **common_kwargs
    )
    + upsampling_idw.plot(
        title="Zoomed Remap (480km to 120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(2).opts(fig_size=200)

Downsampling#

In this example, we use IDW interpolation to downsample data from the 120 km grid to the coarser 480 km grid. Compared to nearest-neighbor remapping, IDW preserves more continuity and mitigates abrupt jumps, though the improvement may be subtle with default parameters. By adjusting:

  • k (number of neighbors)

  • power (distance-decay exponent)

you can control the smoothness and feature retention of the downsampled field. Parameter-tuning strategies will be discussed in the next section.

downsampling_idw = uxds_120["bottomDepth"].remap.inverse_distance_weighted(
    destination_grid=uxds_480.uxgrid, remap_to="faces"
)
(
    uxds_120["bottomDepth"].plot(title="Bottom Depth (120km)", **common_kwargs)
    + downsampling_idw.plot(
        title="Remapped Bottom Depth (120km to 480km)", **common_kwargs
    )
    + uxds_120["bottomDepth"].plot(
        title="Zoomed (120km)", xlim=(-10, 10), ylim=(-5, 5), **common_kwargs
    )
    + downsampling_idw.plot(
        title="Zoomed Remap (120km to 480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(2).opts(fig_size=200)

Comparing k and power Parameters#

  • k
    The number of nearest neighbors included in the interpolation. Larger values draw from more source points, but the impact depends on how those points are weighted.

  • power
    The exponent that governs distance-decay in the weighting function. Higher exponents localize the interpolation by rapidly down-weighting more distant points.

Increasing k alone may have a limited effect if power remains low, and vice versa. To demonstrate their combined influence, we’ll perform two downsampling experiments, one with low k/power values and one with high values—and compare the results side by side.

downsampling_idw_low = uxds_120["bottomDepth"].remap.inverse_distance_weighted(
    uxds_480.uxgrid, remap_to="faces", power=1, k=2
)
downsampling_idw_high = uxds_120["bottomDepth"].remap.inverse_distance_weighted(
    uxds_480.uxgrid, remap_to="faces", power=5, k=128
)
(
    downsampling_idw_low.plot(
        title="Zoomed 480km (power=1, k=2)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
    + downsampling_idw_high.plot(
        title="Zoomed 480km (power=5, k=128)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(1).opts(fig_size=200)
upsampling_idw_low = uxds_480["bottomDepth"].remap.inverse_distance_weighted(
    uxds_120.uxgrid, remap_to="faces", power=1, k=2
)
upsampling_idw_high = uxds_480["bottomDepth"].remap.inverse_distance_weighted(
    uxds_120.uxgrid, remap_to="faces", power=5, k=128
)
(
    upsampling_idw_low.plot(
        title="Zoomed 120km (power=1, k=2)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
    + upsampling_idw_high.plot(
        title="Zoomed 120km (power=5, k=128)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(1).opts(fig_size=200)

When adjusting the k and power parameters during downsampling, the resulting differences are subtle. Downsampling aggregates many fine-scale faces into fewer coarse ones, so whether you include 2 or 128 neighbors, the weighted average over these larger regions changes only marginally.

In contrast, upsampling shows a much stronger response to parameter tuning. Spreading data from a coarse grid (with fewer, larger faces) onto a finer mesh (with more, smaller faces) means that increasing the neighbor count draws in values from a wider spatial context, leading to more pronounced variations in the interpolated field.

Bilinear#

Bilinear remapping breaks down the grid into triangles, and then finds the triangle that contains each point on the destinitation grid. It then uses the values stored at each vertex to bilinearly find a value for the point, depending on it’s postion inside the triangle.

Bilinear remapping finds for each point the polygon that contains it. Then, based on the position of the point within the polygon, it calculates weights. These weights are then applied to the values stored on the polygon. This gives a value to assign to the point. Bilinear remapping tends to give much more smooth results than the previous remapping methods.

Upsampling#

upsampling_bl = uxds_480["bottomDepth"].remap.bilinear(
    destination_grid=uxds_120.uxgrid, remap_to="faces"
)
(
    uxds_480["bottomDepth"].plot(title="Bottom Depth (480km)", **common_kwargs)
    + upsampling_bl.plot(
        title="Remapped Bottom Depth (480km to 120km)", **common_kwargs
    )
    + uxds_480["bottomDepth"].plot(
        title="Zoomed (480km)", xlim=(-10, 10), ylim=(-5, 5), **common_kwargs
    )
    + upsampling_bl.plot(
        title="Zoomed Remap (480km to 120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(2).opts(fig_size=200)

Downsampling#

downsampling_bl = uxds_120["bottomDepth"].remap.bilinear(
    destination_grid=uxds_480.uxgrid, remap_to="faces"
)
(
    uxds_120["bottomDepth"].plot(title="Bottom Depth (120km)", **common_kwargs)
    + downsampling_bl.plot(
        title="Remapped Bottom Depth (120km to 480km)", **common_kwargs
    )
    + uxds_120["bottomDepth"].plot(
        title="Zoomed (120km)", xlim=(-10, 10), ylim=(-5, 5), **common_kwargs
    )
    + downsampling_bl.plot(
        title="Zoomed Remap (120km to 480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        **common_kwargs,
    )
).cols(2).opts(fig_size=200)