Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
21 changes: 17 additions & 4 deletions doc/recipe/preprocessor.rst
Original file line number Diff line number Diff line change
Expand Up @@ -926,10 +926,23 @@ example of its usage in an ESMValTool preprocessor is:
reference: esmf_regrid.schemes:ESMFAreaWeighted
mdtol: 0.7

.. TODO: Remove the following warning once things have settled a bit.
.. warning::
Just as the mesh support in Iris itself, this new regridding package is
still considered experimental.
Additionally, the use of generic schemes that take source and target grid cubes as
arguments is also supported. The call function for such schemes must be defined as
`(src_cube, grid_cube, **kwargs)` and they must return `iris.cube.Cube` objects.
The `regrid` module will automatically pass the source and grid cubes as inputs
of the scheme. An example of this usage is
the :func:`~esmf_regrid.schemes.regrid_rectilinear_to_rectilinear`
scheme available in :doc:`iris-esmf-regrid:index`:git

.. code-block:: yaml

preprocessors:
regrid_preprocessor:
regrid:
target_grid: 2.5x2.5
scheme:
reference: esmf_regrid.schemes:regrid_rectilinear_to_rectilinear
mdtol: 0.7

.. _ensemble statistics:

Expand Down
1 change: 1 addition & 0 deletions environment.yml
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@ dependencies:
- humanfriendly
- importlib_resources
- iris>=3.2.1
- iris-esmf-regrid
- isodate
- jinja2
- nc-time-axis
Expand Down
69 changes: 39 additions & 30 deletions esmvalcore/preprocessor/_regrid.py
Original file line number Diff line number Diff line change
@@ -1,6 +1,7 @@
"""Horizontal and vertical regridding module."""

import importlib
import inspect
import logging
import os
import re
Expand Down Expand Up @@ -536,8 +537,7 @@ def regrid(cube, target_grid, scheme, lat_offset=True, lon_offset=True):
extrapolation_mode: nanmask

To use the area weighted regridder available in
:class:`esmf_regrid.schemes.ESMFAreaWeighted`, make sure that
:doc:`iris-esmf-regrid:index` is installed and use
:class:`esmf_regrid.schemes.ESMFAreaWeighted` use

.. code-block:: yaml

Expand All @@ -547,34 +547,7 @@ def regrid(cube, target_grid, scheme, lat_offset=True, lon_offset=True):
scheme:
reference: esmf_regrid.schemes:ESMFAreaWeighted

.. note::

