spatial
Spatial transformation routines
General Methods
- 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 coordinates
- Returns:
- String denoting input data type
'time series''trajectory''grid'
- pyTMD.spatial.convert_ellipsoid(lat1: 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 [45]
- Parameters:
- lat1: np.ndarray
Latitude of input ellipsoid (degrees)
- h1: np.ndarray
Height above input ellipsoid (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:
- lat2: np.ndarray
Latitude of output ellipsoid (degrees)
- h2: np.ndarray
Height above output ellipsoid (meters)
- pyTMD.spatial.compute_delta_h(lat: ndarray, a1: float, f1: float, a2: float, f2: float)[source]
Compute difference in elevation for two ellipsoids at a given latitude using a simplified empirical relation [45]
- Parameters:
- lat: np.ndarray
Latitudes (degrees north)
- 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
- Returns:
- delta_h: np.ndarray
Difference in elevation for two ellipsoids
- pyTMD.spatial.wrap_longitudes(lon: float | ndarray)[source]
Wraps longitudes to range from -180 to +180
- Parameters:
- lon: float or np.ndarray
Longitude (degrees east)
- pyTMD.spatial.to_dms(d: ndarray)[source]
Convert decimal degrees to degrees, minutes and seconds
- Parameters:
- d: np.ndarray
Angle (decimal degrees)
- Returns:
- degree: np.ndarray
Degrees
- minute: np.ndarray
Minutes (arcminutes)
- second: np.ndarray
Seconds (arcseconds)
- pyTMD.spatial.from_dms(degree: ndarray, minute: ndarray, second: ndarray)[source]
Convert degrees, minutes and seconds to decimal degrees
- Parameters:
- degree: np.ndarray
Degrees
- minute: np.ndarray
Minutes (arcminutes)
- second: np.ndarray
Seconds (arcseconds)
- Returns:
- d: np.ndarray
Angle (decimal degrees)
- pyTMD.spatial.to_cartesian(lon: ndarray, lat: ndarray, h: float | 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
Semi-major 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
- Returns:
- x: np.ndarray
Cartesian x-coordinates (meters)
- y: np.ndarray
Cartesian y-coordinates (meters)
- z: np.ndarray
Cartesian z-coordinates (meters)
- 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 (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- Returns:
- lon: np.ndarray
Longitude (degrees east)
- lat: np.ndarray
Latitude (degrees north)
- rad: np.ndarray
Radius (meters)
- pyTMD.spatial.to_geodetic(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, method: str = 'bowring', eps: float = np.float64(2.220446049250313e-16), iterations: int = 10)[source]
Convert from cartesian coordinates to geodetic coordinates using either iterative or closed-form methods
- Parameters:
- x, np.ndarray
Cartesian x-coordinates (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- a_axis: float, default 6378137.0
Semi-major 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
- Returns:
- lon: np.ndarray
Longitude (degrees east)
- lat: np.ndarray
Latitude (degrees north)
- h: np.ndarray
Height above ellipsoid (meters)
- pyTMD.spatial._moritz_iterative(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, eps: float = np.float64(2.220446049250313e-16), iterations: int = 10)[source]
Convert from cartesian coordinates to geodetic coordinates using the iterative solution of Hofmann-Wellenhof and Moritz [29]
- Parameters:
- x, np.ndarray
Cartesian x-coordinates (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- a_axis: float, default 6378137.0
Semi-major 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
- pyTMD.spatial._bowring_iterative(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, eps: float = np.float64(2.220446049250313e-16), iterations: int = 10)[source]
Convert from cartesian coordinates to geodetic coordinates using the iterative solution of Bowring [5], Bowring [6]
- Parameters:
- x, np.ndarray
Cartesian x-coordinates (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- a_axis: float, default 6378137.0
Semi-major 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
- 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 Zhu [84]
- Parameters:
- x, np.ndarray
Cartesian x-coordinates (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- a_axis: float, default 6378137.0
Semi-major axis of the ellipsoid
- flat: float, default 1.0/298.257223563
Ellipsoidal flattening
- pyTMD.spatial.to_ENU(x: ndarray, y: ndarray, z: ndarray, lon0: float | ndarray = 0.0, lat0: float | ndarray = 0.0, h0: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]
Convert from Earth-Centered Earth-Fixed (ECEF) cartesian coordinates to East-North-Up coordinates (ENU)
- Parameters:
- x, np.ndarray
Cartesian x-coordinates (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- lon0: float or np.ndarray, default 0.0
Reference longitude (degrees east)
- lat0: float or np.ndarray, default 0.0
Reference latitude (degrees north)
- h0: float or np.ndarray, default 0.0
Reference height (meters)
- a_axis: float, default 6378137.0
Semi-major axis of the ellipsoid
- flat: float, default 1.0/298.257223563
Ellipsoidal flattening
- Returns:
- E: np.ndarray
East coordinates (meters)
- N: np.ndarray
North coordinates (meters)
- U: np.ndarray
Up coordinates (meters)
- pyTMD.spatial.from_ENU(E: ndarray, N: ndarray, U: ndarray, lon0: float | ndarray = 0.0, lat0: float | ndarray = 0.0, h0: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]
Convert from East-North-Up coordinates (ENU) to Earth-Centered Earth-Fixed (ECEF) cartesian coordinates
- Parameters:
- E, np.ndarray
East coordinates (meters)
- N, np.ndarray
North coordinates (meters)
- U, np.ndarray
Up coordinates (meters)
- lon0: float or np.ndarray, default 0.0
Reference longitude (degrees east)
- lat0: float or np.ndarray, default 0.0
Reference latitude (degrees north)
- h0: float or np.ndarray, default 0.0
Reference height (meters)
- a_axis: float, default 6378137.0
Semi-major axis of the ellipsoid
- flat: float, default 1.0/298.257223563
Ellipsoidal flattening
- Returns:
- x, float
Cartesian x-coordinates (meters)
- y, float
Cartesian y-coordinates (meters)
- z, float
Cartesian z-coordinates (meters)
- pyTMD.spatial.to_horizontal(E: ndarray, N: ndarray, U: ndarray)[source]
Convert from East-North-Up coordinates (ENU) to a celestial horizontal coordinate system (alt-az)
- Parameters:
- E: np.ndarray
East coordinates (meters)
- N: np.ndarray
North coordinates (meters)
- U: np.ndarray
Up coordinates (meters)
- Returns:
- alpha: np.ndarray
Altitude (elevation) angle (degrees)
- phi: np.ndarray
Azimuth angle (degrees)
- D: np.ndarray
Distance from observer to object (meters)
- pyTMD.spatial.to_zenith(x: ndarray, y: ndarray, z: ndarray, lon0: float | ndarray = 0.0, lat0: float | ndarray = 0.0, h0: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]
Calculate zenith angle of an object from Earth-Centered Earth-Fixed (ECEF) cartesian coordinates
- Parameters:
- x, np.ndarray
Cartesian x-coordinates (meters)
- y, np.ndarray
Cartesian y-coordinates (meters)
- z, np.ndarray
Cartesian z-coordinates (meters)
- lon0: float or np.ndarray, default 0.0
Reference longitude (degrees east)
- lat0: float or np.ndarray, default 0.0
Reference latitude (degrees north)
- h0: float or np.ndarray, default 0.0
Reference height (meters)
- a_axis: float, default 6378137.0
Semi-major axis of the ellipsoid
- flat: float, default 1.0/298.257223563
Ellipsoidal flattening
- Returns:
- zenith: np.ndarray
Zenith angle of object (degrees)
- pyTMD.spatial.geocentric_latitude(lat: ndarray, flat: float = 0.0033528106647474805)[source]
Compute the geocentric latitude from a geodetic latitude using a simplified empirical relation [69]
- Parameters:
- lat: np.ndarray
Geodetic latitudes (degrees north)
- flat: float, default 1.0/298.257223563
Ellipsoidal flattening
- Returns:
- geolat: np.ndarray
Geocentric latitude (degrees)
- 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 [69]
- 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