6#include <isce3/core/forward.h>
13#include <isce3/antenna/Frame.h>
14#include <isce3/core/Ellipsoid.h>
15#include <isce3/geometry/DEMInterpolator.h>
18namespace isce3 {
namespace antenna {
51std::tuple<double, double, bool>
ant2rgdop(
double el_theta,
double az_phi,
52 const isce3::core::Vec3& pos_ecef,
const isce3::core::Vec3& vel_ecef,
53 const isce3::core::Quaternion& quat,
double wavelength,
54 const isce3::geometry::DEMInterpolator& dem_interp = {},
55 double abs_tol = 0.5,
int max_iter = 10,
56 const isce3::antenna::Frame& frame = {},
57 const isce3::core::Ellipsoid& ellips = {});
89std::tuple<Eigen::VectorXd, Eigen::VectorXd, bool>
ant2rgdop(
90 const Eigen::Ref<const Eigen::VectorXd>& el_theta,
double az_phi,
91 const isce3::core::Vec3& pos_ecef,
const isce3::core::Vec3& vel_ecef,
92 const isce3::core::Quaternion& quat,
double wavelength,
93 const isce3::geometry::DEMInterpolator& dem_interp = {},
94 double abs_tol = 0.5,
int max_iter = 10,
95 const isce3::antenna::Frame& frame = {},
96 const isce3::core::Ellipsoid& ellips = {});
127std::tuple<isce3::core::Vec3, bool>
ant2geo(
double el_theta,
double az_phi,
128 const isce3::core::Vec3& pos_ecef,
const isce3::core::Quaternion& quat,
129 const isce3::geometry::DEMInterpolator& dem_interp = {},
130 double abs_tol = 0.5,
int max_iter = 10,
131 const isce3::antenna::Frame& frame = {},
132 const isce3::core::Ellipsoid& ellips = {});
162std::tuple<std::vector<isce3::core::Vec3>,
bool>
ant2geo(
163 const Eigen::Ref<const Eigen::VectorXd>& el_theta,
double az_phi,
164 const isce3::core::Vec3& pos_ecef,
const isce3::core::Quaternion& quat,
165 const isce3::geometry::DEMInterpolator& dem_interp = {},
166 double abs_tol = 0.5,
int max_iter = 10,
167 const isce3::antenna::Frame& frame = {},
168 const isce3::core::Ellipsoid& ellips = {});
187isce3::core::Vec3
rangeAzToXyz(
double slant_range,
double az,
188 const isce3::core::Vec3& pos_ecef,
const isce3::core::Quaternion& quat,
189 const isce3::geometry::DEMInterpolator& dem_interp = {},
190 double el_min = -M_PI / 4,
double el_max = M_PI / 4,
191 double el_tol = 0.0,
const isce3::antenna::Frame& frame = {});
std::tuple< double, double, bool > 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-fixe...
Definition geometryfunc.cpp:70
isce3::core::Vec3 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.
Definition geometryfunc.cpp:164
std::tuple< isce3::core::Vec3, bool > 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-fix...
Definition geometryfunc.cpp:119
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5