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

Optional YAC backend

UXarray uses its native remapping backend by default. If YAC is installed with its yac.core Python bindings, supported remapping methods can use backend="yac".

YAC can be useful for larger remapping targets or workflows that already rely on YAC interpolation methods.

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)
---------------------------------------------------------------------------
ValueError                                Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/core/dimension.py:1381, in Dimensioned._repr_mimebundle_(self, include, exclude)
   1374 def _repr_mimebundle_(self, include=None, exclude=None):
   1375     """Resolves the class hierarchy for the class rendering the
   1376     object using any display hooks registered on Store.display
   1377     hooks.  The output of all registered display_hooks is then
   1378     combined and returned.
   1379 
   1380     """
-> 1381     return Store.render(self)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/core/options.py:1433, in Store.render(cls, obj)
   1431 data, metadata = {}, {}
   1432 for hook in hooks:
-> 1433     ret = hook(obj)
   1434     if ret is None:
   1435         continue

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/ipython/display_hooks.py:340, in pprint_display(obj)
    338 if not ip.display_formatter.formatters['text/plain'].pprint:
    339     return None
--> 340 return display(obj, raw_output=True)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/ipython/display_hooks.py:311, in display(obj, raw_output, **kwargs)
    309 elif isinstance(obj, (Layout, NdLayout, AdjointLayout)):
    310     with option_state(obj):
--> 311         output = layout_display(obj)
    312 elif isinstance(obj, (HoloMap, DynamicMap)):
    313     with option_state(obj):

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/ipython/display_hooks.py:205, in display_hook.<locals>.wrapped(element)
    203 try:
    204     max_frames = OutputSettings.options['max_frames']
--> 205     mimebundle = fn(element, max_frames=max_frames)
    206     if mimebundle is None:
    207         return {}, {}

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/ipython/display_hooks.py:274, in layout_display(layout, max_frames)
    271     max_frame_warning(max_frames)
    272     return None
--> 274 return render(layout)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/ipython/display_hooks.py:76, in render(obj, **kwargs)
     73 if renderer.fig == 'pdf':
     74     renderer = renderer.instance(fig='png')
