User API
Regridder
- class xesmf.frontend.Regridder(ds_in, ds_out, method, locstream_in=False, locstream_out=False, periodic=False, parallel=False, **kwargs)
- __init__(ds_in, ds_out, method, locstream_in=False, locstream_out=False, periodic=False, parallel=False, **kwargs)
Make xESMF regridder
- Parameters:
- ds_in, ds_outxarray Dataset, DataArray, or dictionary
Contain input and output grid coordinates. All variables that the cf-xarray accessor understand are accepted. Otherwise, look for
lon
,lat
, optionallylon_b
,lat_b
for conservative methods, andmask
. Note that for mask, the ESMF convention is used, where masked values are identified by 0, and non-masked values by 1.For conservative methods, if bounds are not present, they will be computed using cf-xarray (only 1D coordinates are currently supported).
Shape can be 1D (n_lon,) and (n_lat,) for rectilinear grids, or 2D (n_y, n_x) for general curvilinear grids. Shape of bounds should be (n+1,) or (n_y+1, n_x+1). CF-bounds (shape (n, 2) or (n, m, 4)) are also accepted if they are accessible through the cf-xarray accessor.
If either dataset includes a 2d mask variable, that will also be used to inform the regridding.
If DataArrays are passed, the are simply converted to Datasets.
- methodstr
Regridding method. Options are
‘bilinear’
‘conservative’, need grid corner information
‘conservative_normed’, need grid corner information
‘patch’
‘nearest_s2d’
‘nearest_d2s’
- periodicbool, optional
Periodic in longitude? Default to False. Only useful for global grids with non-conservative regridding. Will be forced to False for conservative regridding.
- parallelbool, optional
Compute the weights in parallel with Dask. Default to False. If True, weights are computed in parallel with Dask on subsets of the output grid using chunks of the output grid.
- filenamestr, optional
Name for the weight file. The default naming scheme is:
{method}_{Ny_in}x{Nx_in}_{Ny_out}x{Nx_out}.nc
e.g. bilinear_400x600_300x400.nc
- reuse_weightsbool, optional
Whether to read existing weight file to save computing time. False by default (i.e. re-compute, not reuse).
- extrap_methodstr, optional
Extrapolation method. Options are
‘inverse_dist’
‘nearest_s2d’
- extrap_dist_exponentfloat, optional
The exponent to raise the distance to when calculating weights for the extrapolation method. If none are specified, defaults to 2.0
- extrap_num_src_pntsint, optional
The number of source points to use for the extrapolation methods that use more than one source point. If none are specified, defaults to 8
- weightsNone, coo_matrix, dict, str, Dataset, Path,
- Regridding weights, stored as
a scipy.sparse COO matrix,
a dictionary with keys row_dst, col_src and weights,
an xarray Dataset with data variables col, row and S,
or a path to a netCDF file created by ESMF.
If None, compute the weights.
- ignore_degeneratebool, optional
If False (default), raise error if grids contain degenerated cells (i.e. triangles or lines, instead of quadrilaterals)
- unmapped_to_nan: boolean, optional
Set values of unmapped points to np.nan instead of zero (ESMF default). This is useful for target cells lying outside of the source domain when no output mask is defined. If an output mask is defined, or regridding method is nearest_s2d or nearest_d2s, this option has no effect.
- Returns:
- regridderxESMF regridder object
- __call__(indata, keep_attrs=False, skipna=False, na_thres=1.0, output_chunks=None)
Apply regridding to input data.
- Parameters:
- indatanumpy array, dask array, xarray DataArray or Dataset.
If not an xarray object or if input_dìms was not given in the init, the rightmost two dimensions must be the same as
ds_in
. Can have arbitrary additional dimensions.Examples of valid shapes
(n_lat, n_lon), if
ds_in
has shape (n_lat, n_lon)(n_time, n_lev, n_y, n_x), if
ds_in
has shape (Ny, n_x)
Either give input_dims or transpose your input data if the horizontal dimensions are not the rightmost two dimensions
Variables without the regridded dimensions are silently skipped when passing a Dataset.
- keep_attrsbool, optional
Keep attributes for xarray DataArrays or Datasets. Defaults to False.
- skipna: bool, optional
Whether to skip missing values when regridding. When set to False, an output value is masked when a single input value is missing and no grid mask is provided. When set to True, missing values do not contaminate the regridding since only valid values are taken into account. In this case, a given output point is set to NaN only if the ratio of missing values exceeds the level set by na_thres: for instance, when the center of a cell is computed linearly from its four corners, one of which is missing, the output value is set to NaN if na_thres is smaller than 0.25.
- na_thres: float, optional
A value within the [0., 1.] interval that defines the maximum ratio of missing grid points involved in the regrdding over which the output value is set to NaN. For instance, if na_thres is set to 0, the output value is NaN if a single NaN is found in the input values that are used to compute the output value; similarly, if na_thres is set to 1, all input values must be missing to mask the output value.
- output_chunks: dict or tuple, optional
If indata is a dask_array_type, the desired chunks to have on the output data along the spatial axes. Other non-spatial axes inherit the same chunks as indata as those are not affected by the application of the weights. Default behavior is to have the outdata chunks be like the indata chunks. Chunks have to be specified for all spatial dimensions of the output data otherwise regridding will fail. output_chunks can either be a tuple the same size as the spatial axes of outdata or it can be a dict with defined dims. If output_chunks is a dict, the keys must match the dims of the output grid passed when initializing this Regridder.
- Returns:
- outdataData type is the same as input data type, except for datasets.
On the same horizontal grid as
ds_out
, with extra dims indr_in
.Assuming
ds_out
has the shape of (n_y_out, n_x_out), examples of returning shapes are(n_y_out, n_x_out), if
dr_in
is 2D(n_time, n_lev, n_y_out, n_x_out), if
dr_in
has shape (n_time, n_lev, n_y, n_x)
Datasets with dask-backed variables will have modified dtypes. If all input variables are ‘float32’, all output will be ‘float32’, for any other case, all outputs will be ‘float64’.
- regrid_array(indata, weights, skipna=False, na_thres=1.0, output_chunks=None)
See __call__().
- regrid_dataarray(dr_in, keep_attrs=False, skipna=False, na_thres=1.0, output_chunks=None)
See __call__().
- regrid_dataset(ds_in, keep_attrs=False, skipna=False, na_thres=1.0, output_chunks=None)
See __call__().
- to_netcdf(filename=None)
Save weights to disk as a netCDF file.
- property w: xarray.DataArray
Return weights as a 4D DataArray with dimensions (y_out, x_out, y_in, x_in).
ESMF stores the weights in a 2D array with dimensions (out_dim, in_dim), the size of the output and input grids respectively (ny x nx). This property returns the weights reshaped as a 4D array to simplify comparisons with the original grids.
- class xesmf.frontend.SpatialAverager(ds_in, polys, ignore_holes=False, periodic=False, filename=None, reuse_weights=False, weights=None, ignore_degenerate=False, geom_dim_name='geom')
- __init__(ds_in, polys, ignore_holes=False, periodic=False, filename=None, reuse_weights=False, weights=None, ignore_degenerate=False, geom_dim_name='geom')
Compute the exact average of a gridded array over a geometry.
This uses the ESMF conservative regridding method to compute and apply weights mapping a 2D field unto geometries defined by polygons. The conservative method preserves the areal average of the input field. That is, the value at each output grid cell is the average input value over the output grid area. Here, the output grid cells are not rectangles defined by four corners, but polygons defined by multiple vertices (ESMF.Mesh objects). The regridding weights thus compute the areal-average of the input grid over each polygon.
For multi-parts geometries (shapely.MultiPolygon), weights are computed for each geometry, then added, to compute the average over all geometries.
When polygons include holes, the weights over the holes can either be substracted, or ignored.
- Parameters:
- ds_inxr.DataArray or xr.Dataset or dictionary
Contain input and output grid coordinates. Look for variables
lon
,lat
,lon_b
andlat_b
.Optionally looks for
mask
, in which case the ESMF convention is used, where masked values are identified by 0, and non-masked values by 1.Shape can be 1D (n_lon,) and (n_lat,) for rectilinear grids, or 2D (n_y, n_x) for general curvilinear grids. Shape of bounds should be (n+1,) or (n_y+1, n_x+1). DataArrays are converted to Datasets.
- polyssequence of shapely Polygons and MultiPolygons
Sequence of polygons (lon, lat) over which to average ds_in.
- ignore_holesbool
Whether to ignore holes in polygons. Default (True) is to subtract the weight of holes from the weight of the polygon.
- filenamestr, optional
Name for the weight file. The default naming scheme is:
spatialavg_{Ny_in}x{Nx_in}_{Npoly_out}.nc
e.g. spatialavg_400x600_30.nc
- reuse_weightsbool, optional
Whether to read existing weight file to save computing time. False by default (i.e. re-compute, not reuse).
- weightsNone, coo_matrix, dict, str, Dataset, Path,
- Regridding weights, stored as
a scipy.sparse COO matrix,
a dictionary with keys row_dst, col_src and weights,
an xarray Dataset with data variables col, row and S,
or a path to a netCDF file created by ESMF.
If None, compute the weights.
- ignore_degeneratebool, optional
If False (default), raise error if grids contain degenerated cells (i.e. triangles or lines, instead of quadrilaterals)
- self.geom_dim_namestr
Name of dimension along which averages for each polygon are stored.
- Returns:
- xarray.DataArray
Average over polygons along geom_dim_name dimension. The lon and lat coordinates are the polygon centroid coordinates.
References
This approach is inspired by OCGIS.
- property w: xarray.DataArray
Return weights as a 3D DataArray with dimensions (geom, y_in, x_in).
ESMF stores the weights in a 2D array with dimensions (out_dim, in_dim), the size of the output and input grids respectively (ny x nx). This property returns the weights reshaped as a 3D array to simplify comparisons with the original grids.
util
- xesmf.util.grid_2d(lon0_b, lon1_b, d_lon, lat0_b, lat1_b, d_lat)
2D rectilinear grid centers and bounds
- Parameters:
- lon0_b, lon1_bfloat
Longitude bounds
- d_lonfloat
Longitude step size, i.e. grid resolution
- lat0_b, lat1_bfloat
Latitude bounds
- d_latfloat
Latitude step size, i.e. grid resolution
- Returns:
- dsxarray DataSet with coordinate values
- xesmf.util.cf_grid_2d(lon0_b, lon1_b, d_lon, lat0_b, lat1_b, d_lat)
CF compliant 2D rectilinear grid centers and bounds.
- Parameters:
- lon0_b, lon1_bfloat
Longitude bounds
- d_lonfloat
Longitude step size, i.e. grid resolution
- lat0_b, lat1_bfloat
Latitude bounds
- d_latfloat
Latitude step size, i.e. grid resolution
- Returns:
- dsxarray.DataSet with coordinate values
- xesmf.util.grid_global(d_lon, d_lat, cf=False, lon1=180)
Global 2D rectilinear grid centers and bounds
- Parameters:
- d_lonfloat
Longitude step size, i.e. grid resolution
- d_latfloat
Latitude step size, i.e. grid resolution
- cfbool
Return a CF compliant grid.
- lon1{180, 360}
Right longitude bound. According to which convention is used longitudes will vary from -180 to 180 or from 0 to 360.
- Returns:
- dsxarray DataSet with coordinate values
- xesmf.util.split_polygons_and_holes(polys)
Split the exterior boundaries and the holes for a list of polygons.
If MultiPolygons are encountered in the list, they are flattened out in their constituents.
- Parameters:
- polysSequence of shapely Polygons or MultiPolygons
- Returns:
- exteriorslist of Polygons
The polygons without any holes
- holeslist of Polygons
Holes of the polygons as polygons
- i_extlist of integers
The index in polys of each polygon in exteriors.
- i_hollist of integers
The index in polys of the owner of each hole in holes.
- xesmf.util.simple_tripolar_grid(nlons, nlats, lat_cap=60, lon_cut=-300)
Generate a simple tripolar grid, regular under lat_cap.
- Parameters:
- nlons: int
Number of longitude points.
- nlats: int
Number of latitude points.
- lat_cap: float
Latitude of the northern cap.
- lon_cut: float
Longitude of the periodic boundary.
data
Standard test data for regridding benchmark.
- xesmf.data.wave_smooth(lon, lat)
Spherical harmonic with low frequency.
- Parameters:
- lon, lat2D numpy array or xarray DataArray
Longitute/Latitude of cell centers
- Returns:
- f2D numpy array or xarray DataArray depending on input
2D wave field
Notes
\[\begin{split}Y_2^2 = 2 + \cos^2(\\theta) \cos(2 \phi)\end{split}\]References
[1]Jones, P. W. (1999). First-and second-order conservative remapping schemes for grids in spherical coordinates. Monthly Weather Review, 127(9), 2204-2210.
[2]Ullrich, P. A., Lauritzen, P. H., & Jablonowski, C. (2009). Geometrically exact conservative remapping (GECoRe): regular latitude–longitude and cubed-sphere grids. Monthly Weather Review, 137(6), 1721-1741.