spatial

Utilities for reading, writing and operating on spatial data

  • Can read ascii, netCDF4, HDF5 or geotiff files

  • Can output to ascii, netCDF4, HDF5 or geotiff files

Calling Sequence

Reading a netCDF4 file

import pyTMD.spatial
dinput = pyTMD.spatial.from_netCDF4(path_to_netCDF4_file)

Reading a HDF5 file

import pyTMD.spatial
dinput = pyTMD.spatial.from_HDF5(path_to_HDF5_file)

Source code

General Methods

pyTMD.spatial.case_insensitive_filename(filename: str | pathlib.Path)[source]

Searches a directory for a filename without case dependence

Parameters
filename: str

input filename

pyTMD.spatial.data_type(x: ndarray, y: ndarray, t: ndarray) str[source]

Determines input data type based on variable dimensions

Parameters
x: np.ndarray

x-dimension coordinates

y: np.ndarray

y-dimension coordinates

t: np.ndarray

time-dimension coordinates

Returns
string denoting input data type
  • 'time series'

  • 'drift'

  • 'grid'

pyTMD.spatial.from_file(filename: str, format: str, **kwargs)[source]

Wrapper function for reading data from an input format

Parameters
filename: str

full path of input file

format: str

format of input file

**kwargs: dict

Keyword arguments for file reader

pyTMD.spatial.from_ascii(filename: str, **kwargs)[source]

Read data from an ascii file

Parameters
filename: str

full path of input ascii file

compression: str or NoneType, default None

file compression type

columns: list, default [‘time’, ‘y’, ‘x’, ‘data’]

column names of ascii file

delimiter: str, default ‘,’

Delimiter for csv or ascii files

header: int, default 0

header lines to skip from start of file

parse_dates: bool, default False

Try parsing the time column

pyTMD.spatial.from_netCDF4(filename: str, **kwargs)[source]

Read data from a netCDF4 file

Parameters
filename: str

full path of input netCDF4 file

compression: str or NoneType, default None

file compression type

group: str or NoneType, default None

netCDF4 variable group

timename: str, default ‘time’

name for time-dimension variable

xname: str, default ‘lon’

name for x-dimension variable

yname: str, default ‘lat’

name for y-dimension variable

varname: str, default ‘data’

name for data variable

field_mapping: dict, default {}

mapping between output variables and input netCDF4

pyTMD.spatial.from_HDF5(filename: str | pathlib.Path, **kwargs)[source]

Read data from a HDF5 file

Parameters
filename: str

full path of input HDF5 file

compression: str or NoneType, default None

file compression type

group: str or NoneType, default None

netCDF4 variable group

timename: str, default ‘time’

name for time-dimension variable

xname: str, default ‘lon’

name for x-dimension variable

yname: str, default ‘lat’

name for y-dimension variable

varname: str, default ‘data’

name for data variable

field_mapping: dict, default {}

mapping between output variables and input HDF5

pyTMD.spatial.from_geotiff(filename: str, **kwargs)[source]

Read data from a geotiff file

Parameters
filename: str

full path of input geotiff file

compression: str or NoneType, default None

file compression type

bounds: list or NoneType, default bounds

extent of the file to read: [xmin, xmax, ymin, ymax]

pyTMD.spatial.to_ascii(output: dict, attributes: dict, filename: str, **kwargs)[source]

Write data to an ascii file

Parameters
output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

filename: str

full path of output ascii file

delimiter: str, default ‘,’

delimiter for output spatial file

columns: list, default [‘time’, ‘y’, ‘x’, ‘data’]

column names of ascii file

header: bool, default False

create a YAML header with data attributes

pyTMD.spatial.to_netCDF4(output: dict, attributes: dict, filename: str | pathlib.Path, **kwargs)[source]

Wrapper function for writing data to a netCDF4 file

Parameters
output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

filename: str

full path of output netCDF4 file

data_type: str, default ‘drift’

Input data type

  • 'time series'

  • 'drift'

  • 'grid'

pyTMD.spatial._drift_netCDF4(fileID, output: dict, attributes: dict, **kwargs)[source]