---> 76 return renderer.components(obj, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/renderer.py:377, in Renderer.components(self, obj, fmt, comm, **kwargs)
    375     plot = obj
    376 else:
--> 377     plot, fmt = self._validate(obj, fmt)
    379 if not isinstance(plot, Viewable):
    380     html = self._figure_data(plot, fmt, as_script=True, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/renderer.py:308, in Renderer._validate(self, obj, fmt, **kwargs)
    305     plot = HoloViewsPane(obj, center=self.center, backend=self.backend,
    306                          renderer=self)
    307 else:
--> 308     plot = self.get_plot(obj, renderer=self, **kwargs)
    310 all_formats = set(fig_formats + holomap_formats)
    311 if fmt not in all_formats:

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/renderer.py:239, in Renderer.get_plot(self_or_cls, obj, doc, renderer, comm, **kwargs)
    236     defaults = [kd.default for kd in plot.dimensions]
    237     init_key = tuple(v if d is None else d for v, d in
    238                      zip(plot.keys[0], defaults, strict=None))
--> 239     plot.update(init_key)
    240 else:
    241     plot = obj

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/mpl/plot.py:293, in MPLPlot.update(self, key)
    291 def update(self, key):
    292     if len(self) == 1 and key in (0, self.keys[0]) and not self.drawn:
--> 293         return self.initialize_plot()
    294     return self.__getitem__(key)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/mpl/plot.py:62, in mpl_rc_context.<locals>.wrapper(self, *args, **kwargs)
     60 def wrapper(self, *args, **kwargs):
     61     with _rc_context(self.fig_rcparams):
---> 62         return f(self, *args, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/mpl/plot.py:1153, in LayoutPlot.initialize_plot(self)
   1148 traverse_fn = lambda x: x.handles.get('bbox_extra_artists', None)
   1149 extra_artists = list(
   1150     chain.from_iterable(artists for artists in self.traverse(traverse_fn)
   1151                         if artists is not None)
   1152 )
-> 1153 fix_aspect(fig, self.rows, self.cols,
   1154            title_obj, extra_artists,
   1155            vspace=self.vspace*self.fig_scale,
   1156            hspace=self.hspace*self.fig_scale)
   1157 colorbars = self.traverse(specs=[lambda x: hasattr(x, 'colorbar')])
   1158 for cbar_plot in colorbars:

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/holoviews/plotting/mpl/util.py:283, in fix_aspect(fig, nrows, ncols, title, extra_artists, vspace, hspace)
    281 if title and title.get_text():
    282     offset = title.get_window_extent().height/fig.dpi
--> 283 fig.set_size_inches(w, (w*aspect)+offset)
    285 # Redraw and adjust title position if defined
    286 fig.canvas.draw()

File ~/checkouts/readthedocs.org/user_builds/uxarray/conda/v2026.06.0/lib/python3.14/site-packages/matplotlib/figure.py:3167, in Figure.set_size_inches(self, w, h, forward)
   3165 size = np.array([w, h])
   3166 if not np.isfinite(size).all() or (size < 0).any():
-> 3167     raise ValueError(f'figure size must be positive finite not {size}')
   3168 self.bbox_inches.p1 = size
   3169 if forward:

ValueError: figure size must be positive finite not [ 8. nan]
:Layout
   .Overlay.I  :Overlay
      .Image.I     :Image   [Longitude,Latitude]   (Longitude_Latitude bottomDepth)
      .Coastline.I :Feature   [Longitude,Latitude]
   .Overlay.II :Overlay
      .Image.I     :Image   [Longitude,Latitude]   (Longitude_Latitude bottomDepth)
      .Coastline.I :Feature   [Longitude,Latitude]

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

uxds_120.remap

Backend Choice#

UXarray remapping methods use the native UXarray backend by default. If YAC is installed, supported methods can be routed through YAC with backend="yac".

Method

Native backend

YAC backend

nearest_neighbor

Yes

Yes (nnn)

inverse_distance_weighted

Yes

No

bilinear

Yes

Yes (average)

to_rectilinear

Yes

Yes

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.

  • backend
    The remapping backend to use. The default is "uxarray"; set backend="yac" to route the remap through YAC.

  • yac_method
    Required only when backend="yac". Supported values are nnn, average, and conservative; nearest_neighbor() defaults to nnn.

  • yac_options
    Optional dictionary of method-specific keyword arguments forwarded to the selected YAC interpolation routine.

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.

When backend="yac", remap.bilinear() delegates to YAC’s average method. This is the only YAC method exposed through the bilinear convenience accessor; use .remap(..., backend="yac", yac_method=...) if you need to select another YAC method explicitly.

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)

Remapping to a Rectilinear Grid#

The remap.to_rectilinear() accessor remaps unstructured data onto 1D longitude and latitude coordinates and returns a plain xarray.DataArray or xarray.Dataset with regular lat and lon axes.

This is useful when you want to work with UXarray data in tools that expect a structured longitude-latitude grid.

Define the target longitude-latitude grid#

Here, we define a coarse rectilinear grid using 1D longitude and latitude coordinates.

lon = [-180, -120, -60, 0, 60, 120, 180]
lat = [-90, -60, -30, 0, 30, 60, 90]

print(f"Target lon points: {len(lon)}")
print(f"Target lat points: {len(lat)}")

Remap the unstructured field#

Calling to_rectilinear() replaces the unstructured spatial dimension with the new lat and lon dimensions.

rectilinear_bottom_depth = uxds_120["bottomDepth"].remap.to_rectilinear(
    lon=lon,
    lat=lat,
)
rectilinear_bottom_depth

Plot the source and rectilinear output#

First, plot the original unstructured field. Then plot the rectilinear result with xarray directly.

uxds_120["bottomDepth"].plot(title="Bottom Depth (120km)", **common_kwargs)
rectilinear_bottom_depth.plot(
    figsize=(10, 4),
    cmap=cmocean.cm.deep,
    cbar_kwargs={"label": "Bottom Depth"},
)

Alias and YAC backend#

to_structured() and to_lonlat() are aliases for to_rectilinear(). If YAC is installed, pass backend="yac" and choose a yac_method.

Coordinate Handling#

The source data that is being remapped may have spatial coordinate variables present, and they need to be handled properly during the remapping operation. It may include swapping of coordinate values and renaming of some of the coordinates with respect to the dimensions of source data and the remap_to selection. This logic works as follows:

  1. If remap_to matches the source dimension (e.g. source on face centersandremap_to=”faces”` etc.)

    • Swap values of the spatial coordinates with values of the corresponding coordinates from destination_grid.

  2. Else (if remap_to does not match source dimension, e.g. source on face centers but remap_to="nodes" etc.)

    • Swap values of the spatial coordinates with values of the coordinates from destination_grid that are defined on the remap_to dimension.

      • Rename these coordinates to reflect the new element type (e.g. ‘face_x’ → ‘node_x’)

Let us showcase both of these below:

It would be helpful to recall the data array contents again:

Remap with Precomputed Weights#

Use .remap.apply_weights(destination_grid, weight_file) when weights were generated externally with tools such as ESMF or TempestRemap. This path reuses a sparse offline map instead of constructing weights inside UXarray.

remapped = source.remap.apply_weights(destination_grid, "map.nc")

See Remap with Weights for file format details and examples.

uxds_120["bottomDepth"]

Since the source data does not have any coordinate variables, let us add some arbitrarily named lat/lon and x/y coords into it and define it as a new data variable:

uxda_with_coords = ux.core.UxDataArray(
    data=uxds_120["bottomDepth"],
    uxgrid=uxds_120.uxgrid,
    coords={
        "Mesh2_face_lat": uxds_120.uxgrid.face_lat,
        "Mesh_coord_lon": uxds_120.uxgrid.face_lon,
        "Mesh_Faces_x": uxds_120.uxgrid.face_x,
        "Mesh_FACES_y": uxds_120.uxgrid.face_y,
    },
)
uxda_with_coords

Source and remapped output dimensions match#

remapped_same_dims = uxda_with_coords.remap.bilinear(
    destination_grid=uxds_480.uxgrid, remap_to="faces"
)
remapped_same_dims

Note the values of the coordinate variables in the source data have been swapped with the ones from the destination grid, but all the names have been kept as is since this is same-dimension remapping.

Source and remapped output dimensions DO NOT match#

We are now looking into the case of remapping face-centered source data to nodes from the destination grids.

remapped_different_dims = uxda_with_coords.remap.bilinear(
    destination_grid=uxds_480.uxgrid, remap_to="nodes"
)
remapped_different_dims

Note the values of the face-related coordinate variables in the source data have been swapped with the node-related coordinates from the destination grid, and the names of those that had a form of “face” in them have been renamed to have “node” but the one that that did not have any “face” sting has been kept as is.