User API

Regridder

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.

xesmf.util.cell_area(ds, earth_radius=None)

Get cell area of a grid, assuming a sphere.

Parameters:
dsxarray Dataset

Input grid, longitude and latitude required. Curvilinear coordinate system also require cell bounds to be present.

earth_radiusfloat, optional

Earth radius, assuming a sphere, in km.

Returns:
areaxarray DataArray

Cell area. If the earth radius is given, units are km^2, otherwise they are steradian (sr).

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

Equation from [1] [2]:

\[\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.