The isce3::geometry namespace. More...
Classes | |
class | DEMInterpolator |
class | Geocode |
class | Geo2rdr |
Transformer from map coordinates to radar geometry coordinates. More... | |
class | Topo |
Transformer from radar geometry coordinates to map coordinates with DEM / reference altitude. More... | |
class | TopoLayers |
Typedefs | |
typedef OGRLinearRing | Perimeter |
Same as GDAL's OGRLinearRing structure. More... | |
typedef OGREnvelope | BoundingBox |
Same as GDAL's OGREnvelope structure. More... | |
typedef OGRTriangle | Triangle |
Same as GDAL's OGRTriangle structure. More... | |
Enumerations | |
enum | geocodeOutputMode { INTERP = 0, AREA_PROJECTION = 1, AREA_PROJECTION_GAMMA_NAUGHT = 2 } |
Enumeration type to indicate the algorithm used for geocoding. | |
enum | geocodeMemoryMode { AUTO = 0, SINGLE_BLOCK = 1, BLOCKS_GEOGRID = 2, BLOCKS_GEOGRID_AND_RADARGRID = 3 } |
Enumeration type to indicate memory management. | |
enum | rtcInputRadiometry { BETA_NAUGHT = 0, SIGMA_NAUGHT_ELLIPSOID = 1 } |
Enumeration type to indicate input terrain radiometry (for RTC) | |
enum | rtcMemoryMode { RTC_AUTO = 0, RTC_SINGLE_BLOCK = 1, RTC_BLOCKS_GEOGRID = 2 } |
Enumeration type to indicate memory management. | |
enum | rtcAreaMode { AREA = 0, AREA_FACTOR = 1 } |
Enumeration type to indicate RTC area mode (AREA or AREA_FACTOR) | |
enum | rtcAlgorithm { RTC_DAVID_SMALL = 0, RTC_AREA_PROJECTION = 1 } |
Enumeration type to select RTC algorithm (RTC_DAVID_SMALL or RTC_AREA_PROJECTION) | |
Functions | |
Perimeter | getGeoPerimeter (const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::ProjectionBase *proj, const isce3::core::LUT2d< double > &doppler={}, const DEMInterpolator &demInterp=DEMInterpolator(0.), const int pointsPerEdge=11, const double threshold=1.0e-8, const int numiter=15) |
Compute the perimeter of a radar grid in map coordinates. More... | |
BoundingBox | getGeoBoundingBox (const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::ProjectionBase *proj, const isce3::core::LUT2d< double > &doppler={}, const std::vector< double > &hgts={isce3::core::GLOBAL_MIN_HEIGHT, isce3::core::GLOBAL_MAX_HEIGHT}, const double margin=0.0, const int pointsPerEdge=11, const double threshold=1.0e-8, const int numiter=15, bool ignore_out_of_range_exception=false) |
Compute bounding box using min/ max altitude for quick estimates. More... | |
BoundingBox | getGeoBoundingBoxHeightSearch (const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::ProjectionBase *proj, const isce3::core::LUT2d< double > &doppler={}, double min_height=isce3::core::GLOBAL_MIN_HEIGHT, double max_height=isce3::core::GLOBAL_MAX_HEIGHT, const double margin=0.0, const int pointsPerEdge=11, const double threshold=1.0e-8, const int numiter=15, const double height_threshold=100) |
Compute bounding box with auto search within given min/ max height interval. More... | |
template<typename T1 , typename T2 > | |
auto | operator* (const std::complex< T1 > &lhs, const T2 &rhs) |
template<typename T1 , typename T2 > | |
auto | operator* (const T1 &lhs, const std::complex< T2 > &rhs) |
template<typename T , typename T_out > | |
void | _convertToOutputType (T a, T_out &b) |
template<typename T , typename T_out > | |
void | _convertToOutputType (std::complex< T > a, T_out &b) |
template<typename T , typename T_out > | |
void | _convertToOutputType (std::complex< T > a, std::complex< T_out > &b) |
template<typename T , typename T_out > | |
void | _accumulate (T_out &band_value, T a, double b) |
template<typename T > | |
constexpr bool | is_complex () |
std::vector< float > | getGeoAreaElementMean (const std::vector< double > &x_vect, const std::vector< double > &y_vect, const isce3::product::RadarGridParameters &radar_grid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &input_dop, isce3::io::Raster &input_raster, isce3::io::Raster &dem_raster, isce3::geometry::rtcInputRadiometry input_radiometry, int exponent, geocodeOutputMode output_mode, double geogrid_upsampling, float rtc_min_value_db, double abs_cal_factor, float radar_grid_nlooks, float *out_nlooks, isce3::core::dataInterpMethod interp_method, double threshold, int num_iter, double delta_range) |
template<typename T > | |
std::vector< float > | _getGeoAreaElementMean (const std::vector< double > &r_vect, const std::vector< double > &a_vect, int x_min, int y_min, isce3::core::Matrix< float > &rtc_area, const isce3::product::RadarGridParameters &radar_grid, isce3::io::Raster &input_raster, geocodeOutputMode output_mode, float rtc_min_value, float *out_nlooks, double abs_cal_factor, float radar_grid_nlooks) |
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. More... | |
int | rdr2geo (const isce3::core::Pixel &pixel, const isce3::core::Basis &TCNbasis, const isce3::core::Vec3 &pos, const isce3::core::Vec3 &vel, const isce3::core::Ellipsoid &ellipsoid, const DEMInterpolator &demInterp, isce3::core::Vec3 &targetLLH, isce3::core::LookSide side, double threshold, int maxIter, int extraIter) |
Radar geometry coordinates to map coordinates transformer. More... | |
int | rdr2geo (const isce3::core::Vec3 &radarXYZ, const isce3::core::Vec3 &axis, double angle, double range, const DEMInterpolator &dem, isce3::core::Vec3 &targetXYZ, isce3::core::LookSide side, double threshold, int maxIter, int extraIter) |
"Cone" interface to rdr2geo. More... | |
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. More... | |
int | geo2rdr (const isce3::core::Vec3 &inputLLH, const isce3::core::Ellipsoid &ellipsoid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &doppler, double &aztime, double &slantRange, double wavelength, isce3::core::LookSide side, double threshold, int maxIter, double deltaRange) |
Map coordinates to radar geometry coordinates transformer. More... | |
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. More... | |
template<class T > | |
double | _compute_doppler_aztime_diff (isce3::core::Vec3 dr, isce3::core::Vec3 satvel, T &doppler, double wavelength, double aztime, double slantRange, double deltaRange) |
int | _omp_thread_count () |
template<typename T > | |
void | _applyRTC (isce3::io::Raster &input_raster, isce3::io::Raster &input_rtc, isce3::io::Raster &output_raster, float rtc_min_value, double abs_cal_factor, pyre::journal::info_t &info, bool flag_complex_to_real_squared) |
void | applyRTC (const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &dop, isce3::io::Raster &input_raster, isce3::io::Raster &dem_raster, isce3::io::Raster &output_raster, rtcInputRadiometry inputRadiometry=rtcInputRadiometry::BETA_NAUGHT, int exponent=0, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), double abs_cal_factor=1, float radar_grid_nlooks=1, isce3::io::Raster *out_nlooks=nullptr, isce3::io::Raster *input_rtc=nullptr, isce3::io::Raster *output_rtc=nullptr, rtcMemoryMode rtc_memory_mode=rtcMemoryMode::RTC_AUTO) |
Apply radiometric terrain correction (RTC) over an input raster. More... | |
double | computeUpsamplingFactor (const DEMInterpolator &dem_interp, const isce3::product::RadarGridParameters &radar_grid, const isce3::core::Ellipsoid &ellps) |
int | areaProjGetNBlocks (int array_length, pyre::journal::info_t *channel, int upsampling, int *block_length_with_upsampling, int *block_length, int min_block_length, int max_block_length) |
void | facetRTC (isce3::product::Product &product, isce3::io::Raster &dem_raster, isce3::io::Raster &output_raster, char frequency= 'A', bool native_doppler=false, rtcInputRadiometry inputRadiometry=rtcInputRadiometry::BETA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), size_t nlooks_az=1, size_t nlooks_rg=1, isce3::io::Raster *out_nlooks=nullptr, rtcMemoryMode rtc_memory_mode=rtcMemoryMode::RTC_AUTO) |
Generate radiometric terrain correction (RTC) area or area factor. More... | |
void | facetRTC (const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &dop, isce3::io::Raster &dem, isce3::io::Raster &output_raster, rtcInputRadiometry inputRadiometry=rtcInputRadiometry::BETA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), float radar_grid_nlooks=1, isce3::io::Raster *out_nlooks=nullptr, rtcMemoryMode rtc_memory_mode=rtcMemoryMode::RTC_AUTO, isce3::core::dataInterpMethod interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, double threshold=1e-4, int num_iter=100, double delta_range=1e-4) |
Generate radiometric terrain correction (RTC) area or area factor. More... | |
void | facetRTC (isce3::io::Raster &dem_raster, isce3::io::Raster &output_raster, const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &dop, const double y0, const double dy, const double x0, const double dx, const int geogrid_length, const int geogrid_width, const int epsg, rtcInputRadiometry inputRadiometry=rtcInputRadiometry::BETA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), float radar_grid_nlooks=1, isce3::io::Raster *out_geo_vertices=nullptr, isce3::io::Raster *out_geo_grid=nullptr, isce3::io::Raster *out_nlooks=nullptr, rtcMemoryMode rtc_memory_mode=rtcMemoryMode::RTC_AUTO, isce3::core::dataInterpMethod interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, double threshold=1e-4, int num_iter=100, double delta_range=1e-4) |
Generate radiometric terrain correction (RTC) area or area factor. More... | |
void | areaProjIntegrateSegment (double y1, double y2, double x1, double x2, int length, int width, isce3::core::Matrix< double > &w_arr, double &nlooks, int plane_orientation) |
void | _addArea (double area, isce3::core::Matrix< float > &out_array, float radar_grid_nlooks, isce3::core::Matrix< float > &out_nlooks_array, int length, int width, int x_min, int y_min, int size_x, int size_y, isce3::core::Matrix< double > &w_arr, double nlooks, isce3::core::Matrix< double > &w_arr_out, double &nlooks_out, double x_center, double x_left, double x_right, double y_center, double y_left, double y_right, int plane_orientation) |
double | computeFacet (Vec3 xyz_center, Vec3 xyz_left, Vec3 xyz_right, Vec3 lookXYZ, double p1, double &p3, double divisor, bool clockwise_direction) |
void | facetRTCDavidSmall (isce3::io::Raster &dem_raster, isce3::io::Raster &output_raster, const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &dop, const double y0, const double dy, const double x0, const double dx, const int geogrid_length, const int geogrid_width, const int epsg, rtcInputRadiometry inputRadiometry=rtcInputRadiometry::BETA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN()) |
Generate radiometric terrain correction (RTC) area or area factor using the David Small algorithm. More... | |
void | _RunBlock (const int jmax, int block_size, int block_size_with_upsampling, int block, int &numdone, int progress_block, double geogrid_upsampling, isce3::core::dataInterpMethod interp_method, isce3::io::Raster &dem_raster, isce3::io::Raster *out_geo_vertices, isce3::io::Raster *out_geo_grid, const double start, const double pixazm, const double dr, double r0, int xbound, int ybound, const double y0, const double dy, const double x0, const double dx, const int geogrid_length, const int geogrid_width, const isce3::product::RadarGridParameters &radar_grid, const isce3::core::LUT2d< double > &dop, const isce3::core::Ellipsoid &ellipsoid, const isce3::core::Orbit &orbit, double threshold, int num_iter, double delta_range, isce3::core::Matrix< float > &out_array, isce3::core::Matrix< float > &out_nlooks_array, isce3::core::ProjectionBase *proj, rtcAreaMode rtc_area_mode, rtcInputRadiometry input_radiometry, float radar_grid_nlooks) |
void | facetRTCAreaProj (isce3::io::Raster &dem, isce3::io::Raster &output_raster, const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &dop, const double y0, const double dy, const double x0, const double dx, const int geogrid_length, const int geogrid_width, const int epsg, rtcInputRadiometry inputRadiometry=rtcInputRadiometry::BETA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), float radar_grid_nlooks=1, isce3::io::Raster *out_geo_vertices=nullptr, isce3::io::Raster *out_geo_grid=nullptr, isce3::io::Raster *out_nlooks=nullptr, rtcMemoryMode rtc_memory_mode=rtcMemoryMode::RTC_AUTO, isce3::core::dataInterpMethod interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, double threshold=1e-4, int num_iter=100, double delta_range=1e-4) |
Generate radiometric terrain correction (RTC) area or area factor using the area projection algorithms. More... | |
std::string | get_input_radiometry_str (rtcInputRadiometry input_radiometry) |
Convert enum input_radiometry to string. | |
std::string | get_rtc_area_mode_str (rtcAreaMode rtc_area_mode) |
Convert enum output_mode to string. | |
std::string | get_rtc_algorithm_str (rtcAlgorithm rtc_algorithm) |
Convert enum output_mode to string. | |
void | print_parameters (pyre::journal::info_t &channel, const isce3::product::RadarGridParameters &radar_grid, const double y0, const double dy, const double x0, const double dx, const int geogrid_length, const int geogrid_width, rtcInputRadiometry input_radiometry, rtcAreaMode rtc_area_mode, double geogrid_upsampling, float rtc_min_value_db) |
template<typename T > | |
void | load_archive (std::string metadata, char *objectTag, T *object) |
template<class Archive > | |
void | save (Archive &, const Topo &) |
template<class Archive > | |
void | load (Archive &archive, Topo &topo) |
template<class Archive > | |
void | save (Archive &, const Geo2rdr &) |
template<class Archive > | |
void | load (Archive &archive, Geo2rdr &geo) |
The isce3::geometry namespace.
typedef OGREnvelope isce3::geometry::BoundingBox |
Same as GDAL's OGREnvelope structure.
typedef OGRLinearRing isce3::geometry::Perimeter |
Same as GDAL's OGRLinearRing structure.
typedef OGRTriangle isce3::geometry::Triangle |
Same as GDAL's OGRTriangle structure.
void isce3::geometry::applyRTC | ( | const isce3::product::RadarGridParameters & | radarGrid, |
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | dop, | ||
isce3::io::Raster & | input_raster, | ||
isce3::io::Raster & | dem_raster, | ||
isce3::io::Raster & | output_raster, | ||
rtcInputRadiometry | inputRadiometry = rtcInputRadiometry::BETA_NAUGHT , |
||
int | exponent = 0 , |
||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR , |
||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION , |
||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN() , |
||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN() , |
||
double | abs_cal_factor = 1 , |
||
float | radar_grid_nlooks = 1 , |
||
isce3::io::Raster * | out_nlooks = nullptr , |
||
isce3::io::Raster * | input_rtc = nullptr , |
||
isce3::io::Raster * | output_rtc = nullptr , |
||
rtcMemoryMode | rtc_memory_mode = rtcMemoryMode::RTC_AUTO |
||
) |
Apply radiometric terrain correction (RTC) over an input raster.
[in] | radarGrid | Radar Grid |
[in] | orbit | Orbit |
[in] | input_dop | Doppler LUT |
[in] | input_raster | Input raster |
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | input_radiometry | Terrain radiometry of the input raster |
[in] | exponent | Exponent to be applied to the input data. The value 0 indicates that the the exponent is based on the data type of the input raster (1 for real and 2 for complex rasters). |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_algorithm | RTC algorithm (RTC_DAVID_SMALL or RTC_AREA_PROJECTION) |
[in] | geogrid_upsampling | Geogrid upsampling (in each direction) |
[in] | rtc_min_value_db | Minimum value for the RTC area factor. Radar data with RTC area factor below this limit are ignored. |
[in] | abs_cal_factor | Absolute calibration factor. |
[in] | radar_grid_nlooks | Radar grid number of looks. This parameters determines the multilooking factor used to compute out_nlooks. |
[out] | out_nlooks | Raster to which the number of radar-grid looks associated with the geogrid will be saved |
[in] | input_rtc | Raster containing pre-computed RTC area factor |
[out] | output_rtc | Output RTC area factor |
[in] | rtc_memory_mode | Select memory mode |
void isce3::geometry::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.
[in] | orbit | Orbit object |
[in] | ellipsoid | Ellipsoid object |
[in] | doppler | LUT2d doppler object |
[in] | lookSide | Left or Right |
[in] | radarGrid | RadarGridParameters object |
[in] | xoff | Column index of radar subwindow |
[in] | yoff | Row index of radar subwindow |
[in] | xsize | Number of columns of radar subwindow |
[in] | ysize | Number of rows of radar subwindiw |
[in] | margin | Padding of extracted DEM (radians) |
[out] | min_lon | Minimum longitude of geographic region (radians) |
[out] | min_lat | Minimum latitude of geographic region (radians) |
[out] | max_lon | Maximum longitude of geographic region (radians) |
[out] | max_lat | Maximum latitude of geographic region (radians) |
void isce3::geometry::facetRTC | ( | isce3::product::Product & | product, |
isce3::io::Raster & | dem_raster, | ||
isce3::io::Raster & | output_raster, | ||
char | frequency = 'A' , |
||
bool | native_doppler = false , |
||
rtcInputRadiometry | inputRadiometry = rtcInputRadiometry::BETA_NAUGHT , |
||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR , |
||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION , |
||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN() , |
||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN() , |
||
size_t | nlooks_az = 1 , |
||
size_t | nlooks_rg = 1 , |
||
isce3::io::Raster * | out_nlooks = nullptr , |
||
rtcMemoryMode | rtc_memory_mode = rtcMemoryMode::RTC_AUTO |
||
) |
Generate radiometric terrain correction (RTC) area or area factor.
[in] | product | Product |
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | frequency | Product frequency |
[in] | native_doppler | Use native doppler |
[in] | input_radiometry | Terrain radiometry of the input raster |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_algorithm | RTC algorithm (RTC_DAVID_SMALL or RTC_AREA_PROJECTION) |
[in] | geogrid_upsampling | Geogrid upsampling (in each direction) |
[in] | rtc_min_value_db | Minimum value for the RTC area factor. Radar data with RTC area factor below this limit are ignored. |
[in] | nlooks_az | Number of azimuth looks. |
[in] | nlooks_rg | Number of range looks. |
[out] | out_nlooks | Raster to which the number of radar-grid |
[in] | rtc_memory_mode | Select memory mode looks associated with the geogrid will be saved |
void isce3::geometry::facetRTC | ( | const isce3::product::RadarGridParameters & | radarGrid, |
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | dop, | ||
isce3::io::Raster & | dem, | ||
isce3::io::Raster & | output_raster, | ||
rtcInputRadiometry | inputRadiometry = rtcInputRadiometry::BETA_NAUGHT , |
||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR , |
||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION , |
||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN() , |
||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN() , |
||
float | radar_grid_nlooks = 1 , |
||
isce3::io::Raster * | out_nlooks = nullptr , |
||
rtcMemoryMode | rtc_memory_mode = rtcMemoryMode::RTC_AUTO , |
||
isce3::core::dataInterpMethod | interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD , |
||
double | threshold = 1e-4 , |
||
int | num_iter = 100 , |
||
double | delta_range = 1e-4 |
||
) |
Generate radiometric terrain correction (RTC) area or area factor.
[in] | product | Product |
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | frequency | Product frequency |
[in] | input_radiometry | Terrain radiometry of the input raster |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_algorithm | RTC algorithm (RTC_DAVID_SMALL or RTC_AREA_PROJECTION) |
[in] | geogrid_upsampling | Geogrid upsampling (in each direction) |
[in] | rtc_min_value_db | Minimum value for the RTC area factor. Radar data with RTC area factor below this limit are ignored. |
[in] | radar_grid_nlooks | Radar grid number of looks. This parameters determines the multilooking factor used to compute out_nlooks. |
[out] | out_nlooks | Raster to which the number of radar-grid looks associated with the geogrid will be saved |
[in] | rtc_memory_mode | Select memory mode |
[in] | interp_method | Interpolation Method |
[in] | threshold | Distance threshold for convergence |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
void isce3::geometry::facetRTC | ( | isce3::io::Raster & | dem_raster, |
isce3::io::Raster & | output_raster, | ||
const isce3::product::RadarGridParameters & | radarGrid, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | dop, | ||
const double | y0, | ||
const double | dy, | ||
const double | x0, | ||
const double | dx, | ||
const int | geogrid_length, | ||
const int | geogrid_width, | ||
const int | epsg, | ||
rtcInputRadiometry | inputRadiometry = rtcInputRadiometry::BETA_NAUGHT , |
||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR , |
||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION , |
||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN() , |
||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN() , |
||
float | radar_grid_nlooks = 1 , |
||
isce3::io::Raster * | out_geo_vertices = nullptr , |
||
isce3::io::Raster * | out_geo_grid = nullptr , |
||
isce3::io::Raster * | out_nlooks = nullptr , |
||
rtcMemoryMode | rtc_memory_mode = rtcMemoryMode::RTC_AUTO , |
||
isce3::core::dataInterpMethod | interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD , |
||
double | threshold = 1e-4 , |
||
int | num_iter = 100 , |
||
double | delta_range = 1e-4 |
||
) |
Generate radiometric terrain correction (RTC) area or area factor.
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | radarGrid | Radar Grid |
[in] | orbit | Orbit |
[in] | input_dop | Doppler LUT |
[in] | y0 | Starting northing position |
[in] | dy | Northing step size |
[in] | x0 | Starting easting position |
[in] | dx | Easting step size |
[in] | geogrid_length | Geographic length (number of pixels) in the Northing direction |
[in] | geogrid_width | Geographic width (number of pixels) in the Easting direction |
[in] | epsg | Output geographic grid EPSG |
[in] | input_radiometry | Terrain radiometry of the input raster |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_algorithm | RTC algorithm (RTC_DAVID_SMALL or RTC_AREA_PROJECTION) |
[in] | geogrid_upsampling | Geogrid upsampling (in each direction) |
[in] | rtc_min_value_db | Minimum value for the RTC area factor. Radar data with RTC area factor below this limit are ignored. |
[in] | radar_grid_nlooks | Radar grid number of looks. This parameters determines the multilooking factor used to compute out_nlooks. |
[out] | out_geo_vertices | Raster to which the radar-grid positions (range and azimuth) of the geogrid pixels vertices will be saved. |
[out] | out_geo_grid | Raster to which the radar-grid positions (range and azimuth) of the geogrid pixels center will be saved. |
[out] | out_nlooks | Raster to which the number of radar-grid looks associated with the geogrid will be saved |
[in] | rtc_memory_mode | Select memory mode |
[in] | interp_method | Interpolation Method |
[in] | threshold | Distance threshold for convergence |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
void isce3::geometry::facetRTCAreaProj | ( | isce3::io::Raster & | dem, |
isce3::io::Raster & | output_raster, | ||
const isce3::product::RadarGridParameters & | radarGrid, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | dop, | ||
const double | y0, | ||
const double | dy, | ||
const double | x0, | ||
const double | dx, | ||
const int | geogrid_length, | ||
const int | geogrid_width, | ||
const int | epsg, | ||
rtcInputRadiometry | inputRadiometry = rtcInputRadiometry::BETA_NAUGHT , |
||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR , |
||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN() , |
||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN() , |
||
float | radar_grid_nlooks = 1 , |
||
isce3::io::Raster * | out_geo_vertices = nullptr , |
||
isce3::io::Raster * | out_geo_grid = nullptr , |
||
isce3::io::Raster * | out_nlooks = nullptr , |
||
rtcMemoryMode | rtc_memory_mode = rtcMemoryMode::RTC_AUTO , |
||
isce3::core::dataInterpMethod | interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD , |
||
double | threshold = 1e-4 , |
||
int | num_iter = 100 , |
||
double | delta_range = 1e-4 |
||
) |
Generate radiometric terrain correction (RTC) area or area factor using the area projection algorithms.
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | radarGrid | Radar Grid |
[in] | orbit | Orbit |
[in] | input_dop | Doppler LUT |
[in] | y0 | Starting easting position |
[in] | dy | |
[in] | input_dop | Doppler LUT |
[in] | y0 | Starting northing position |
[in] | dy | Northing step size |
[in] | x0 | Starting easting position |
[in] | dx | Easting step size |
[in] | geogrid_length | Geographic length (number of pixels) in the Northing direction |
[in] | geogrid_width | Geographic width (number of pixels) in the Easting direction |
[in] | epsg | Output geographic grid EPSG |
[in] | input_radiometry | Terrain radiometry of the input raster |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | geogrid_upsampling | Geogrid upsampling (in each direction) |
[in] | rtc_min_value_db | Minimum value for the RTC area factor. Radar data with RTC area factor below this limit are ignored. |
[in] | radar_grid_nlooks | Radar grid number of looks. This parameters determines the multilooking factor used to compute out_nlooks. |
[out] | out_geo_vertices | Raster to which the radar-grid positions (range and azimuth) of the geogrid pixels vertices will be saved. |
[out] | out_geo_grid | Raster to which the radar-grid positions (range and azimuth) of the geogrid pixels center will be saved. |
[out] | out_nlooks | Raster to which the number of radar-grid looks associated with the geogrid will be saved |
[in] | rtc_memory_mode | Select memory mode |
[in] | interp_method | Interpolation Method |
[in] | threshold | Distance threshold for convergence |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
void isce3::geometry::facetRTCDavidSmall | ( | isce3::io::Raster & | dem_raster, |
isce3::io::Raster & | output_raster, | ||
const isce3::product::RadarGridParameters & | radarGrid, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | dop, | ||
const double | y0, | ||
const double | dy, | ||
const double | x0, | ||
const double | dx, | ||
const int | geogrid_length, | ||
const int | geogrid_width, | ||
const int | epsg, | ||
rtcInputRadiometry | inputRadiometry = rtcInputRadiometry::BETA_NAUGHT , |
||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR , |
||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN() |
||
) |
Generate radiometric terrain correction (RTC) area or area factor using the David Small algorithm.
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | radarGrid | Radar Grid |
[in] | orbit | Orbit |
[in] | input_dop | Doppler LUT |
[in] | y0 | Starting northing position |
[in] | dy | Northing step size |
[in] | x0 | Starting easting position |
[in] | dx | Easting step size |
[in] | geogrid_length | Geographic length (number of pixels) in the Northing direction |
[in] | geogrid_width | Geographic width (number of pixels) in the Easting direction |
[in] | epsg | Output geographic grid EPSG |
[in] | input_radiometry | Terrain radiometry of the input raster |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | geogrid_upsampling | Geogrid upsampling (in each direction) |
int isce3::geometry::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.
This is the elementary transformation from map geometry to radar geometry. The transformation is applicable for a single lon/lat/h coordinate (i.e., a single point target). For algorithmic details, see geometry overview.
[in] | inputLLH | Lon/Lat/Hae of target of interest |
[in] | ellipsoid | Ellipsoid object |
[in] | orbit | Orbit object |
[in] | doppler | Poly2D Doppler model |
[out] | aztime | azimuth time of inputLLH w.r.t reference epoch of the orbit |
[out] | slantRange | slant range to inputLLH |
[in] | wavelength | Radar wavelength |
[in] | startingRange | Starting slant range of reference image |
[in] | rangePixelSpacing | Slant range pixel spacing |
[in] | rwidth | Width (number of samples) of reference image |
[in] | side | Left or Right |
[in] | threshold | azimuth time convergence threshold in seconds |
[in] | maxIter | Maximum number of Newton-Raphson iterations |
[in] | deltaRange | step size used for computing derivative of doppler |
int isce3::geometry::geo2rdr | ( | const isce3::core::Vec3 & | inputLLH, |
const isce3::core::Ellipsoid & | ellipsoid, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | doppler, | ||
double & | aztime, | ||
double & | slantRange, | ||
double | wavelength, | ||
isce3::core::LookSide | side, | ||
double | threshold, | ||
int | maxIter, | ||
double | deltaRange | ||
) |
Map coordinates to radar geometry coordinates transformer.
This is the elementary transformation from map geometry to radar geometry. The transformation is applicable for a single lon/lat/h coordinate (i.e., a single point target). For algorithmic details, see geometry overview.
[in] | inputLLH | Lon/Lat/Hae of target of interest |
[in] | ellipsoid | Ellipsoid object |
[in] | orbit | Orbit object |
[in] | doppler | LUT2d Doppler model |
[out] | aztime | azimuth time of inputLLH w.r.t reference epoch of the orbit |
[out] | slantRange | slant range to inputLLH |
[in] | wavelength | Radar wavelength |
[in] | side | Left or Right |
[in] | threshold | azimuth time convergence threshold in seconds |
[in] | maxIter | Maximum number of Newton-Raphson iterations |
[in] | deltaRange | step size used for computing derivative of doppler |
isce3::geometry::BoundingBox isce3::geometry::getGeoBoundingBox | ( | const isce3::product::RadarGridParameters & | radarGrid, |
const isce3::core::Orbit & | orbit, | ||
const isce3::core::ProjectionBase * | proj, | ||
const isce3::core::LUT2d< double > & | doppler = {} , |
||
const std::vector< double > & | hgts = {isce3::core::GLOBAL_MIN_HEIGHT, isce3::core::GLOBAL_MAX_HEIGHT} , |
||
const double | margin = 0.0 , |
||
const int | pointsPerEdge = 11 , |
||
const double | threshold = 1.0e-8 , |
||
const int | numiter = 15 , |
||
bool | ignore_out_of_range_exception = false |
||
) |
Compute bounding box using min/ max altitude for quick estimates.
[in] | radarGrid | RadarGridParameters object |
[in] | orbit | Orbit object |
[in] | proj | ProjectionBase object indicating desired projection of output. |
[in] | doppler | LUT2d doppler model |
[in] | hgts | Vector of heights to use for the scene |
[in] | margin | Marging to add to estimated bounding box in decimal degrees |
[in] | pointsPerEge | Number of points to use on each edge of radar grid |
[in] | threshold | Slant range threshold for convergence |
[in] | numiter | Max number of iterations for convergence |
The output of this method is an OGREnvelope.
isce3::geometry::BoundingBox isce3::geometry::getGeoBoundingBoxHeightSearch | ( | const isce3::product::RadarGridParameters & | radarGrid, |
const isce3::core::Orbit & | orbit, | ||
const isce3::core::ProjectionBase * | proj, | ||
const isce3::core::LUT2d< double > & | doppler = {} , |
||
double | min_height = isce3::core::GLOBAL_MIN_HEIGHT , |
||
double | max_height = isce3::core::GLOBAL_MAX_HEIGHT , |
||
const double | margin = 0.0 , |
||
const int | pointsPerEdge = 11 , |
||
const double | threshold = 1.0e-8 , |
||
const int | numiter = 15 , |
||
const double | height_threshold = 100 |
||
) |
Compute bounding box with auto search within given min/ max height interval.
[in] | radarGrid | RadarGridParameters object |
[in] | orbit | Orbit object |
[in] | proj | ProjectionBase object indicating desired projection of output. |
[in] | doppler | LUT2d doppler model |
[in] | minHeight | Height lower bound |
[in] | maxHeight | Height upper bound |
[in] | margin | Marging to add to estimated bounding box in decimal degrees |
[in] | pointsPerEge | Number of points to use on each edge of radar grid |
[in] | threshold | Slant range threshold for convergence |
[in] | numiter | Max number of iterations for convergence |
[in] | height_threshold | Height threshold for convergence The output of this method is an OGREnvelope. |
isce3::geometry::Perimeter isce3::geometry::getGeoPerimeter | ( | const isce3::product::RadarGridParameters & | radarGrid, |
const isce3::core::Orbit & | orbit, | ||
const isce3::core::ProjectionBase * | proj, | ||
const isce3::core::LUT2d< double > & | doppler = {} , |
||
const DEMInterpolator & | demInterp = DEMInterpolator(0.) , |
||
const int | pointsPerEdge = 11 , |
||
const double | threshold = 1.0e-8 , |
||
const int | numiter = 15 |
||
) |
Compute the perimeter of a radar grid in map coordinates.
[in] | radarGrid | RadarGridParameters object |
[in] | orbit | Orbit object |
[in] | proj | ProjectionBase object indicating desired projection of output. |
[in] | doppler | LUT2d doppler model |
[in] | demInterp | DEM Interpolator |
[in] | pointsPerEge | Number of points to use on each edge of radar grid |
[in] | threshold | Slant range threshold for convergence |
[in] | numiter | Max number of iterations for convergence |
The outputs of this method is an OGRLinearRing. Transformer from radar geometry coordinates to map coordiantes with a DEM The sequence of walking the perimeter is always in the following order :
int isce3::geometry::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.
This is meant to be the light version of isce3::geometry::Topo and not meant to be used for processing large number of targets of interest. Note that doppler and wavelength are meant for completeness and this method can be used with both Native and Zero Doppler geometries. For details of the algorithm, see the geometry overview.
[in] | aztime | azimuth time corresponding to line of interest |
[in] | slantRange | slant range corresponding to pixel of interest |
[in] | doppler | doppler model value corresponding to line,pixel |
[in] | orbit | Orbit object |
[in] | ellipsoid | Ellipsoid object |
[in] | demInterp | DEMInterpolator object |
[out] | targetLLH | output Lon/Lat/Hae corresponding to aztime and slantRange |
[in] | wvl | imaging wavelength |
[in] | side | Left or Right. |
[in] | threshold | Distance threshold for convergence |
[in] | maxIter | Number of primary iterations |
[in] | extraIter | Number of secondary iterations |
int isce3::geometry::rdr2geo | ( | const isce3::core::Pixel & | pixel, |
const isce3::core::Basis & | TCNbasis, | ||
const isce3::core::Vec3 & | pos, | ||
const isce3::core::Vec3 & | vel, | ||
const isce3::core::Ellipsoid & | ellipsoid, | ||
const DEMInterpolator & | demInterp, | ||
isce3::core::Vec3 & | targetLLH, | ||
isce3::core::LookSide | side, | ||
double | threshold, | ||
int | maxIter, | ||
int | extraIter | ||
) |
Radar geometry coordinates to map coordinates transformer.
This is the elementary transformation from radar geometry to map geometry. The transformation is applicable for a single slant range and azimuth time (i.e., a single point target). The slant range and Doppler information are encapsulated in the Pixel object, so this function can work for both zero and native Doppler geometries. The azimuth time information is encapsulated in the TCNbasis and StateVector of the platform. For algorithmic details, see geometry overview.
[in] | pixel | Pixel object |
[in] | TCNbasis | Geocentric TCN basis corresponding to pixel |
[in] | pos/vel | position and velocity as Vec3 objects |
[in] | ellipsoid | Ellipsoid object |
[in] | demInterp | DEMInterpolator object |
[out] | targetLLH | output Lon/Lat/Hae corresponding to pixel |
[in] | side | Left or Right |
[in] | threshold | Distance threshold for convergence |
[in] | maxIter | Number of primary iterations |
[in] | extraIter | Number of secondary iterations |
int isce3::geometry::rdr2geo | ( | const isce3::core::Vec3 & | radarXYZ, |
const isce3::core::Vec3 & | axis, | ||
double | angle, | ||
double | range, | ||
const DEMInterpolator & | dem, | ||
isce3::core::Vec3 & | targetXYZ, | ||
isce3::core::LookSide | side, | ||
double | threshold, | ||
int | maxIter, | ||
int | extraIter | ||
) |
"Cone" interface to rdr2geo.
Solve for target position given radar position, range, and cone angle. The cone is described by a generating axis and the complement of the angle to that axis (e.g., angle=0 means a plane perpendicular to the axis). The vertex of the cone is at the radar position, as is the center of the range sphere.
Typically axis
is the velocity vector and angle
is the squint angle. However, with this interface you can also set axis
equal to the long axis of the antenna, in which case angle
is an azimuth angle. In this manner one can determine where the antenna boresight intersects the ground at a given range and therefore determine the Doppler from pointing.
[in] | radarXYZ | Position of antenna phase center, meters ECEF XYZ. |
[in] | axis | Cone generating axis (typically velocity), ECEF XYZ. |
[in] | angle | Complement of cone angle, radians. |
[in] | range | Range to target, meters. |
[in] | dem | Digital elevation model, meters above ellipsoid, |
[out] | targetXYZ | Target position, ECEF XYZ. |
[in] | side | Radar look direction |
[in] | threshold | Range convergence threshold, meters. |
[in] | maxIter | Maximum iterations. |
[in] | extraIter | Additional iterations. |