Write drift data variables to a netCDF4 file object

Parameters
fileID: obj

open netCDF4 file object

output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

pyTMD.spatial._grid_netCDF4(fileID, output: dict, attributes: dict, **kwargs)[source]

Write gridded data variables to a netCDF4 file object

Parameters
fileID: obj

open netCDF4 file object

output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

pyTMD.spatial._time_series_netCDF4(fileID, output: dict, attributes: dict, **kwargs)[source]

Write time series data variables to a netCDF4 file object

Parameters
fileID: obj

open netCDF4 file object

output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

pyTMD.spatial.to_HDF5(output: dict, attributes: dict, filename: str, **kwargs)[source]

Write data to a HDF5 file

Parameters
output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

filename: str

full path of output HDF5 file

pyTMD.spatial.to_geotiff(output: dict, attributes: dict, filename: str, **kwargs)[source]

Write data to a geotiff file

Parameters
output: dict

python dictionary of output data

attributes: dict

python dictionary of output attributes

filename: str

full path of output geotiff file

varname: str, default ‘data’

output variable name

driver: str, default ‘cog’

GDAL driver

  • 'GTiff': GeoTIFF

  • 'cog': Cloud Optimized GeoTIFF

dtype: obj, default osgeo.gdal.GDT_Float64

GDAL data type

options: list, default [‘COMPRESS=LZW’]

GDAL driver creation options

pyTMD.spatial.expand_dims(obj: dict, varname: str = 'data')[source]

Add a singleton dimension to a spatial dictionary if non-existent

Parameters
obj: dict

python dictionary of data

varname: str, default data

variable name to expand

pyTMD.spatial.default_field_mapping(variables: list | numpy.ndarray)[source]

Builds field mappings from a variable list

Parameters
variables: list

Variables from argument parser

  • time

  • yname

  • xname

  • varname

Returns
field_mapping: dict

Field mappings for netCDF4/HDF5 read functions

pyTMD.spatial.inverse_mapping(field_mapping)[source]

Reverses the field mappings of a dictionary

Parameters
field_mapping: dict

Field mappings for netCDF4/HDF5 read functions

pyTMD.spatial.convert_ellipsoid(phi1: ndarray, h1: ndarray, a1: float, f1: float, a2: float, f2: float, eps: float = 1e-12, itmax: int = 10)[source]

Convert latitudes and heights to a different ellipsoid using Newton-Raphson

Parameters
phi1: np.ndarray

latitude of input ellipsoid in degrees

h1: np.ndarray

height above input ellipsoid in meters

a1: float

semi-major axis of input ellipsoid

f1: float

flattening of input ellipsoid

a2: float

semi-major axis of output ellipsoid

f2: float

flattening of output ellipsoid

eps: float, default 1e-12

tolerance to prevent division by small numbers and to determine convergence

itmax: int, default 10

maximum number of iterations to use in Newton-Raphson

Returns
phi2: np.ndarray

latitude of output ellipsoid in degrees

h2: np.ndarray

height above output ellipsoid in meters

References

1
  1. Meeus, Astronomical Algorithms, 2nd edition, 477 pp., (1998).

pyTMD.spatial.compute_delta_h(a1: float, f1: float, a2: float, f2: float, lat: ndarray)[source]
Compute difference in elevation for two ellipsoids at a given

latitude using a simplified empirical equation

Parameters
a1: float

semi-major axis of input ellipsoid

f1: float

flattening of input ellipsoid

a2: float

semi-major axis of output ellipsoid

f2: float

flattening of output ellipsoid

lat: np.ndarray

latitudes (degrees north)

Returns
delta_h: np.ndarray

difference in elevation for two ellipsoids

References

1

J Meeus, Astronomical Algorithms, pp. 77–82, (1991).

pyTMD.spatial.wrap_longitudes(lon: float | numpy.ndarray)[source]

Wraps longitudes to range from -180 to +180

Parameters
lon: float or np.ndarray

longitude (degrees east)

