isce3 0.25.0
Loading...
Searching...
No Matches
geometryfunc.h File Reference

A collection of antenna-related geometry functions. More...

#include <isce3/core/forward.h>
#include <tuple>
#include <vector>
#include <Eigen/Dense>
#include <isce3/antenna/Frame.h>
#include <isce3/core/Ellipsoid.h>
#include <isce3/geometry/DEMInterpolator.h>

Go to the source code of this file.

Namespaces

namespace  isce3
 base interpolator is an abstract base class
 

Functions

std::tuple< double, double, bool > isce3::antenna::ant2rgdop (double el_theta, double az_phi, const isce3::core::Vec3 &pos_ecef, const isce3::core::Vec3 &vel_ecef, const isce3::core::Quaternion &quat, double wavelength, const isce3::geometry::DEMInterpolator &dem_interp={}, double abs_tol=0.5, int max_iter=10, const isce3::antenna::Frame &frame={}, const isce3::core::Ellipsoid &ellips={})
 Estimate Radar products, Slant range and Doppler centroid, from spherical angles in antenna body-fixed domain for a certain spacecraft position, velocity,and attitude at a certain height w.r.t.
 
std::tuple< Eigen::VectorXd, Eigen::VectorXd, bool > isce3::antenna::ant2rgdop (const Eigen::Ref< const Eigen::VectorXd > &el_theta, double az_phi, const isce3::core::Vec3 &pos_ecef, const isce3::core::Vec3 &vel_ecef, const isce3::core::Quaternion &quat, double wavelength, const isce3::geometry::DEMInterpolator &dem_interp={}, double abs_tol=0.5, int max_iter=10, const isce3::antenna::Frame &frame={}, const isce3::core::Ellipsoid &ellips={})
 Overloaded function to estimate Radar products, Slant ranges and Doppler centroids, from spherical angles in antenna body-fixed domain for a certain spacecraft position, velocity,and attitude at a certain height w.r.t.
 
std::tuple< isce3::core::Vec3, bool > isce3::antenna::ant2geo (double el_theta, double az_phi, const isce3::core::Vec3 &pos_ecef, const isce3::core::Quaternion &quat, const isce3::geometry::DEMInterpolator &dem_interp={}, double abs_tol=0.5, int max_iter=10, const isce3::antenna::Frame &frame={}, const isce3::core::Ellipsoid &ellips={})
 Estimate geodetic geolocation (longitude, latitude, height) from spherical angles in antenna body-fixed domain for a certain spacecraft position and attitude at a certain height w.r.t.
 
std::tuple< std::vector< isce3::core::Vec3 >, bool > isce3::antenna::ant2geo (const Eigen::Ref< const Eigen::VectorXd > &el_theta, double az_phi, const isce3::core::Vec3 &pos_ecef, const isce3::core::Quaternion &quat, const isce3::geometry::DEMInterpolator &dem_interp={}, double abs_tol=0.5, int max_iter=10, const isce3::antenna::Frame &frame={}, const isce3::core::Ellipsoid &ellips={})
 Overloaded function to estimate geodetic geolocation (longitude, latitude, height) from spherical angles in antenna body-fixed domain for a certain spacecraft position and attitude at a certain height w.r.t.
 
isce3::core::Vec3 isce3::antenna::rangeAzToXyz (double slant_range, double az, const isce3::core::Vec3 &pos_ecef, const isce3::core::Quaternion &quat, const isce3::geometry::DEMInterpolator &dem_interp={}, double el_min=-M_PI/4, double el_max=M_PI/4, double el_tol=0.0, const isce3::antenna::Frame &frame={})
 Compute target position given range and AZ angle by varying EL until height matches DEM.
 

Detailed Description

A collection of antenna-related geometry functions.

Function Documentation

◆ ant2geo() [1/2]

std::tuple< std::vector< Vec3 >, bool > isce3::antenna::ant2geo ( const Eigen::Ref< const Eigen::VectorXd > & el_theta,
double az_phi,
const isce3::core::Vec3 & pos_ecef,
const isce3::core::Quaternion & quat,
const isce3::geometry::DEMInterpolator & dem_interp = {},
double abs_tol = 0.5,
int max_iter = 10,
const isce3::antenna::Frame & frame = {},
const isce3::core::Ellipsoid & ellips = {} )

