API#
Grids#
HEALPix#
From this notebook: https://colab.research.google.com/drive/1MzTyeNFiy-7RNY6UtGKsmDavX5dk6epU
- class earth2grid.healpix.Compass(*values)#
Cardinal directions in counter clockwise order
- class earth2grid.healpix.Grid(level: int, pixel_order: PixelOrder | XY = PixelOrder.RING)#
A Healpix Grid
- Attrs:
level: 2^level = nside pixel_order: the ordering convection of the data
- ang2pix(lon: Tensor, lat: Tensor) Tensor#
Get the pixel index containing the given longitude and latitude
- get_bilinear_regridder_to(lat: ndarray, lon: ndarray)#
Get regridder to the specified lat and lon points
- pix2ang(pix: Tensor, dtype=torch.float64) tuple[Tensor, Tensor]#
Pixel to (lon, lat)
- reorder(order: PixelOrder | XY, x: Tensor) Tensor#
Rorder the pixels of
xto haveorder
- to_image(x: Tensor, fill_value=nan) Tensor#
Use the 45 degree rotated grid pixelation i points to SE, j point to NE
- class earth2grid.healpix.PaddingBackends(*values)#
- class earth2grid.healpix.PixelOrder(*values)#
Healpy has two indexing conventions NEST and RING. But for convolutions we want 2D array indexing in row or column major order. Here are some vectorized routines nest2xy and x2nest for going in between these conventions. The previous code shared by Dale used string computations to handle these operations, which was probably quite slow. Here we use vectorized bit-shifting.
XY orientation#
For array-like indexing can have a different origin and orientation. For example, the default is the origin is S and the data arr[f, y, x] follows the right hand rule. In other words, (x + 1, y) being counterclockwise from (x, y) when looking down on the face.
HPXPad#
For the missing diamonds, N and S of the equatorial faces, the method is not what I expected. Instead of traversing the edge graph as I would expect, it fills in the missing diamond first, and then treats this as if it were an actual tile. For, the NW tile, the rows (along x) are simply shifted along x into the missing face until they hit the diagonal until they intersect the diagonal.
- class earth2grid.healpix.XY(origin: Compass = Compass.S, clockwise: bool = False)#
- Assumes
i = n * n * f + n * y + x
the origin (x,y)=(0,0) is South
if clockwise=False follows the hand rule:
Space | | | / y | / |/______ x
(Thumb points towards Space, index finger towards x, middle finger towards y)
- earth2grid.healpix.angular_to_global(lon: Tensor, lat: Tensor) tuple[Tensor, Tensor]#
Convert angular coordinates (lon, lat) to global HEALPix coordinates (xₛ, yₛ).
- Parameters:
lon – Longitude in degrees East
lat – Latitude in degrees North
- Returns:
tuple of (xs, ys) coordinates
- earth2grid.healpix.conv2d(input, weight, bias=None, stride=1, padding=0, dilation=1, groups=1) Tensor#
Applies a 2D convolution over an input image composed of several input planes.
This operator supports TensorFloat32.
See
Conv2dfor details and output shape.
- earth2grid.healpix.face_to_global(x: Tensor, y: Tensor, face: Tensor) tuple[Tensor, Tensor]#
Convert (x, y, f) to (xₛ, yₛ)
Inverse of
global_to_face- Parameters:
x – x coordinate [0, 1) within face
y – y coordinate [0, 1) within face
face – face index [0, 11]
- Returns:
tuple of (xs, ys) coordinates
- earth2grid.healpix.global_to_angular(xs: Tensor, ys: Tensor) tuple[Tensor, Tensor]#
Inverse HEALPix projection equations.
- Returns:
tuple of (longitude, latitude) in degrees
- earth2grid.healpix.global_to_face(xs: Tensor, ys: Tensor) tuple[Tensor, Tensor, Tensor]#
Convert (xₛ, yₛ) to (x, y, f)
Continuous version of
earth2grid.healpix.double2xy- Parameters:
xs – xs coordinate [0, 2 pi)
ys – ys coordinate [-pi / 2, pi / 2]
- Returns:
tuple of (x, y, f)
- earth2grid.healpix.local2xy(nside: int, x: Tensor, y: Tensor, face: Tensor) Tensor#
Convert a local x, y coordinate in a given face (S origin) to a global pixel index
The local coordinates can be < 0 or > nside
This can be used to implement padding
- Parameters:
nside – int
x – local coordinate [-nside, 2 * nside), origin=S
y – local coordinate [-nside, 2 * nside), origin=S
face – index of the face [0, 12)
right_first – if True then traverse to the face in the x-direction first
- Returns:
- 0 <= x,y< nside. f>12 are the missing faces.
See
_xy_with_filled_tileandpadfor the original hpxpad methods for filling them in.
- Return type:
x, y, f
- earth2grid.healpix.nest2xy(nside, i, pixel_order: ~earth2grid.healpix.core.XY = XY(origin=<Compass.S: 0>, clockwise=False))#
convert NEST to XY index
- earth2grid.healpix.pad(x: Tensor, padding: int) Tensor#
Pad each face consistently with its according neighbors in the HEALPix
- Parameters:
x – The input tensor of shape [N, F, H, W] or [N, F, C, H, W]. Must be ordered in HEALPIX_PAD_XY pixel ordering.
padding – the amount of padding
- Returns:
The tensor padded along the spatial dimensions. For 4D input [N, F, H, W], returns shape [N, F, H+2*padding, W+2*padding]. For 5D input, returns [N, F, C, H+2*padding, W+2*padding],
Examples
Ths example show to pad data described by a
Gridobject.>>> grid = Grid(level=4, pixel_order=PixelOrder.RING) >>> lon = torch.from_numpy(grid.lon) >>> faces = grid.reorder(HEALPIX_PAD_XY, lon) >>> faces = faces.view(1, 12, grid._nside(), grid._nside()) >>> faces.shape torch.Size([1, 12, 16, 16]) >>> padded = pad(faces, padding=1) >>> padded.shape torch.Size([1, 12, 18, 18])
- earth2grid.healpix.pad_backend(backend: PaddingBackends)#
Select the backend for padding
- earth2grid.healpix.pad_with_dim(x, padding, dim=1, pixel_order=XY(origin=<Compass.S: 0>, clockwise=False))#
x[dim] is the spatial dim
- earth2grid.healpix.pcolormesh(z, ax=None, show_axes: bool = False, **kwargs)#
Plot a RING ordered HEALPix map as a pcolormesh.
- Parameters:
z – 1D
ax – An optional matplotlib axis object.
- earth2grid.healpix.reorder(x: Tensor, src_pixel_order: PixelOrder | XY, dest_pixel_order: PixelOrder | XY)#
Reorder x from one pixel order to another
- earth2grid.healpix.ring2double(nside: int, p: ArrayT)#
Compute the (i,j) index in the double pixelization scheme of Calabretta (2007)
This is a visually appealing way to visualize healpix data without any interpolation.
See Fig 5
Calabretta, M. R., & Roukema, B. F. (2007). Mapping on the HEALPix grid. Monthly Notices of the Royal Astronomical Society, 381(2), 865–872. https://doi.org/10.1111/j.1365-2966.2007.12297.x
- earth2grid.healpix.to_double_pixelization(x: ArrayT, fill_value=0) ArrayT#
Convert the array x to 2D-image w/ the double pixelization
xmust be in RING pixel order
- earth2grid.healpix.to_rotated_pixelization(x, fill_value=nan)#
Convert an array to a 2D-iamge w/ the rotated pixelization
- earth2grid.healpix.xy2nest(nside, i, pixel_order: ~earth2grid.healpix.core.XY = XY(origin=<Compass.S: 0>, clockwise=False))#
convert XY index to NEST
- earth2grid.healpix.xy2xy(nside: int, src: XY, dest: XY, i: Tensor)#
Convert flat index between pixel ordering conventions`
- Parameters:
i – int64
- earth2grid.healpix.zonal_average(x: ArrayT, dim=-1) ArrayT#
Compute the zonal average of a map in ring format
Regridding#
- earth2grid.get_regridder(src: Grid, dest: Grid) Module#
Get a regridder from src to dest
- earth2grid.KNNS2Interpolator(src_lon: Tensor, src_lat: Tensor, dest_lon: Tensor, dest_lat: Tensor, k: int = 1, eps=1e-07) Regridder#
K-nearest neighbor interpolator with inverse distance weighting
- Parameters:
src_lon – (m,) source longitude in degrees E
src_lat – (m,) source latitude in degrees N
dest_lon – (n,) output longitude in degrees E
dest_lat – (n,) output latitude in degrees N
k – number of neighbors, default: 1
eps – regularization factor for inverse distance weighting. Only used if k > 1.
- earth2grid.BilinearInterpolator(x_coords: Tensor, y_coords: Tensor, x_query: Tensor, y_query: Tensor, fill_value=nan) None#
Bilinear interpolation for a non-uniform grid
Other utilities#
- earth2grid.healpix.reorder(x: Tensor, src_pixel_order: PixelOrder | XY, dest_pixel_order: PixelOrder | XY)#
Reorder x from one pixel order to another
- earth2grid.healpix.pad(x: Tensor, padding: int) Tensor#
Pad each face consistently with its according neighbors in the HEALPix
- Parameters:
x – The input tensor of shape [N, F, H, W] or [N, F, C, H, W]. Must be ordered in HEALPIX_PAD_XY pixel ordering.
padding – the amount of padding
- Returns:
The tensor padded along the spatial dimensions. For 4D input [N, F, H, W], returns shape [N, F, H+2*padding, W+2*padding]. For 5D input, returns [N, F, C, H+2*padding, W+2*padding],
Examples
Ths example show to pad data described by a
Gridobject.>>> grid = Grid(level=4, pixel_order=PixelOrder.RING) >>> lon = torch.from_numpy(grid.lon) >>> faces = grid.reorder(HEALPIX_PAD_XY, lon) >>> faces = faces.view(1, 12, grid._nside(), grid._nside()) >>> faces.shape torch.Size([1, 12, 16, 16]) >>> padded = pad(faces, padding=1) >>> padded.shape torch.Size([1, 12, 18, 18])