pyTMD.spatial.to_cartesian(lon: ndarray, lat: ndarray, h: float | numpy.ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Converts geodetic coordinates to Cartesian coordinates

Parameters
lon: np.ndarray

longitude (degrees east)

lat: np.ndarray

latitude (degrees north)

h: float or np.ndarray, default 0.0

height above ellipsoid (or sphere)

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

for spherical coordinates set to radius of the Earth

flat: float, default 1.0/298.257223563

ellipsoidal flattening

for spherical coordinates set to 0

pyTMD.spatial.to_sphere(x: ndarray, y: ndarray, z: ndarray)[source]

Convert from cartesian coordinates to spherical coordinates

Parameters
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

pyTMD.spatial.to_geodetic(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, method: str = 'bowring', eps: float = 2.220446049250313e-16, iterations: int = 10)[source]

Convert from cartesian coordinates to geodetic coordinates using either iterative or closed-form methods

Parameters
x, float

cartesian x-coordinates

y, float

cartesian y-coordinates

z, float

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

method: str, default ‘bowring’

method to use for conversion

  • 'moritz': iterative solution

  • 'bowring': iterative solution

  • 'zhu': closed-form solution

eps: float, default np.finfo(np.float64).eps

tolerance for iterative methods

iterations: int, default 10

maximum number of iterations

pyTMD.spatial._moritz_iterative(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, eps: float = 2.220446049250313e-16, iterations: int = 10)[source]

Convert from cartesian coordinates to geodetic coordinates using the iterative solution of [1]

Parameters
x, float

cartesian x-coordinates

y, float

cartesian y-coordinates

z, float

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

eps: float, default np.finfo(np.float64).eps

tolerance for iterative method

iterations: int, default 10

maximum number of iterations

References

1

B. Hofmann-Wellenhof and H. Moritz, Physical Geodesy, 2nd Edition, 403 pp., (2006). doi: 10.1007/978-3-211-33545-1

pyTMD.spatial._bowring_iterative(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, eps: float = 2.220446049250313e-16, iterations: int = 10)[source]

Convert from cartesian coordinates to geodetic coordinates using the iterative solution of [1] [2]

Parameters
x, float

cartesian x-coordinates

y, float

cartesian y-coordinates

z, float

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

eps: float, default np.finfo(np.float64).eps

tolerance for iterative method

iterations: int, default 10

maximum number of iterations

References

1

B. R. Bowring, “Transformation from spatial to geodetic coordinates,” Survey Review, 23(181), 323–327, (1976). doi: 10.1179/sre.1976.23.181.323

2

B. R. Bowring, “The Accuracy Of Geodetic Latitude and Height Equations,” Survey Review, 28(218), 202–206, (1985). doi: 10.1179/sre.1985.28.218.202

pyTMD.spatial._zhu_closed_form(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Convert from cartesian coordinates to geodetic coordinates using the closed-form solution of [1]

Parameters
x, float

cartesian x-coordinates

y, float

cartesian y-coordinates

z, float

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

References

1

J. Zhu, “Exact conversion of Earth-centered, Earth-fixed coordinates to geodetic coordinates,” Journal of Guidance, Control, and Dynamics, 16(2), 389–391, (1993). doi: 10.2514/3.21016

pyTMD.spatial.scale_factors(lat: ndarray, flat: float = 0.0033528106647474805, reference_latitude: float = 70.0, metric: str = 'area')[source]

Calculates scaling factors to account for polar stereographic distortion including special case of at the exact pole [1] [2]

Parameters
lat: np.ndarray

latitude (degrees north)

flat: float, default 1.0/298.257223563

ellipsoidal flattening

reference_latitude: float, default 70.0

reference latitude (true scale latitude)

metric: str, default ‘area’

metric to calculate scaling factors

  • 'distance': scale factors for distance

  • 'area': scale factors for area

Returns
scale: np.ndarray

scaling factors at input latitudes

References

1

J. P. Snyder, Map Projections used by the U.S. Geological Survey, Geological Survey Bulletin 1532, U.S. Government Printing Office, (1982).

2

JPL Technical Memorandum 3349-85-101