Overloaded function to estimate geodetic geolocation (longitude, latitude, height) from spherical angles in antenna body-fixed domain for a certain spacecraft position and attitude at a certain height w.r.t.

an ellipsoid.

Parameters
[in]el_theta: a vector of either elevation or theta angle in radians depending on the "frame" object.
[in]az_phi: either azimuth or phi angle in radians depending on the "frame" object.
[in]pos_ecef: antenna/spacecraft position in ECEF (m,m,m)
[in]quat: isce3 quaternion object for transformation from antenna body-fixed to ECEF
[in]dem_interp(optional): isce3 DEMInterpolator object w.r.t ellipsoid. Default is zero height.
[in]abs_tol(optional): Abs error/tolerance in height estimation (m) between desired input height and final output height. Default is 0.5.
[in]max_iter(optional): Max number of iterations in height estimation. Default is 10.
[in]frame(optional): isce3 Frame object to define antenna spherical coordinate system. Default is based on "EL_AND_AZ" spherical grid.
[in]ellips(optional): isce3 Ellipsoid object defining the ellipsoidal planet. Default is WGS84 ellipsoid.
Returns
a vector of isce3::core::Vec3 of geodetic LLH(longitude,latitude,height) in (rad,rad,m)
a bool which is true if height tolerance is met, false otherwise.
Exceptions
InvalidArgument,RuntimeError[ReeTechDesDoc]

◆ ant2geo() [2/2]

std::tuple< Vec3, bool > isce3::antenna::ant2geo ( double el_theta,
double az_phi,
const isce3::core::Vec3 & pos_ecef,
const isce3::core::Quaternion & quat,
const isce3::geometry::DEMInterpolator & dem_interp = {},
double abs_tol = 0.5,
int max_iter = 10,
const isce3::antenna::Frame & frame = {},
const isce3::core::Ellipsoid & ellips = {} )

Estimate geodetic geolocation (longitude, latitude, height) from spherical angles in antenna body-fixed domain for a certain spacecraft position and attitude at a certain height w.r.t.

an ellipsoid.

Parameters
[in]el_theta: either elevation or theta angle in radians depending on the "frame" object.
[in]az_phi: either azimuth or phi angle in radians depending on the "frame" object.
[in]pos_ecef: antenna/spacecraft position in ECEF (m,m,m)
[in]quat: isce3 quaternion object for transformation from antenna body-fixed to ECEF
[in]dem_interp(optional): isce3 DEMInterpolator object w.r.t ellipsoid. Default is zero height.
[in]abs_tol(optional): Abs error/tolerance in height estimation (m) between desired input height and final output height. Default is 0.5.
[in]max_iter(optional): Max number of iterations in height estimation. Default is 10.
[in]frame(optional): isce3 Frame object to define antenna spherical coordinate system. Default is based on "EL_AND_AZ" spherical grid.
[in]ellips(optional): isce3 Ellipsoid object defining the ellipsoidal planet. Default is WGS84 ellipsoid.
Returns
an isce3::core::Vec3 of geodetic LLH(longitude,latitude,height) in (rad,rad,m)
a bool which is true if height tolerance is met, false otherwise.
Exceptions
InvalidArgument,RuntimeError[ReeTechDesDoc]

◆ ant2rgdop() [1/2]

std::tuple< VecXd, VecXd, bool > isce3::antenna::ant2rgdop ( const Eigen::Ref< const Eigen::VectorXd > & el_theta,
double az_phi,
const isce3::core::Vec3 & pos_ecef,
const isce3::core::Vec3 & vel_ecef,
const isce3::core::Quaternion & quat,
double wavelength,
const isce3::geometry::DEMInterpolator & dem_interp = {},
double abs_tol = 0.5,
int max_iter = 10,
const isce3::antenna::Frame & frame = {},
const isce3::core::Ellipsoid & ellips = {} )

Overloaded function to estimate Radar products, Slant ranges and Doppler centroids, from spherical angles in antenna body-fixed domain for a certain spacecraft position, velocity,and attitude at a certain height w.r.t.

an ellipsoid.

