pygmt.sphinterpolate

pygmt.sphinterpolate(data, outgrid=None, **kwargs)[source]

Create spherical grid files in tension of data.

Reads a table containing lon, lat, z columns and performs a Delaunay triangulation to set up a spherical interpolation in tension. Several options may be used to affect the outcome, such as choosing local versus global gradient estimation or optimize the tension selection to satisfy one of four criteria.

Full option list at https://docs.generic-mapping-tools.org/6.5/sphinterpolate.html

Aliases:

  • I = spacing

  • R = region

  • V = verbose

Parameters:
  • data (str, numpy.ndarray, pandas.DataFrame, xarray.Dataset, or geopandas.GeoDataFrame) – Pass in (x, y, z) or (longitude, latitude, elevation) values by providing a file name to an ASCII data table, a 2-D numpy.ndarray, a pandas.DataFrame, an xarray.Dataset made up of 1-D xarray.DataArray data variables, or a geopandas.GeoDataFrame containing the tabular data.

  • outgrid (str | None, default: None) – Name of the output netCDF grid file. If not specified, will return an xarray.DataArray object. For writing a specific grid file format or applying basic data operations to the output grid, see https://docs.generic-mapping-tools.org/6.5/gmt.html#grd-inout-full for the available modifiers.

  • spacing (float, str, or list) –

    x_inc[+e|n][/y_inc[+e|n]]. x_inc [and optionally y_inc] is the grid spacing.

    • Geographical (degrees) coordinates: Optionally, append an increment unit. Choose among m to indicate arc-minutes or s to indicate arc-seconds. If one of the units e, f, k, M, n or u is appended instead, the increment is assumed to be given in meter, foot, km, mile, nautical mile or US survey foot, respectively, and will be converted to the equivalent degrees longitude at the middle latitude of the region (the conversion depends on PROJ_ELLIPSOID). If y_inc is given but set to 0 it will be reset equal to x_inc; otherwise it will be converted to degrees latitude.

    • All coordinates: If +e is appended then the corresponding max x (east) or y (north) may be slightly adjusted to fit exactly the given increment [by default the increment may be adjusted slightly to fit the given domain]. Finally, instead of giving an increment you may specify the number of nodes desired by appending +n to the supplied integer argument; the increment is then recalculated from the number of nodes, the registration, and the domain. The resulting increment value depends on whether you have selected a gridline-registered or pixel-registered grid; see GMT File Formats for details.

    Note: If region=grdfile is used then the grid spacing and the registration have already been initialized; use spacing and registration to override these values.

  • region (str or list) – xmin/xmax/ymin/ymax[+r][+uunit]. Specify the region of interest.

  • verbose (bool or str) –

    Select verbosity level [Default is w], which modulates the messages written to stderr. Choose among 7 levels of verbosity:

    • q - Quiet, not even fatal error messages are produced

    • e - Error messages only

    • w - Warnings [Default]

    • t - Timings (report runtimes for time-intensive algorithms)

    • i - Informational messages (same as verbose=True)

    • c - Compatibility warnings

    • d - Debugging messages

Returns:

ret (xarray.DataArray or None) – Return type depends on whether the outgrid parameter is set:

  • xarray.DataArray if outgrid is not set

  • None if outgrid is set (grid output will be stored in file set by outgrid)

Example

>>> import pygmt
>>> # Load a table of Mars with longitude/latitude/radius columns
>>> mars_shape = pygmt.datasets.load_sample_data(name="mars_shape")
>>> # Perform Delaunay triangulation on the table data
>>> # to produce a grid with a 1 arc-degree spacing
>>> grid = pygmt.sphinterpolate(data=mars_shape, spacing=1, region="g")