Class inertialsim::geodesy::Coordinates¶
ClassList > inertialsim > geodesy > Coordinates
3-dimensional coordinates with respect to Earth. More...
#include <coordinates.h>
Public Functions¶
| Type | Name |
|---|---|
| Vector3D | AsEcef () const Return ECEF coordinates. |
| Array3D | AsEllipsoidalHarmonic () const Return ellipsoidal harmonic coordinates. |
| Array3D | AsGeodetic () const Return geodetic coordinates. |
| Array3D | AsSpherical () const Return spherical coordinates. |
| Vector3D | AsTopocentric (const std::optional< TopocentricDatum > & datum=std::nullopt) const Return topocentric coordinates. |
| Coordinates | Interpolate (const Timestamps & time, geometry::Interpolator method=geometry::Interpolator::LINEAR) const Interpolate coordinates at requested times. |
| geometry::RigidTransform | RigidTransform (CoordinateSystem from_coordinates, CoordinateSystem to_coordinates, const std::optional< TopocentricDatum > & from_topocentric_datum=std::nullopt, const std::optional< TopocentricDatum > & to_topocentric_datum=std::nullopt) const Return rigid transform between coordinate systems. |
| Coordinates | Slice (int start, int end) const Get a slice of coordinates. |
| Vector3D | TimeDerivative (int order=1, int accuracy=4) const Numerical time derivative in current coordinate system. |
| CoordinateSystem | coordinate_system () const The coordinate system of the stored coordinates. |
| const Array3D & | coordinates () const Array of raw coordinates. |
| std::optional< GeodeticDatum > | geodetic_datum () const The geodetic datum associated with the stored coordinates. |
| int | num_coordinates () const The number of coordinates stored in the instance. |
| Coordinates | operator[] (int index) const Get a single coordinate by index. |
| std::optional< Timestamps > | time () const Array of times corresponding to each coordinate (may be empty). |
| std::optional< TopocentricDatum > | topocentric_datum () const The topocentric datum associated with the stored coordinates. |
Public Static Functions¶
| Type | Name |
|---|---|
| Coordinates | FromEcef (const Vector3D & xyz, const GeodeticDatum & datum=WGS84, const std::optional< Timestamps > & time=std::nullopt) Construct from ECEF coordinates. |
| Coordinates | FromEllipsoidalHarmonic (const Array3D & upl, const GeodeticDatum & datum=WGS84, const std::optional< Timestamps > & time=std::nullopt) Construct from ellipsoidal harmonic coordinates. |
| Coordinates | FromGeodetic (const Array3D & lla, const GeodeticDatum & datum=WGS84, const std::optional< Timestamps > & time=std::nullopt) Construct from geodetic coordinates. |
| Coordinates | FromInertial (const Vector3D & xyz, const std::optional< Timestamps > & time=std::nullopt) Construct from inertial coordinates. |
| Coordinates | FromLocal (const Vector3D & position, const std::optional< Timestamps > & time=std::nullopt) Construct from simple local coordinates. |
| Coordinates | FromLocal (const Vector3D & position, const Array3D & frame_origin, const geometry::Rotation & frame_rotation, const LocalDatum & local_datum, const std::optional< Timestamps > & time=std::nullopt) Construct from local coordinates with geodetic frame. |
| Coordinates | FromSpherical (const Array3D & rpl, const GeodeticDatum & datum=WGS84, const std::optional< Timestamps > & time=std::nullopt) Construct from spherical coordinates. |
| Coordinates | FromTopocentric (const Vector3D & ned, const TopocentricDatum & datum, const std::optional< Timestamps > & time=std::nullopt) Construct from topocentric coordinates. |
Protected Functions¶
| Type | Name |
|---|---|
| Coordinates () = default Default constructor. |
|
| Coordinates (const Array3D & coordinates, CoordinateSystem coordinate_system, const std::optional< GeodeticDatum > & geodetic_datum=std::nullopt, const std::optional< TopocentricDatum > & topocentric_datum=std::nullopt, const std::optional< Timestamps > & time=std::nullopt) Construct with raw data. |
Detailed Description¶
Geodetic coordinates describe the location of points relative to the Earth. In contrast, see geometry::Vector for basic geometric coordinates.
Use the From* methods to construct an instance. Do not initialize an instance directly using the protected constructor.
Warning:
Coordinate conventions: There are a variety of coordinate conventions and they are a common source of errors. Because many of these conventions cannot be automatically validated, users must be aware of them and manage interfaces with other software accordingly. InertialSim conventions are documented in each method.
** **
Consistent with inertial navigation practice, a coordinate system is defined as an axis system (an orthogonal Cartesian basis) and an origin point relative to Earth.
Note:
This class stores an internal ECEF representation for conversions between coordinate systems.
Note:
C++ Implementation Uses 3xN Matrix Layout: Unlike the Python implementation which uses Nx3 layout, the C++ implementation stores coordinates as 3xN Eigen matrices where each column represents a coordinate point. This is consistent with the C++ geometry module convention and optimizes for Eigen's column-major storage.
See also: geometry::Vector for basic geometric coordinates
See also: GlobalPose for full pose (attitude + position + velocity)
See also: FromGeodetic(), FromEcef(), FromTopocentric() for factory methods
Public Functions Documentation¶
function AsEcef¶
Return ECEF coordinates.
Returns:
ECEF coordinates as 3xN matrix.
function AsEllipsoidalHarmonic¶
Return ellipsoidal harmonic coordinates.
Returns:
Ellipsoidal harmonic coordinates (semi-minor axis, polar angle, longitude) as 3xN matrix.
function AsGeodetic¶
Return geodetic coordinates.
Returns:
Geodetic coordinates (latitude, longitude, altitude) as 3xN matrix.
function AsSpherical¶
Return spherical coordinates.
Returns:
Spherical coordinates (radius, polar angle, longitude) as 3xN matrix.
function AsTopocentric¶
Return topocentric coordinates.
Vector3D inertialsim::geodesy::Coordinates::AsTopocentric (
const std::optional< TopocentricDatum > & datum=std::nullopt
) const
Parameters:
datumOptional topocentric datum. If not provided, uses the datum from construction.
Returns:
Topocentric coordinates (NED or ENU) as 3xN matrix.
function Interpolate¶
Interpolate coordinates at requested times.
Coordinates inertialsim::geodesy::Coordinates::Interpolate (
const Timestamps & time,
geometry::Interpolator method=geometry::Interpolator::LINEAR
) const
function RigidTransform¶
Return rigid transform between coordinate systems.
geometry::RigidTransform inertialsim::geodesy::Coordinates::RigidTransform (
CoordinateSystem from_coordinates,
CoordinateSystem to_coordinates,
const std::optional< TopocentricDatum > & from_topocentric_datum=std::nullopt,
const std::optional< TopocentricDatum > & to_topocentric_datum=std::nullopt
) const
function Slice¶
Get a slice of coordinates.
Parameters:
startStart index (inclusive).endEnd index (exclusive).
Returns:
Coordinates containing the sliced elements.
function TimeDerivative¶
Numerical time derivative in current coordinate system.
function coordinate_system¶
The coordinate system of the stored coordinates.
function coordinates¶
Array of raw coordinates.
function geodetic_datum¶
The geodetic datum associated with the stored coordinates.
function num_coordinates¶
The number of coordinates stored in the instance.
function operator[]¶
Get a single coordinate by index.
Parameters:
indexIndex of the coordinate to retrieve.
Returns:
Coordinates containing a single coordinate.
function time¶
Array of times corresponding to each coordinate (may be empty).
function topocentric_datum¶
The topocentric datum associated with the stored coordinates.
inline std::optional< TopocentricDatum > inertialsim::geodesy::Coordinates::topocentric_datum () const
Public Static Functions Documentation¶
function FromEcef¶
Construct from ECEF coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromEcef (
const Vector3D & xyz,
const GeodeticDatum & datum=WGS84,
const std::optional< Timestamps > & time=std::nullopt
)
The z-axis points along Earth's rotation axis (through the IERS Reference Pole), the x-axis is the line between the equatorial plane and the zero-meridian plane (IERS Reference Meridian) and the y-axis forms a right-hand orthonormal basis.
** **
Virtually all references use the convention described above, but a notable exception is Reference [01] which places the y-axis along Earth's rotation axis.
Parameters:
xyzECEF coordinates as 3xN matrix (meters).datumA datum describing Earth's shape and gravity field.timeOptional time array of unique, strictly increasing times.
Returns:
Coordinates instance.
function FromEllipsoidalHarmonic¶
Construct from ellipsoidal harmonic coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromEllipsoidalHarmonic (
const Array3D & upl,
const GeodeticDatum & datum=WGS84,
const std::optional< Timestamps > & time=std::nullopt
)
Parameters:
uplEllipsoidal harmonic coordinates (semi-minor axis, polar angle, longitude) as 3xN matrix.datumA datum describing Earth's shape and gravity field.timeOptional time array.
Returns:
Coordinates instance.
function FromGeodetic¶
Construct from geodetic coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromGeodetic (
const Array3D & lla,
const GeodeticDatum & datum=WGS84,
const std::optional< Timestamps > & time=std::nullopt
)
** **
Geodetic coordinates are traditionally ordered latitude, longitude, altitude. This is the convention that InertialSim adopts. Users should note that the local Cartesian basis of these coordinates is a left-handed system.
** **
Latitudes north of the equator are positive, latitudes south of the equator are negative. Longitude is measured relative to the zero-meridian (the IERS Reference Meridian, also known as the Greenwich meridian). Longitudes east of the zero-meridian are positive, longitudes west of the zero-meridian are negative. To be valid, the latitude and longitude coordinates must meet the following requirements: \(\[ -\pi/2 \leq latitude \leq \pi/2 \\ -\pi < longitude \leq \pi \text{ or } 0 \leq longitude < 2\pi \]\) To standardize, the longitude will be adjusted such that: \(-\pi < longitude \leq \pi\)
Altitude is measured relative to the ellipsoid. Positive altitude is above the ellipsoid, negative altitude is below the ellipsoid.
Parameters:
llaGeodetic coordinates (latitude, longitude, altitude) as 3xN matrix. Latitude and longitude in radians, altitude in meters.datumA datum describing Earth's shape and gravity field.timeOptional time array of unique, strictly increasing times.
Returns:
Coordinates instance.
function FromInertial¶
Construct from inertial coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromInertial (
const Vector3D & xyz,
const std::optional< Timestamps > & time=std::nullopt
)
Parameters:
xyzInertial coordinates (3xN matrix).timeOptional time array.
Returns:
Coordinates instance.
function FromLocal [½]¶
Construct from simple local coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromLocal (
const Vector3D & position,
const std::optional< Timestamps > & time=std::nullopt
)
Creates local coordinates without a geodetic datum or transformation.
Parameters:
positionLocal coordinates (3xN matrix).timeOptional time array.
Returns:
Coordinates instance.
function FromLocal [2/2]¶
Construct from local coordinates with geodetic frame.
static Coordinates inertialsim::geodesy::Coordinates::FromLocal (
const Vector3D & position,
const Array3D & frame_origin,
const geometry::Rotation & frame_rotation,
const LocalDatum & local_datum,
const std::optional< Timestamps > & time=std::nullopt
)
Creates local coordinates relative to a geodetic origin with rotation.
Parameters:
positionLocal coordinates (3xN matrix).frame_originGeodetic origin for the local frame.frame_rotationRotation from topocentric to local axes.local_datumAxis ordering for the local frame.timeOptional time array.
Returns:
Coordinates instance.
function FromSpherical¶
Construct from spherical coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromSpherical (
const Array3D & rpl,
const GeodeticDatum & datum=WGS84,
const std::optional< Timestamps > & time=std::nullopt
)
Parameters:
rplSpherical coordinates (radius, polar angle, longitude) as 3xN matrix.datumA datum describing Earth's shape and gravity field.timeOptional time array.
Returns:
Coordinates instance.
function FromTopocentric¶
Construct from topocentric coordinates.
static Coordinates inertialsim::geodesy::Coordinates::FromTopocentric (
const Vector3D & ned,
const TopocentricDatum & datum,
const std::optional< Timestamps > & time=std::nullopt
)
Parameters:
nedTopocentric coordinates (NED or ENU) as 3xN matrix.datumA datum describing the origin and axis order.timeOptional time array.
Returns:
Coordinates instance.
Protected Functions Documentation¶
function Coordinates [½]¶
Default constructor.
function Coordinates [2/2]¶
Construct with raw data.
inline inertialsim::geodesy::Coordinates::Coordinates (
const Array3D & coordinates,
CoordinateSystem coordinate_system,
const std::optional< GeodeticDatum > & geodetic_datum=std::nullopt,
const std::optional< TopocentricDatum > & topocentric_datum=std::nullopt,
const std::optional< Timestamps > & time=std::nullopt
)
Use the factory methods (FromGeodetic, FromEcef, etc.) instead of calling this constructor directly.
Parameters:
coordinatesRaw coordinate data (3xN matrix).coordinate_systemThe coordinate system of the data.geodetic_datumOptional geodetic datum.topocentric_datumOptional topocentric datum.timeOptional time array.
The documentation for this class was generated from the following file cpp/include/inertialsim/geodesy/coordinates.h