Parameters
[in]el_theta: a vector of either elevation or theta angle in radians depending on the "frame" object.
[in]az_phi: either azimuth or phi angle in radians depending on the "frame" object.
[in]pos_ecef: antenna/spacecraft position in ECEF (m,m,m)
[in]vel_ecef: spacecraft velocity in ECEF (m/s,m/s,m/s)
[in]quat: isce3 quaternion object for transformation from antenna body-fixed to ECEF
[in]wavelength: Radar wavelength in (m).
[in]dem_interp(optional): isce3 DEMInterpolator object w.r.t ellipsoid. Default is zero height.
[in]abs_tol(optional): Abs error/tolerance in height estimation (m) between desired input height and final output height. Default is 0.5.
[in]max_iter(optional): Max number of iterations in height estimation. Default is 10.
[in]frame(optional): isce3 Frame object to define antenna spherical coordinate system. Default is based on "EL_AND_AZ" spherical grid.
[in]ellips(optional): isce3 Ellipsoid object defining the ellipsoidal planet. Default is WGS84 ellipsoid.
Returns
an Eigen::VectorXd of slant ranges (m)
an Eigen::VectorXd of Doppler values (Hz)
a bool which is true if height tolerance is met, false otherwise.
Exceptions
InvalidArgument,RuntimeError[ReeTechDesDoc]

◆ ant2rgdop() [2/2]

std::tuple< double, double, bool > isce3::antenna::ant2rgdop ( double el_theta,
double az_phi,
const isce3::core::Vec3 & pos_ecef,
const isce3::core::Vec3 & vel_ecef,
const isce3::core::Quaternion & quat,
double wavelength,
const isce3::geometry::DEMInterpolator & dem_interp = {},
double abs_tol = 0.5,
int max_iter = 10,
const isce3::antenna::Frame & frame = {},
const isce3::core::Ellipsoid & ellips = {} )

Estimate Radar products, Slant range and Doppler centroid, from spherical angles in antenna body-fixed domain for a certain spacecraft position, velocity,and attitude at a certain height w.r.t.

an ellipsoid.

Parameters
[in]el_theta: either elevation or theta angle in radians depending on the "frame" object.
[in]az_phi: either azimuth or phi angle in radians depending on the "frame" object.
[in]pos_ecef: antenna/spacecraft position in ECEF (m,m,m)
[in]vel_ecef: spacecraft velocity in ECEF (m/s,m/s,m/s)
[in]quat: isce3 quaternion object for transformation from antenna body-fixed to ECEF
[in]wavelength: Radar wavelength in (m).
[in]dem_interp(optional): isce3 DEMInterpolator object w.r.t ellipsoid. Default is zero height.
[in]abs_tol(optional): Abs error/tolerance in height estimation (m) between desired input height and final output height. Default is 0.5.
[in]max_iter(optional): Max number of iterations in height estimation. Default is 10.
[in]frame(optional): isce3 Frame object to define antenna spherical coordinate system. Default is based on "EL_AND_AZ" spherical grid.
[in]ellips(optional): isce3 Ellipsoid object defining the ellipsoidal planet. Default is WGS84 ellipsoid.
Returns
slantrange (m)
Doppler (Hz)
a bool which is true if height tolerance is met, false otherwise.
Exceptions
InvalidArgument,RuntimeError[ReeTechDesDoc]

◆ rangeAzToXyz()

Vec3 isce3::antenna::rangeAzToXyz ( double slant_range,
double az,
const isce3::core::Vec3 & pos_ecef,
const isce3::core::Quaternion & quat,
const isce3::geometry::DEMInterpolator & dem_interp = {},
double el_min = -M_PI / 4,
double el_max = M_PI / 4,
double el_tol = 0.0,
const isce3::antenna::Frame & frame = {} )

Compute target position given range and AZ angle by varying EL until height matches DEM.

Parameters
[in]slant_rangeRange to target in m
[in]azAZ angle in rad
[in]pos_ecefECEF XYZ position of radar in m
[in]quatOrientation of the antenna (RCS to ECEF quaternion)
[in]dem_interpDigital elevation model in m above ellipsoid
[in]el_minLower bound for EL solution in rad (default=-45 deg)
[in]el_maxUpper bound for EL solution in rad (default=+45 deg)
[in]el_tolAllowable absolute error in EL solution in rad Zero for maximum possible precision. (default=0)
[in]frameCoordinate convention for (EL, AZ) to cartesian transformation. (default=EL_AND_AZ)
Returns
Target position in ECEF in m

Generated for ISCE3.0 by doxygen 1.13.2.