Remapping#

Remapping, or commonly referred to as Regridding, is the process of taking data that resides on one grid and mapping it to another. Details on various remapping methods can be found here. This user guide section will cover the two native remapping methods that are supported by UXarray:

  • Nearest Neighbor

  • Inverse Distance Weighted

import uxarray as ux
import geoviews.feature as gf
import cartopy.crs as ccrs
import holoviews as hv
import os
import urllib.request
from pathlib import Path
import warnings

warnings.filterwarnings("ignore")

hv.extension("bokeh")

features = gf.coastline(projection=ccrs.PlateCarree(), scale="50m")

Simple Remapping Example#

Provided below is an example of remapping using a simple grid. A has data, B does not. We will remap to B from A.

grid_path = "../../test/meshfiles/ugrid/quad-hexagon/grid.nc"
data_path = "../../test/meshfiles/ugrid/quad-hexagon/data.nc"
destination_grid = "../../test/meshfiles/ugrid/quad-hexagon/triangulated-grid.nc"
grid = ux.open_grid(destination_grid)
uxds = ux.open_dataset(grid_path, data_path)
(
    uxds["t2m"].plot(
        colorbar=False,
        cmap=ux.cmaps.sequential_green_blue,
        title="Mesh with Data",
    )
    * uxds.uxgrid.plot.edges(color="black")
    + grid.plot(colorbar=False, title="Mesh without Data", color="black")
).cols(1)
remapped_grid = uxds["t2m"].remap.inverse_distance_weighted(
    grid, k=3, remap_to="face centers"
)
remapped_grid.plot(
    colorbar=True,
    cmap=ux.cmaps.sequential_green_blue,
    title="Mesh With Remapped Data",
) * remapped_grid.uxgrid.plot.mesh(color="black")

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.

hv.extension("matplotlib")

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)",
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + uxds_120["bottomDepth"].plot(
        title="Bottom Depth (120km)",
        cmap=ux.cmaps.sequential_blue,
    )
    * features
).cols(1).opts(fig_size=125)

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

uxds_120["bottomDepth"].remap
<uxarray.UxDataArray.remap>
Supported Methods:
  * nearest_neighbor(destination_obj, remap_to, coord_type)
  * inverse_distance_weighted(destination_obj, remap_to, coord_type, power, k)

Nearest Neighbor#

Nearest Neighbor Remapping is a point-based method that uses the nearest neighbor from the source grid when assigning data to the destination grid. It is a distance-based approach that uses kd_tree or ball_tree to determine the distance between points. We can use the UxDataArray.remap.nearest_neighbor() method, which takes in the following parameters:

  • destination_grid is the grid object that is being remapped to.

  • destination_obj is being deprecated and soon will no longer be used, it allows remapping to data arrays, grids, and datasets.

  • remap_to specifies the location of the remapping, either to nodes, face centers, or edge centers.

  • coord_type refers to what coordinate system to use, either spherical or cartesian.

Upsampling#

We can remap from the 480km grid to the 120km one, which would perform an upsampling operation. Our destination_obj will be our 120km mesh, and we will make sure to specify to do the remap to face centers. View the plots below to see a comparison of the original 120km mesh compared to the remapped one.

upsampling = uxds_480["bottomDepth"].remap.nearest_neighbor(
    uxds_120.uxgrid, remap_to="face centers"
)
(
    uxds_480["bottomDepth"].plot(
        title="Bottom Depth (480km)", cmap=ux.cmaps.sequential_blue
    )
    * features
    + upsampling.plot(
        title="Remapped Bottom Depth (480km to 120km)",
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + uxds_480["bottomDepth"].plot(
        title="Zoomed (480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + upsampling.plot(
        title="Zoomed Remap (480km to 120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
).cols(1).opts(fig_size=125)

As we can see, the nearest neighbor remap has successfully upsampled, going from a lower resolution to a higher resolution. Of course, since all the data is coming from the lower-resolution grid, it should be noted that no new information is being added, we are simply remapping data from a lower-resolution grid to a higher one. This can be observed in the result since multiple data points from the source grid are repeated in the destination grid.

Downsampling#

Now we can do the reverse, which is downsampling. We we go from a higher resolution (120km) to a lower one (480km). Once again comparision plots are found below, with a zoomed in version provided.

downsampling = uxds_120["bottomDepth"].remap.nearest_neighbor(
    uxds_480.uxgrid, remap_to="face centers"
)
(
    uxds_120["bottomDepth"].plot(
        title="Bottom Depth (120km)", cmap=ux.cmaps.sequential_blue
    )
    * features
    + downsampling.plot(
        title="Remapped Bottom Depth (120km to 480km)",
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + uxds_120["bottomDepth"].plot(
        title="Zoomed (120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + downsampling.plot(
        title="Zoomed Remap (120km to 480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
).cols(1).opts(fig_size=125)

As you can see, the two datasets here don’t look as similar as with the upsampling one. This is an important note, when you are downsampling, there will always be a loss of information. This is because there are many more data points in the source grid than in the destination grid.

Inverse Distance Weighted#

Inverse distance weighted remapping assigns a value based on the weighted average of a specified number of nearby points. This gives a more smooth remapping than the nearest neighbor and helps decrease the chance of outliers. Unlike the nearest neighbor remapping, it constructs new values based on the values around it.

For inverse distance weighted remapping the parameters are the same as nearest neighbor with the addition of two extra parameters we can use to change the remapping parameters.

  • power controls how local or global the remapping is, the larger this value the less influence points that are further away have.

  • k is the number of neighbors to use in the weighted average.

Using the same examples as before, we can see the differences between nearest neighbor and inverse distance weighted remapping.

Upsampling#

We will upsample from the 480km grid to 120km grid again. This time we will see that the results are a lot smoother, and that the 120km looks a lot more like it’s original dataset compared to when the nearest neighbor remap was used.

upsampling_idw = uxds_480["bottomDepth"].remap.inverse_distance_weighted(
    uxds_120.uxgrid, remap_to="face centers"
)
(
    uxds_480["bottomDepth"].plot(
        title="Bottom Depth (480km)", cmap=ux.cmaps.sequential_blue
    )
    * features
    + upsampling_idw.plot(
        title="Remapped Bottom Depth (480km to 120km)",
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + uxds_480["bottomDepth"].plot(
        title="Zoomed (480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + upsampling_idw.plot(
        title="Zoomed Remap (480km to 120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
).cols(1).opts(fig_size=125)

Downsampling#

Here we can see an example of downsampling with inverse distance weighted. It has less of a noticable difference compared to the nearest neighbor, however this may be due to the k number and the default power value used in the calculation. Depending on your graph and the scale difference, you may need to use different k and power values to achieve a smoother result. That will be discussed in the next secition.

downsampling_idw = uxds_120["bottomDepth"].remap.inverse_distance_weighted(
    uxds_480.uxgrid, remap_to="face centers"
)
(
    uxds_120["bottomDepth"].plot(
        title="Bottom Depth (120km)", cmap=ux.cmaps.sequential_blue
    )
    * features
    + downsampling_idw.plot(
        title="Remapped Bottom Depth (120km to 480km)",
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + uxds_120["bottomDepth"].plot(
        title="Zoomed (120km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
    + downsampling_idw.plot(
        title="Zoomed Remap (120km to 480km)",
        xlim=(-10, 10),
        ylim=(-5, 5),
        cmap=ux.cmaps.sequential_blue,
    )
    * features
).cols(1).opts(fig_size=125)