shxarray.core.xr_accessor#

class shxarray.core.xr_accessor.SHDaAccessor(xarray_obj)#

Bases: ShXrBase

analysis(nmax=100, engine='shlib', **kwargs)#

Apply spherical harmonic analysis from the given ints :param nmax : Spherical harmonic truncation degree of output :return: A datarray with spherical harmonic coefficients derived from the input dataArray

Depending on the method applied different scenarios can be handled input is given on an equidistant longitude, latitude grid (but may be different step size in lon and lat direction)

basinav(dabasins, filtername=None, leakage_corr=None, **kwargs)#

Compute the basin averages of the spherical harmonic coefficients using basin coefficients and possibly a filter and leakage correction :param dabasins: spherical harmonic basin coefficients (basins are uniform masks with 1 inside and 0 outside the region) :type dabasins: xr.DataArray :param filtername: Name of the filter to apply to the coefficients (see sh.filter) :type filtername: str or None :param leakage_corr: Name of the leakage correction to apply to the coefficients. :type leakage_corr: str or None

convolve(kernel)#

Execute a convolution of the dataarray with a given kernel/filter: :param kernel: Kernel object (see e.g shxarray.kernels)

degvar(mean=False)#

Compute the degree variances of spherical harmonic data :param mean: Take the average power per degree instead of the sum (divide by 2n+1) :return: A Dataarray with the degree variances

dvplot(ax=None, mean=False, **kwargs)#

Plot the degree variance of spherical harmonic coefficients

Parameters:
  • ax – Matplotlib axis object to plot on. If None, a new axis will be created

  • **kwargs – Additional arguments passed to the matplotlib function

  • Returns

  • -------- – ax: Matplotlib axis object

filter(filtername, **kwargs)#

Apply well known filters to Spherical harmonic data :param filtername: currently ‘DDKX’ or ‘Gauss’ :param **kwargs:

transpose (default=False): apply he transpose of the filter (only makes sense for anisotropic filters halfwidth (int): specify the halfwidth in km’s for the Gaussian filter, alternatively specify in the format Gauss300 truncate (bool) (default=True): Truncate low degree coefficients which fall outside the filter. Set to False to keep unfiltered input

static from_geoseries(gseries, nmax: int, auxcoord=None, engine='shlib', **kwargs)#

Convert a GeoSeries (from geopandas) to spherical harmonic coefficients

Parameters:
  • gseries (geopandas.GeoSeries) – A GeoSeries Instance of points or polygons

  • nmax (int) – maximum spherical harmonic degree and order to resolve

  • auxcoord (named Pandas.Series or dict(dimname=coordvalues)) – Auxiliary coordinate to map to the dimension of gseries. The default will construct a coordinate with an sequential numerical index and index “id”

  • engine (str, default: 'shlib') – Compute engine to use for the computation. Other options could be ‘shtns’ (when installed). This option has no effect when the input is a GeoSeries of points

  • **kwargs (dict) – Optional arguments which will be passed to either polygon2sh, or point2sh

Returns:

  • xr.DataArray

  • A DataArray holding the spherical harmonic coefficients up to maximum degree specified

static gaunt(n2, n3, m2, m3, engine='shlib')#
static gauntReal(n2, n3, m2, m3, engine='shlib')#
geoid(**kwargs)#
geoplot(ax=None, engine='shlib', add_colorbar=True, **kwargs)#

Plot the spherical harmonic coefficients on a global map using cartopy

Parameters:
  • ax – Cartopy axis object to plot on. If None, a new axis will be created

  • engine (str, default: 'shlib') – Compute engine to use for the synthesis of the spherical harmonic coefficients. Other options could be ‘shtns’ (when installed).

  • add_colorbar (bool, default: True) – Whether to add a colorbar to the plot

  • **kwargs – Additional arguments passed to the pcolormesh plot function

Returns:

ax – The axis on which the spherical harmonic coefficients are plotted

Return type:

cartopy axis object

gravfunctional(outgravtype, ingravtype=None, **kwargs)#

Compute a gravity functional from the spherical harmonic coefficients :param outgravtype: The gravity functional type to compute (see shxarray.kernels.gravfunctionals.gtypes) :param ingravtype: The input gravity functional type to use (default is the one set in the DataArray) :return: A DataArray with the computed gravity functional

horzdef(**kwargs)#
multiply(daother, engine='exp', truncate=False, **kwargs)#
static ones(nmax, nmin=0, name='cnm', auxcoords={}, order='C', nshdims=1)#

1-Initialize an spherical harmonic DataArray based on nmax and nmin

p2s(engine='shlib')#

Returns the product to sum matrix of the spherical harmonic coefficients in the dataarray The output matrix will be symmetric with sides spanning up to nmax/2 of the input

stokes(**kwargs)#
synthesis(lon=None, lat=None, engine='shlib', gtype=None, **kwargs)#

Apply spherical harmonic synthesis on a set of longitude, latitude points gen :param lon: Longitude in degrees East :param lat: Latitude in degrees North :param gtype: Set to false if lon,lat pairs represent individual points :param engine: Spherical harmonic compute engine to use for the computation :return: A datarray for which the spherical harmonic coefficient dimension is mapped to set of points The following scenarios can be handled:

1: lon, lat are Xarray coordinate variables sharing the same dimension mension. Map to a list of points (SH dimension is mapped to a single dimension) 2: lon, lat are Xarray coordinate variables with different dimensions: Map tot a grid spanned by lon,lat 3 lon, lat are list-like objects with the same length: Map to a grid unless (grid=False) 4. lon,lat are list-like objects of different lengths: Map to a grid

synthesis_like(dslonlat=None, engine='shlib', **kwargs)#
to_ascii(out_obj=None)#

Return a string representing the ascii file content of the spherical harmonic coefficients

triplot(ax=None, **kwargs)#

Plot the spherical harmonic coefficients as a triangular plot

Parameters:
  • ax – Matplotlib axis object to plot on. If None, a new axis will be created

  • **kwargs – Additional arguments passed to the pcolormesh plot function

  • Returns

  • -------- – ax: Matplotlib axis object

tws(**kwargs)#
uplift(**kwargs)#
static wigner3j(j2, j3, m2, m3, engine='shlib')#
static zeros(nmax, nmin=0, name='cnm', auxcoords={}, order='C', nshdims=1)#

0-Initialize an spherical harmonic DataArray based on nmax and nmin

class shxarray.core.xr_accessor.SHDsAccessor(xarray_obj)#

Bases: ShXrBase

static ones(nmax, nmin=0, squeeze=True, name='cnm', auxcoords={}, order='C', nshdims=1)#

1-Initialize an spherical harmonic Dataset based on nmax and nmin

synthesis(lon=None, lat=None, engine='shlib', **kwargs)#

Calls the spherical harmonic synthesis operation on all DataArrays which have a ‘nm’ index

static zeros(nmax, nmin=0, squeeze=True, name='cnm', auxcoords={}, order='C', nshdims=1)#

0-Initialize an spherical harmonic Dataset based on nmax and nmin