Note that :doc:`iris-esmf-regrid:index` is still experimental.
"""
if isinstance(scheme, dict):
try:
object_ref = scheme.pop("reference")
except KeyError as key_err:
raise ValueError(
"No reference specified for generic regridding.") from key_err
module_name, separator, scheme_name = object_ref.partition(":")
try:
obj = importlib.import_module(module_name)
except ImportError as import_err:
raise ValueError(
"Could not import specified generic regridding module. "
"Please double check spelling and that the required module is "
"installed.") from import_err
if separator:
for attr in scheme_name.split('.'):
obj = getattr(obj, attr)
loaded_scheme = obj(**scheme)
else:
loaded_scheme = HORIZONTAL_SCHEMES.get(scheme.lower())
if loaded_scheme is None:
emsg = 'Unknown regridding scheme, got {!r}.'
raise ValueError(emsg.format(scheme))

if isinstance(target_grid, str):
if os.path.isfile(target_grid):
target_grid = iris.load_cube(target_grid)
Expand All @@ -592,14 +565,46 @@ def regrid(cube, target_grid, scheme, lat_offset=True, lon_offset=True):
ycoord = target_grid.coord(axis='y', dim_coords=True)
xcoord.coord_system = src_cs
ycoord.coord_system = src_cs

elif isinstance(target_grid, dict):
# Generate a target grid from the provided specification,
target_grid = _regional_stock_cube(target_grid)

if not isinstance(target_grid, iris.cube.Cube):
raise ValueError(f'Expecting a cube, got {target_grid}.')

if isinstance(scheme, dict):
try:
object_ref = scheme.pop("reference")
except KeyError as key_err:
raise ValueError(
"No reference specified for generic regridding.") from key_err
module_name, separator, scheme_name = object_ref.partition(":")
try:
obj = importlib.import_module(module_name)
except ImportError as import_err:
raise ValueError(
"Could not import specified generic regridding module. "
"Please double check spelling and that the required module is "
"installed.") from import_err
if separator:
for attr in scheme_name.split('.'):
obj = getattr(obj, attr)

scheme_args = inspect.getfullargspec(obj).args
# Add source and target cubes as arguments if required
if 'src_cube' in scheme_args:
scheme['src_cube'] = cube
if 'grid_cube' in scheme_args:
scheme['grid_cube'] = target_grid

loaded_scheme = obj(**scheme)
else:
loaded_scheme = HORIZONTAL_SCHEMES.get(scheme.lower())

if loaded_scheme is None:
emsg = 'Unknown regridding scheme, got {!r}.'
raise ValueError(emsg.format(scheme))

# Unstructured regridding requires x2 2d spatial coordinates,
# so ensure to purge any 1d native spatial dimension coordinates
# for the regridder.
Expand Down Expand Up @@ -629,6 +634,10 @@ def regrid(cube, target_grid, scheme, lat_offset=True, lon_offset=True):
# Perform the horizontal regridding
if _attempt_irregular_regridding(cube, scheme):
cube = esmpy_regrid(cube, target_grid, scheme)
elif isinstance(loaded_scheme, iris.cube.Cube):
# Return regridded cube in cases in which the
# scheme is a function f(src_cube, grid_cube) -> Cube
return loaded_scheme
else:
cube = cube.regrid(target_grid, loaded_scheme)

Expand Down
1 change: 1 addition & 0 deletions setup.py
Original file line number Diff line number Diff line change
Expand Up @@ -32,6 +32,7 @@
'cf-units',
'dask[array]',
'esgf-pyclient>=0.3.1',
'esmf-regrid',
'esmpy!=8.1.0', # see github.com/ESMValGroup/ESMValCore/issues/1208
'fiona',
'fire',
Expand Down
33 changes: 33 additions & 0 deletions tests/integration/preprocessor/_regrid/test_regrid.py
Original file line number Diff line number Diff line change
Expand Up @@ -92,6 +92,18 @@ def test_regrid__linear(self):
expected = np.array([[[1.5]], [[5.5]], [[9.5]]])
self.assert_array_equal(result.data, expected)

def test_regrid__esmf_rectilinear(self):
scheme_name = 'esmf_regrid.schemes:regrid_rectilinear_to_rectilinear'
scheme = {
'reference': scheme_name
}
result = regrid(
self.cube,
self.grid_for_linear,
scheme)
expected = np.array([[[1.5]], [[5.5]], [[9.5]]])
np.testing.assert_array_almost_equal(result.data, expected, decimal=1)

def test_regrid__regular_coordinates(self):
data = np.ones((1, 1))
lons = iris.coords.DimCoord([1.50000000000001],
Expand Down Expand Up @@ -214,6 +226,27 @@ def test_regrid__area_weighted(self):
expected = np.array([1.499886, 5.499886, 9.499886])
np.testing.assert_array_almost_equal(result.data, expected, decimal=6)

def test_regrid__esmf_area_weighted(self):
data = np.empty((1, 1))
lons = iris.coords.DimCoord([1.6],
standard_name='longitude',
bounds=[[1, 2]],
units='degrees_east',
coord_system=self.cs)
lats = iris.coords.DimCoord([1.6],
standard_name='latitude',
bounds=[[1, 2]],
units='degrees_north',
coord_system=self.cs)
coords_spec = [(lats, 0), (lons, 1)]
grid = iris.cube.Cube(data, dim_coords_and_dims=coords_spec)
scheme = {
'reference': 'esmf_regrid.schemes:ESMFAreaWeighted'
}
result = regrid(self.cube, grid, scheme)
expected = np.array([[[1.499886]], [[5.499886]], [[9.499886]]])
np.testing.assert_array_almost_equal(result.data, expected, decimal=6)

def test_regrid__unstructured_nearest_float(self):
"""Test unstructured_nearest regridding with cube of floats."""
result = regrid(self.unstructured_grid_cube,
Expand Down
3 changes: 1 addition & 2 deletions tests/unit/preprocessor/_regrid/test_regrid.py
Original file line number Diff line number Diff line change
Expand Up @@ -103,10 +103,9 @@ def test_invalid_tgt_grid__unknown(self):
regrid(self.src_cube, dummy, scheme)

def test_invalid_scheme__unknown(self):
dummy = mock.sentinel.dummy
emsg = 'Unknown regridding scheme'
with self.assertRaisesRegex(ValueError, emsg):
regrid(dummy, dummy, 'wibble')
regrid(self.src_cube, self.src_cube, 'wibble')

def test_horizontal_schemes(self):
self.assertEqual(
Expand Down