3#include <isce3/core/forward.h>
7#include <isce3/core/Common.h>
8#include <isce3/core/LookSide.h>
9#include <isce3/error/ErrorCode.h>
16 double threshold = 1e-8;
47template<
class Orbit,
class DEMInterpolator>
48CUDA_HOSTDEV isce3::error::ErrorCode
49rdr2geo(isce3::core::Vec3* llh,
double t,
double r,
double fd,
52 isce3::core::LookSide side,
double h0 = 0.,
75template<
class DEMInterpolator>
76CUDA_HOSTDEV isce3::error::ErrorCode
81 double h0 = 0.,
const Rdr2GeoParams& params = {});
85inline constexpr double DEFAULT_TOL_HEIGHT = 1e-5;
90 double tol_height = DEFAULT_TOL_HEIGHT;
93 double look_min = 0.0;
96 double look_max = M_PI / 2;
125template<
class Orbit,
class DEMInterpolator>
126CUDA_HOSTDEV isce3::error::ErrorCode
127rdr2geo_bracket(isce3::core::Vec3* xyz,
128 double aztime,
double slantRange,
double doppler,
const Orbit& orbit,
130 double wavelength, isce3::core::LookSide side,
135#include "Rdr2Geo.icc"
Definition DEMInterpolator.h:25
Sequence of platform ephemeris samples (state vectors) with uniform temporal spacing,...
Definition Orbit.h:44
Simple class to store three-dimensional basis vectors.
Definition Basis.h:20
Data structure to store Ellipsoid information.
Definition Ellipsoid.h:20
Helper datastructure to handle slant range information for a pixel.
Definition Pixel.h:13
Definition DEMInterpolator.h:25
The isce3::geometry namespace.
Definition boundingbox.h:15
int rdr2geo(double aztime, double slantRange, double doppler, const isce3::core::Orbit &orbit, const isce3::core::Ellipsoid &ellipsoid, const DEMInterpolator &demInterp, isce3::core::Vec3 &targetLLH, double wvl, isce3::core::LookSide side, double threshold, int maxIter, int extraIter)
Radar geometry coordinates to map coordinates transformer.
Definition geometry.cpp:39
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5