isce3 0.25.0
Loading...
Searching...
No Matches
geometry.h
Go to the documentation of this file.
1// -*- C++ -*-
2// -*- coding: utf-8 -*-
3//
4// Author: Bryan Riel , Hirad Ghaemi
5// Copyright 2017-2018 , 2020-2021
6//
14
15#pragma once
16
17#include "forward.h"
18#include <isce3/core/forward.h>
19#include <isce3/product/forward.h>
20
21#include <optional>
22#include <tuple>
23
24#include <Eigen/Dense>
25
26#include <isce3/core/Constants.h>
27#include <isce3/core/Ellipsoid.h>
28#include <isce3/geometry/DEMInterpolator.h>
29
30// Declaration
31namespace isce3 {
33namespace geometry {
34
58int rdr2geo(double aztime, double slantRange, double doppler,
59 const isce3::core::Orbit& orbit,
60 const isce3::core::Ellipsoid& ellipsoid,
61 const DEMInterpolator& demInterp, isce3::core::Vec3& targetLLH,
62 double wvl, isce3::core::LookSide side, double threshold, int maxIter,
63 int extraIter);
64
88int rdr2geo(const isce3::core::Pixel& pixel, const isce3::core::Basis& TCNbasis,
89 const isce3::core::Vec3& pos, const isce3::core::Vec3& vel,
90 const isce3::core::Ellipsoid& ellipsoid,
91 const DEMInterpolator& demInterp, isce3::core::Vec3& targetLLH,
92 isce3::core::LookSide side, double threshold, int maxIter,
93 int extraIter);
94
122int rdr2geo(const isce3::core::Vec3& radarXYZ, const isce3::core::Vec3& axis,
123 double angle, double range, const DEMInterpolator& dem,
124 isce3::core::Vec3& targetXYZ, isce3::core::LookSide side,
125 double threshold, int maxIter, int extraIter);
126
151int geo2rdr(const isce3::core::Vec3& inputLLH,
152 const isce3::core::Ellipsoid& ellipsoid,
153 const isce3::core::Orbit& orbit, const isce3::core::Poly2d& doppler,
154 double& aztime, double& slantRange, double wavelength,
155 double startingRange, double rangePixelSpacing, size_t rwidth,
156 isce3::core::LookSide side, double threshold, int maxIter,
157 double deltaRange);
158
180int geo2rdr(const isce3::core::Vec3& inputLLH,
181 const isce3::core::Ellipsoid& ellipsoid,
182 const isce3::core::Orbit& orbit,
183 const isce3::core::LUT2d<double>& doppler, double& aztime,
184 double& slantRange, double wavelength, isce3::core::LookSide side,
185 double threshold, int maxIter, double deltaRange);
186
205void computeDEMBounds(const isce3::core::Orbit& orbit,
206 const isce3::core::Ellipsoid& ellipsoid,
207 const isce3::core::LUT2d<double>& doppler,
208 const isce3::product::RadarGridParameters& radarGrid, size_t xoff,
209 size_t yoff, size_t xsize, size_t ysize, double margin, double& min_lon,
210 double& min_lat, double& max_lon, double& max_lat);
211
212template<class T>
213double _compute_doppler_aztime_diff(isce3::core::Vec3 dr,
214 isce3::core::Vec3 satvel, T& doppler, double wavelength, double aztime,
215 double slantRange, double deltaRange);
216
227isce3::core::Vec3 nedVector(
228 double lon, double lat, const isce3::core::Vec3& vel);
229
240isce3::core::Vec3 nwuVector(
241 double lon, double lat, const isce3::core::Vec3& vel);
242
253isce3::core::Vec3 enuVector(
254 double lon, double lat, const isce3::core::Vec3& vel);
255
265double heading(double lon, double lat, const isce3::core::Vec3& vel);
266
281double slantRangeFromLookVec(const isce3::core::Vec3& pos,
282 const isce3::core::Vec3& lkvec,
283 const isce3::core::Ellipsoid& ellips = {});
284
316std::pair<int, double> srPosFromLookVecDem(double& sr,
317 isce3::core::Vec3& tg_pos, isce3::core::Vec3& llh,
318 const isce3::core::Vec3& sc_pos, const isce3::core::Vec3& lkvec,
319 const DEMInterpolator& dem_interp = {}, double hgt_err = 0.5,
320 int num_iter = 10, const isce3::core::Ellipsoid& ellips = {},
321 std::optional<double> initial_height = {});
322
357std::tuple<double, double> lookIncAngFromSlantRange(double slant_range,
358 const isce3::core::Orbit& orbit, std::optional<double> az_time = {},
359 const DEMInterpolator& dem_interp = {},
360 const isce3::core::Ellipsoid& ellips = {});
361
396std::tuple<Eigen::ArrayXd, Eigen::ArrayXd> lookIncAngFromSlantRange(
397 const Eigen::Ref<const Eigen::ArrayXd>& slant_range,
398 const isce3::core::Orbit& orbit, std::optional<double> az_time = {},
399 const DEMInterpolator& dem_interp = {},
400 const isce3::core::Ellipsoid& ellips = {});
401
406double compute_mean_dem(const DEMInterpolator& dem);
407
408} // namespace geometry
409} // namespace isce3
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
int geo2rdr(const isce3::core::Vec3 &inputLLH, const isce3::core::Ellipsoid &ellipsoid, const isce3::core::Orbit &orbit, const isce3::core::Poly2d &doppler, double &aztime, double &slantRange, double wavelength, double startingRange, double rangePixelSpacing, size_t rwidth, isce3::core::LookSide side, double threshold, int maxIter, double deltaRange)
Map coordinates to radar geometry coordinates transformer.
Definition geometry.cpp:176
double slantRangeFromLookVec(const isce3::core::Vec3 &pos, const isce3::core::Vec3 &lkvec, const isce3::core::Ellipsoid &ellips={})
Get slant range (m) from platform/antenna position in ECEF (x,y,z) to Reference Ellipsoid given unit ...
Definition geometry.cpp:430
double compute_mean_dem(const DEMInterpolator &dem)
Definition geometry.cpp:619
isce3::core::Vec3 nwuVector(double lon, double lat, const isce3::core::Vec3 &vel)
Get NWU(north,west,up) velocity or unit vector from ECEF velocity or unit vector at a certain geodeti...
Definition geometry.cpp:384
void computeDEMBounds(const isce3::core::Orbit &orbit, const isce3::core::Ellipsoid &ellipsoid, const isce3::core::LUT2d< double > &doppler, const isce3::product::RadarGridParameters &radarGrid, size_t xoff, size_t yoff, size_t xsize, size_t ysize, double margin, double &min_lon, double &min_lat, double &max_lon, double &max_lat)
Utility function to compute geographic bounds for a radar grid.
Definition geometry.cpp:254
std::tuple< double, double > lookIncAngFromSlantRange(double slant_range, const isce3::core::Orbit &orbit, std::optional< double > az_time={}, const DEMInterpolator &dem_interp={}, const isce3::core::Ellipsoid &ellips={})
Estimate look angle (off-nadir angle) and ellipsoidal incidence angle at a desired slant range from o...
Definition geometry.cpp:537
double heading(double lon, double lat, const isce3::core::Vec3 &vel)
Get spacecraft heading/track angle from velocity vector at a certain geodetic location of Spacecraft.
Definition geometry.cpp:399
isce3::core::Vec3 nedVector(double lon, double lat, const isce3::core::Vec3 &vel)
Get unit NED(north,east,down) velocity or unit vector from ECEF velocity or unit vector at a certain ...
Definition geometry.cpp:369
isce3::core::Vec3 enuVector(double lon, double lat, const isce3::core::Vec3 &vel)
Get unit ENU(east,north,up) velocity or unit vector from ECEF velocity or unit vector at a certain ge...
Definition geometry.cpp:391
std::pair< int, double > srPosFromLookVecDem(double &sr, isce3::core::Vec3 &tg_pos, isce3::core::Vec3 &llh, const isce3::core::Vec3 &sc_pos, const isce3::core::Vec3 &lkvec, const DEMInterpolator &dem_interp={}, double hgt_err=0.5, int num_iter=10, const isce3::core::Ellipsoid &ellips={}, std::optional< double > initial_height={})
Get an approximatre ECEF, LLH position and respective Slant range at a certain height above the refer...
Definition geometry.cpp:456
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.