isce3 0.25.0
|
The isce3::geometry namespace. More...
Classes | |
class | DEMInterpolator |
class | Geo2rdr |
Transformer from map coordinates to radar geometry coordinates. More... | |
class | Geocode |
struct | RadarGridBoundingBox |
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. | |
typedef OGREnvelope | BoundingBox |
Same as GDAL's OGREnvelope structure. | |
typedef OGRTriangle | Triangle |
Same as GDAL's OGRTriangle structure. | |
Enumerations | |
enum | rtcInputTerrainRadiometry { BETA_NAUGHT = 0 , SIGMA_NAUGHT_ELLIPSOID = 1 } |
Enumeration type to indicate input terrain radiometry (for RTC) | |
enum | rtcOutputTerrainRadiometry { SIGMA_NAUGHT = 0 , GAMMA_NAUGHT = 1 } |
Enumeration type to indicate output terrain radiometry (for RTC) | |
enum | rtcAreaMode { AREA = 0 , AREA_FACTOR = 1 } |
Enumeration type to indicate RTC area mode (AREA or AREA_FACTOR) | |
enum | rtcAlgorithm { RTC_BILINEAR_DISTRIBUTION = 0 , RTC_AREA_PROJECTION = 1 } |
Enumeration type to select RTC algorithm (RTC_BILINEAR_DISTRIBUTION or RTC_AREA_PROJECTION) | |
enum | rtcAreaBetaMode { AUTO = 0 , PIXEL_AREA = 1 , PROJECTION_ANGLE = 2 } |
Enumeration type to indicate RTC area beta mode (option only available for rtcAlgorithm.RTC_AREA_PROJECTION) More... | |
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=detail::DEFAULT_TOL_HEIGHT) |
Compute the perimeter of a radar grid in map coordinates. | |
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=detail::DEFAULT_TOL_HEIGHT, bool ignore_out_of_range_exception=false) |
Compute bounding box using min/ max altitude for quick estimates. | |
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=detail::DEFAULT_TOL_HEIGHT, const double height_threshold=100) |
Compute bounding box with auto search within given min/ max height interval. | |
RadarGridBoundingBox | getRadarBoundingBox (const isce3::product::GeoGridParameters &geoGrid, const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, const float min_height=isce3::core::GLOBAL_MIN_HEIGHT, const float max_height=isce3::core::GLOBAL_MAX_HEIGHT, const isce3::core::LUT2d< double > &doppler={}, const int margin=50, const isce3::geometry::detail::Geo2RdrParams &g2r_params={}, const int geogrid_expansion_threshold=100) |
Compute bounding box of a geocoded grid within radar grid. | |
RadarGridBoundingBox | getRadarBoundingBox (const isce3::product::GeoGridParameters &geoGrid, const isce3::product::RadarGridParameters &radarGrid, const isce3::core::Orbit &orbit, isce3::geometry::DEMInterpolator &dem_interp, const isce3::core::LUT2d< double > &doppler={}, const int margin=50, const isce3::geometry::detail::Geo2RdrParams &g2r_params={}, const int geogrid_expansion_threshold=100) |
Compute bounding box of a geocoded grid within radar grid. | |
int | geo2rdr_bracket (const isce3::core::Vec3 &x, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &doppler, double &aztime, double &range, double wavelength, isce3::core::LookSide side, double tolAzTime=isce3::geometry::detail::DEFAULT_TOL_AZ_TIME, std::optional< double > timeStart=std::nullopt, std::optional< double > timeEnd=std::nullopt) |
Solve the inverse mapping problem using a derivative-free method. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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. | |
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) |
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 geodetic location of spacecraft. | |
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 geodetic location of spacecraft. | |
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 geodetic location of spacecraft. | |
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. | |
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 look vector (poitning) in ECEF. | |
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 reference ellipsoid of planet for a look vector looking from a certain spacecraft position in ECEF towards the planet. | |
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 orbit (spacecraft/antenna statevector) and at a certain relative azimuth time. | |
std::tuple< Eigen::ArrayXd, Eigen::ArrayXd > | lookIncAngFromSlantRange (const Eigen::Ref< const Eigen::ArrayXd > &slant_range, const isce3::core::Orbit &orbit, std::optional< double > az_time={}, const DEMInterpolator &dem_interp={}, const isce3::core::Ellipsoid &ellips={}) |
Overloaded vectorized version of estimating look angle (off-nadir angle) and ellipsoidal incidence angle at a desired slant range from orbit (spacecraft/antenna statevector) and at a certain relative azimuth time. | |
double | compute_mean_dem (const DEMInterpolator &dem) |
void | getGeolocationGrid (isce3::io::Raster &dem_raster, const isce3::product::RadarGridParameters &radar_grid, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &native_doppler, const isce3::core::LUT2d< double > &grid_doppler, const int epsg, isce3::core::dataInterpMethod dem_interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, const isce3::geometry::detail::Rdr2GeoParams &rdr2geo_params={}, const isce3::geometry::detail::Geo2RdrParams &geo2rdr_params={}, isce3::io::Raster *interpolated_dem_raster=nullptr, isce3::io::Raster *coordinate_x_raster=nullptr, isce3::io::Raster *coordinate_y_raster=nullptr, isce3::io::Raster *incidence_angle_raster=nullptr, isce3::io::Raster *los_unit_vector_x_raster=nullptr, isce3::io::Raster *los_unit_vector_y_raster=nullptr, isce3::io::Raster *along_track_unit_vector_x_raster=nullptr, isce3::io::Raster *along_track_unit_vector_y_raster=nullptr, isce3::io::Raster *elevation_angle_raster=nullptr, isce3::io::Raster *ground_track_velocity_raster=nullptr) |
Get geolocation grid from L1 products. | |
DEMInterpolator | DEMRasterToInterpolator (isce3::io::Raster &demRaster, const isce3::product::GeoGridParameters &geoGrid, const int demMarginInPixels=50, const isce3::core::dataInterpMethod demInterpMethod=isce3::core::BIQUINTIC_METHOD) |
returns a DEM interpolator for a geocoded grid. | |
DEMInterpolator | DEMRasterToInterpolator (isce3::io::Raster &demRaster, const isce3::product::GeoGridParameters &geoGrid, const int lineStart, const int blockLength, const int blockWidth, const int demMarginInPixels=50, const isce3::core::dataInterpMethod demInterpMethod=isce3::core::BIQUINTIC_METHOD) |
returns a DEM interpolator for a block of geocoded grid. | |
isce3::error::ErrorCode | loadDemFromProj (isce3::io::Raster &dem_raster, const double minX, const double maxX, const double minY, const double maxY, isce3::geometry::DEMInterpolator *dem_interp, isce3::core::ProjectionBase *proj=nullptr, const int dem_margin_x_in_pixels=100, const int dem_margin_y_in_pixels=100, const int dem_raster_band=1) |
Load DEM raster into a DEMInterpolator object around a given bounding box in the same or different coordinate system as the DEM raster. | |
Vec3 | getDemCoordsSameEpsg (double x, double y, const DEMInterpolator &dem_interp, isce3::core::ProjectionBase *) |
Vec3 | getDemCoordsDiffEpsg (double x, double y, const DEMInterpolator &dem_interp, isce3::core::ProjectionBase *input_proj) |
void | writeVectorDerivedCubes (const int array_pos_i, const int array_pos_j, const double native_azimuth_time, const isce3::core::Vec3 &target_llh, const isce3::core::Orbit &orbit, const isce3::core::Ellipsoid &ellipsoid, isce3::io::Raster *incidence_angle_raster, isce3::core::Matrix< float > &incidence_angle_array, isce3::io::Raster *los_unit_vector_x_raster, isce3::core::Matrix< float > &los_unit_vector_x_array, isce3::io::Raster *los_unit_vector_y_raster, isce3::core::Matrix< float > &los_unit_vector_y_array, isce3::io::Raster *along_track_unit_vector_x_raster, isce3::core::Matrix< float > &along_track_unit_vector_x_array, isce3::io::Raster *along_track_unit_vector_y_raster, isce3::core::Matrix< float > &along_track_unit_vector_y_array, isce3::io::Raster *elevation_angle_raster, isce3::core::Matrix< float > &elevation_angle_array, isce3::io::Raster *ground_track_velocity_raster, isce3::core::Matrix< double > &ground_track_velocity_array, isce3::io::Raster *local_incidence_angle_raster, isce3::core::Matrix< float > &local_incidence_angle_array, isce3::io::Raster *projection_angle_raster, isce3::core::Matrix< float > &projection_angle_array, isce3::io::Raster *simulated_radar_brightness_raster, isce3::core::Matrix< float > &simulated_radar_brightness_array, isce3::core::Vec3 *terrain_normal_unit_vec_enu=nullptr, isce3::core::LookSide *lookside=nullptr) |
Compute geometry vectors and metadata cube values for a given target point and write to the corresponding arrays. | |
void | makeRadarGridCubes (const isce3::product::RadarGridParameters &radar_grid, const isce3::product::GeoGridParameters &geogrid, const std::vector< double > &heights, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &native_doppler, const isce3::core::LUT2d< double > &grid_doppler, isce3::io::Raster *slant_range_raster=nullptr, isce3::io::Raster *azimuth_time_raster=nullptr, isce3::io::Raster *incidence_angle_raster=nullptr, isce3::io::Raster *los_unit_vector_x_raster=nullptr, isce3::io::Raster *los_unit_vector_y_raster=nullptr, isce3::io::Raster *along_track_unit_vector_x_raster=nullptr, isce3::io::Raster *along_track_unit_vector_y_raster=nullptr, isce3::io::Raster *elevation_angle_raster=nullptr, isce3::io::Raster *ground_track_velocity_raster=nullptr, const double threshold_geo2rdr=1e-8, const int numiter_geo2rdr=100, const double delta_range=1e-8, bool flag_set_output_rasters_geolocation=false) |
Make metadata radar grid cubes. | |
void | makeGeolocationGridCubes (const isce3::product::RadarGridParameters &radar_grid, const std::vector< double > &heights, const isce3::core::Orbit &orbit, const isce3::core::LUT2d< double > &native_doppler, const isce3::core::LUT2d< double > &grid_doppler, const int epsg, isce3::io::Raster *coordinate_x_raster=nullptr, isce3::io::Raster *coordinate_y_raster=nullptr, isce3::io::Raster *incidence_angle_raster=nullptr, isce3::io::Raster *los_unit_vector_x_raster=nullptr, isce3::io::Raster *los_unit_vector_y_raster=nullptr, isce3::io::Raster *along_track_unit_vector_x_raster=nullptr, isce3::io::Raster *along_track_unit_vector_y_raster=nullptr, isce3::io::Raster *elevation_angle_raster=nullptr, isce3::io::Raster *ground_track_velocity_raster=nullptr, const double threshold_geo2rdr=1e-8, const int numiter_geo2rdr=100, const double delta_range=1e-8) |
Make metadata geolocation grid cubes. | |
int | rdr2geo_bracket (double aztime, double slantRange, double doppler, const isce3::core::Orbit &orbit, const isce3::geometry::DEMInterpolator &dem, isce3::core::Vec3 &targetXYZ, double wavelength, isce3::core::LookSide side, double tolHeight=isce3::geometry::detail::DEFAULT_TOL_HEIGHT, double lookMin=0.0, double lookMax=M_PI/2) |
Transform from radar coordinates (range, azimuth) to ECEF XYZ coordinates. | |
template<typename T> | |
void | _clip_min_max (T &radar_value, float clip_min, float clip_max) |
template<typename T> | |
void | _clip_min_max (std::complex< T > &radar_value, float clip_min, float clip_max) |
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, float clip_min, float clip_max, pyre::journal::info_t &info, bool flag_complex_to_real_squared) |
void | _applyRtcMinValueDb (isce3::core::Matrix< float > &out_array, float rtc_min_value_db, rtcAreaMode rtc_area_mode, pyre::journal::info_t &info) |
void | _normalizeRtcArea (isce3::core::Matrix< float > &numerator_array, const isce3::core::Matrix< float > &denominator_array, pyre::journal::info_t &info) |
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, rtcInputTerrainRadiometry input_terrain_radiometry=rtcInputTerrainRadiometry::BETA_NAUGHT, rtcOutputTerrainRadiometry output_terrain_radiometry=rtcOutputTerrainRadiometry::GAMMA_NAUGHT, int exponent=0, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, rtcAreaBetaMode rtc_area_beta_mode=rtcAreaBetaMode::AUTO, 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 clip_min=std::numeric_limits< float >::quiet_NaN(), float clip_max=std::numeric_limits< float >::quiet_NaN(), isce3::io::Raster *out_sigma=nullptr, const isce3::core::LUT2d< double > &az_time_correction={}, const isce3::core::LUT2d< double > &slant_range_correction={}, isce3::io::Raster *input_rtc=nullptr, isce3::io::Raster *output_rtc=nullptr, isce3::core::MemoryModeBlocksY rtc_memory_mode=isce3::core::MemoryModeBlocksY::AutoBlocksY) |
Apply radiometric terrain correction (RTC) over an input raster. | |
double | computeUpsamplingFactor (const DEMInterpolator &dem_interp, const isce3::product::RadarGridParameters &radar_grid, const isce3::core::Ellipsoid &ellps) |
void | computeRtc (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, rtcInputTerrainRadiometry inputTerrainRadiometry=rtcInputTerrainRadiometry::BETA_NAUGHT, rtcOutputTerrainRadiometry outputTerrainRadiometry=rtcOutputTerrainRadiometry::GAMMA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, rtcAreaBetaMode rtc_area_beta_mode=rtcAreaBetaMode::AUTO, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), isce3::io::Raster *out_sigma=nullptr, const isce3::core::LUT2d< double > &az_time_correction={}, const isce3::core::LUT2d< double > &slant_range_correction={}, isce3::core::MemoryModeBlocksY rtc_memory_mode=isce3::core::MemoryModeBlocksY::AutoBlocksY, isce3::core::dataInterpMethod interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, double threshold=1e-8, int num_iter=100, double delta_range=1e-8, const long long min_block_size=isce3::core::DEFAULT_MIN_BLOCK_SIZE, const long long max_block_size=isce3::core::DEFAULT_MAX_BLOCK_SIZE) |
Generate radiometric terrain correction (RTC) area or area normalization factor. | |
void | computeRtc (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, rtcInputTerrainRadiometry inputTerrainRadiometry=rtcInputTerrainRadiometry::BETA_NAUGHT, rtcOutputTerrainRadiometry outputTerrainRadiometry=rtcOutputTerrainRadiometry::GAMMA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAlgorithm rtc_algorithm=rtcAlgorithm::RTC_AREA_PROJECTION, rtcAreaBetaMode rtc_area_beta_mode=rtcAreaBetaMode::AUTO, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), isce3::io::Raster *out_geo_rdr=nullptr, isce3::io::Raster *out_geo_grid=nullptr, isce3::io::Raster *out_sigma=nullptr, const isce3::core::LUT2d< double > &az_time_correction={}, const isce3::core::LUT2d< double > &slant_range_correction={}, isce3::core::MemoryModeBlocksY rtc_memory_mode=isce3::core::MemoryModeBlocksY::AutoBlocksY, isce3::core::dataInterpMethod interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, double threshold=1e-8, int num_iter=100, double delta_range=1e-8, const long long min_block_size=isce3::core::DEFAULT_MIN_BLOCK_SIZE, const long long max_block_size=isce3::core::DEFAULT_MAX_BLOCK_SIZE) |
Generate radiometric terrain correction (RTC) area or area normalization factor. | |
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 gamma_naught_area, double sigma_naught_area, double beta_naught_area, isce3::core::Matrix< float > &out_gamma_array, isce3::core::Matrix< float > &out_beta_array, isce3::core::Matrix< float > &out_sigma_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 target_to_sensor_xyz, Vec3 image_normal_xyz, rtcAreaMode rtc_area_mode, rtcAreaBetaMode rtc_area_beta_mode, double p1, double &p3, double divisor, rtcOutputTerrainRadiometry output_terrain_radiometry, double &sigma_naught_area, double &beta_naught_area) |
void | computeRtcBilinearDistribution (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 isce3::product::GeoGridParameters &geogrid, rtcInputTerrainRadiometry input_terrain_radiometry=rtcInputTerrainRadiometry::BETA_NAUGHT, rtcOutputTerrainRadiometry output_terrain_radiometry=rtcOutputTerrainRadiometry::GAMMA_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(), isce3::io::Raster *out_sigma=nullptr, double threshold=1e-8, int num_iter=100, double delta_range=1e-8, const isce3::core::LUT2d< double > &az_time_correction={}, const isce3::core::LUT2d< double > &slant_range_correction={}) |
Generate radiometric terrain correction (RTC) area or area normalization factor using the Bilinear Distribution (D. | |
void | _RunBlock (const int jmax, const int block_size, const int block_size_with_upsampling, const int block, long long &numdone, const long long progress_block, const double geogrid_upsampling, isce3::core::dataInterpMethod interp_method, isce3::io::Raster &dem_raster, isce3::io::Raster *out_geo_rdr, isce3::io::Raster *out_geo_grid, const double start, const double pixazm, const double dr, double r0, int xbound, int ybound, const isce3::product::GeoGridParameters &geogrid, 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_gamma_array, isce3::core::Matrix< float > &out_beta_array, isce3::core::Matrix< float > &out_sigma_array, const isce3::core::LUT2d< double > &az_time_correction, const isce3::core::LUT2d< double > &slant_range_correction, isce3::core::ProjectionBase *proj, rtcAreaMode rtc_area_mode, rtcAreaBetaMode rtc_area_beta_mode, rtcInputTerrainRadiometry input_terrain_radiometry, rtcOutputTerrainRadiometry output_terrain_radiometry) |
void | computeRtcAreaProj (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 isce3::product::GeoGridParameters &geogrid, rtcInputTerrainRadiometry input_terrain_radiometry=rtcInputTerrainRadiometry::BETA_NAUGHT, rtcOutputTerrainRadiometry output_terrain_radiometry=rtcOutputTerrainRadiometry::GAMMA_NAUGHT, rtcAreaMode rtc_area_mode=rtcAreaMode::AREA_FACTOR, rtcAreaBetaMode rtc_area_beta_mode=rtcAreaBetaMode::AUTO, double geogrid_upsampling=std::numeric_limits< double >::quiet_NaN(), float rtc_min_value_db=std::numeric_limits< float >::quiet_NaN(), isce3::io::Raster *out_geo_rdr=nullptr, isce3::io::Raster *out_geo_grid=nullptr, isce3::io::Raster *out_sigma=nullptr, const isce3::core::LUT2d< double > &az_time_correction={}, const isce3::core::LUT2d< double > &slant_range_correction={}, isce3::core::MemoryModeBlocksY rtc_memory_mode=isce3::core::MemoryModeBlocksY::AutoBlocksY, isce3::core::dataInterpMethod interp_method=isce3::core::dataInterpMethod::BIQUINTIC_METHOD, double threshold=1e-8, int num_iter=100, double delta_range=1e-8, const long long min_block_size=isce3::core::DEFAULT_MIN_BLOCK_SIZE, const long long max_block_size=isce3::core::DEFAULT_MAX_BLOCK_SIZE) |
Generate radiometric terrain correction (RTC) area or area normalization factor using the area projection algorithms. | |
std::string | get_input_terrain_radiometry_str (rtcInputTerrainRadiometry input_terrain_radiometry) |
Convert enum input terrain radiometry to string. | |
std::string | get_output_terrain_radiometry_str (rtcOutputTerrainRadiometry output_terrain_radiometry) |
Convert enum output terrain radiometry to string. | |
std::string | get_rtc_area_mode_str (rtcAreaMode rtc_area_mode) |
Convert enum rtc_area_mode to string. | |
std::string | get_rtc_area_beta_mode_str (rtcAreaBetaMode rtc_area_beta_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, rtcInputTerrainRadiometry input_terrain_radiometry, rtcOutputTerrainRadiometry output_terrain_radiometry, rtcAreaMode rtc_area_mode, rtcAreaBetaMode rtc_area_beta_mode, double geogrid_upsampling, float rtc_min_value_db) |
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.
Enumeration type to indicate RTC area beta mode (option only available for rtcAlgorithm.RTC_AREA_PROJECTION)
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, | ||
rtcInputTerrainRadiometry | input_terrain_radiometry = rtcInputTerrainRadiometry::BETA_NAUGHT, | ||
rtcOutputTerrainRadiometry | output_terrain_radiometry = rtcOutputTerrainRadiometry::GAMMA_NAUGHT, | ||
int | exponent = 0, | ||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR, | ||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION, | ||
rtcAreaBetaMode | rtc_area_beta_mode = rtcAreaBetaMode::AUTO, | ||
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 | clip_min = std::numeric_limits< float >::quiet_NaN(), | ||
float | clip_max = std::numeric_limits< float >::quiet_NaN(), | ||
isce3::io::Raster * | out_sigma = nullptr, | ||
const isce3::core::LUT2d< double > & | az_time_correction = {}, | ||
const isce3::core::LUT2d< double > & | slant_range_correction = {}, | ||
isce3::io::Raster * | input_rtc = nullptr, | ||
isce3::io::Raster * | output_rtc = nullptr, | ||
isce3::core::MemoryModeBlocksY | rtc_memory_mode = isce3::core::MemoryModeBlocksY::AutoBlocksY ) |
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_terrain_radiometry | Input terrain radiometry |
[in] | output_terrain_radiometry | Output terrain radiometry |
[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_BILINEAR_DISTRIBUTION or RTC_AREA_PROJECTION) |
[in] | rtc_area_beta_mode | RTC area beta mode (AUTO, PIXEL_AREA, PROJECTION_ANGLE) |
[in] | geogrid_upsampling | Geogrid upsampling |
[in] | rtc_min_value_db | Minimum value for the RTC area normalization factor. Radar data with RTC area normalization factor below this limit will be set to NaN. |
[in] | abs_cal_factor | Absolute calibration factor. |
[in] | clip_min | Clip minimum output values |
[in] | clip_max | Clip maximum output values |
[out] | out_sigma | Output sigma surface area (rtc_area_mode = AREA) or area factor (rtc_area_mode = AREA_FACTOR) raster |
[in] | az_time_correction | Azimuth additive correction, in seconds, as a function of azimuth and range |
[in] | slant_range_correction | Slant range additive correction, in meters, as a function of azimuth and range |
[in] | input_rtc | Raster containing pre-computed RTC area factor |
[out] | output_rtc | Output RTC area normalization factor |
[in] | rtc_memory_mode | Select memory mode |
double isce3::geometry::compute_mean_dem | ( | const DEMInterpolator & | dem | ) |
[in] | dem | DEMInterpolator object. |
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::computeRtc | ( | 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, | ||
rtcInputTerrainRadiometry | inputTerrainRadiometry = rtcInputTerrainRadiometry::BETA_NAUGHT, | ||
rtcOutputTerrainRadiometry | outputTerrainRadiometry = rtcOutputTerrainRadiometry::GAMMA_NAUGHT, | ||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR, | ||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION, | ||
rtcAreaBetaMode | rtc_area_beta_mode = rtcAreaBetaMode::AUTO, | ||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN(), | ||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN(), | ||
isce3::io::Raster * | out_sigma = nullptr, | ||
const isce3::core::LUT2d< double > & | az_time_correction = {}, | ||
const isce3::core::LUT2d< double > & | slant_range_correction = {}, | ||
isce3::core::MemoryModeBlocksY | rtc_memory_mode = isce3::core::MemoryModeBlocksY::AutoBlocksY, | ||
isce3::core::dataInterpMethod | interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD, | ||
double | threshold = 1e-8, | ||
int | num_iter = 100, | ||
double | delta_range = 1e-8, | ||
const long long | min_block_size = isce3::core::DEFAULT_MIN_BLOCK_SIZE, | ||
const long long | max_block_size = isce3::core::DEFAULT_MAX_BLOCK_SIZE ) |
Generate radiometric terrain correction (RTC) area or area normalization factor.
[in] | radarGrid | Radar Grid |
[in] | orbit | Orbit |
[in] | input_dop | Doppler LUT |
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | input_terrain_radiometry | Input terrain radiometry |
[in] | output_terrain_radiometry | Output terrain radiometry |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_algorithm | RTC algorithm (RTC_BILINEAR_DISTRIBUTION or RTC_AREA_PROJECTION) |
[in] | rtc_area_beta_mode | RTC area beta mode (AUTO, PIXEL_AREA, PROJECTION_ANGLE) |
[in] | geogrid_upsampling | Geogrid upsampling |
[in] | rtc_min_value_db | Minimum value for the RTC area normalization factor. Radar data with RTC area normalization factor below this limit will be set to NaN.. |
[out] | out_sigma | Output sigma surface area (rtc_area_mode = AREA) or area factor (rtc_area_mode = AREA_FACTOR) raster |
[in] | az_time_correction | Azimuth additive correction, in seconds, as a function of azimuth and range |
[in] | slant_range_correction | Slant range additive correction, in meters, as a function of azimuth and range |
[in] | rtc_memory_mode | Select memory mode |
[in] | interp_method | Interpolation Method |
[in] | threshold | Azimuth time threshold for convergence (s) |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
[in] | min_block_size | Minimum block size (per thread) |
[in] | max_block_size | Maximum block size (per thread) |
void isce3::geometry::computeRtc | ( | 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, | ||
rtcInputTerrainRadiometry | inputTerrainRadiometry = rtcInputTerrainRadiometry::BETA_NAUGHT, | ||
rtcOutputTerrainRadiometry | outputTerrainRadiometry = rtcOutputTerrainRadiometry::GAMMA_NAUGHT, | ||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR, | ||
rtcAlgorithm | rtc_algorithm = rtcAlgorithm::RTC_AREA_PROJECTION, | ||
rtcAreaBetaMode | rtc_area_beta_mode = rtcAreaBetaMode::AUTO, | ||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN(), | ||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN(), | ||
isce3::io::Raster * | out_geo_rdr = nullptr, | ||
isce3::io::Raster * | out_geo_grid = nullptr, | ||
isce3::io::Raster * | out_sigma = nullptr, | ||
const isce3::core::LUT2d< double > & | az_time_correction = {}, | ||
const isce3::core::LUT2d< double > & | slant_range_correction = {}, | ||
isce3::core::MemoryModeBlocksY | rtc_memory_mode = isce3::core::MemoryModeBlocksY::AutoBlocksY, | ||
isce3::core::dataInterpMethod | interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD, | ||
double | threshold = 1e-8, | ||
int | num_iter = 100, | ||
double | delta_range = 1e-8, | ||
const long long | min_block_size = isce3::core::DEFAULT_MIN_BLOCK_SIZE, | ||
const long long | max_block_size = isce3::core::DEFAULT_MAX_BLOCK_SIZE ) |
Generate radiometric terrain correction (RTC) area or area normalization 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_terrain_radiometry | Input terrain radiometry |
[in] | output_terrain_radiometry | Output terrain radiometry |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_algorithm | RTC algorithm (RTC_BILINEAR_DISTRIBUTION or RTC_AREA_PROJECTION) |
[in] | rtc_area_beta_mode | RTC area beta mode (AUTO, PIXEL_AREA, PROJECTION_ANGLE) |
[in] | geogrid_upsampling | Geogrid upsampling |
[in] | rtc_min_value_db | Minimum value for the RTC area normalization factor. Radar data with RTC area normalization factor below this limit will be set to NaN |
[out] | out_geo_rdr | 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_sigma | Output sigma surface area (rtc_area_mode = AREA) or area factor (rtc_area_mode = AREA_FACTOR) raster |
[in] | az_time_correction | Azimuth additive correction, in seconds, as a function of azimuth and range |
[in] | slant_range_correction | Slant range additive correction, in meters, as a function of azimuth and range |
[in] | rtc_memory_mode | Select memory mode |
[in] | interp_method | Interpolation Method |
[in] | threshold | Azimuth time threshold for convergence (s) |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
[in] | min_block_size | Minimum block size (per thread) |
[in] | max_block_size | Maximum block size (per thread) |
void isce3::geometry::computeRtcAreaProj | ( | 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 isce3::product::GeoGridParameters & | geogrid, | ||
rtcInputTerrainRadiometry | input_terrain_radiometry = rtcInputTerrainRadiometry::BETA_NAUGHT, | ||
rtcOutputTerrainRadiometry | output_terrain_radiometry = rtcOutputTerrainRadiometry::GAMMA_NAUGHT, | ||
rtcAreaMode | rtc_area_mode = rtcAreaMode::AREA_FACTOR, | ||
rtcAreaBetaMode | rtc_area_beta_mode = rtcAreaBetaMode::AUTO, | ||
double | geogrid_upsampling = std::numeric_limits< double >::quiet_NaN(), | ||
float | rtc_min_value_db = std::numeric_limits< float >::quiet_NaN(), | ||
isce3::io::Raster * | out_geo_rdr = nullptr, | ||
isce3::io::Raster * | out_geo_grid = nullptr, | ||
isce3::io::Raster * | out_sigma = nullptr, | ||
const isce3::core::LUT2d< double > & | az_time_correction = {}, | ||
const isce3::core::LUT2d< double > & | slant_range_correction = {}, | ||
isce3::core::MemoryModeBlocksY | rtc_memory_mode = isce3::core::MemoryModeBlocksY::AutoBlocksY, | ||
isce3::core::dataInterpMethod | interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD, | ||
double | threshold = 1e-8, | ||
int | num_iter = 100, | ||
double | delta_range = 1e-8, | ||
const long long | min_block_size = isce3::core::DEFAULT_MIN_BLOCK_SIZE, | ||
const long long | max_block_size = isce3::core::DEFAULT_MAX_BLOCK_SIZE ) |
Generate radiometric terrain correction (RTC) area or area normalization 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] | geogrid | Geogrid parameters |
[in] | input_terrain_radiometry | Input terrain radiometry |
[in] | output_terrain_radiometry | Output terrain radiometry |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | rtc_area_beta_mode | RTC area beta mode (AUTO, PIXEL_AREA, PROJECTION_ANGLE) |
[in] | geogrid_upsampling | Geogrid upsampling |
[in] | rtc_min_value_db | Minimum value for the RTC area normalization factor. Radar data with RTC area normalization factor below this limit will be set to NaN.. |
[out] | out_geo_rdr | 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_sigma | Output sigma surface area (rtc_area_mode = AREA) or area factor (rtc_area_mode = AREA_FACTOR) raster |
[in] | az_time_correction | Azimuth additive correction, in seconds, as a function of azimuth and range |
[in] | slant_range_correction | Slant range additive correction, in meters, as a function of azimuth and range |
[in] | rtc_memory_mode | Select memory mode |
[in] | interp_method | Interpolation Method |
[in] | threshold | Azimuth time threshold for convergence (s) |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
[in] | min_block_size | Minimum block size (per thread) |
[in] | max_block_size | Maximum block size (per thread) |
void isce3::geometry::computeRtcBilinearDistribution | ( | 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 isce3::product::GeoGridParameters & | geogrid, | ||
rtcInputTerrainRadiometry | input_terrain_radiometry = rtcInputTerrainRadiometry::BETA_NAUGHT, | ||
rtcOutputTerrainRadiometry | output_terrain_radiometry = rtcOutputTerrainRadiometry::GAMMA_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(), | ||
isce3::io::Raster * | out_sigma = nullptr, | ||
double | threshold = 1e-8, | ||
int | num_iter = 100, | ||
double | delta_range = 1e-8, | ||
const isce3::core::LUT2d< double > & | az_time_correction = {}, | ||
const isce3::core::LUT2d< double > & | slant_range_correction = {} ) |
Generate radiometric terrain correction (RTC) area or area normalization factor using the Bilinear Distribution (D.
Small) algorithm [small2011].
[in] | dem_raster | Input DEM raster |
[out] | output_raster | Output raster |
[in] | radarGrid | Radar Grid |
[in] | orbit | Orbit |
[in] | input_dop | Doppler LUT |
[in] | geogrid | Geogrid parameters |
[in] | epsg | Output geographic grid EPSG |
[in] | input_terrain_radiometry | Input terrain radiometry |
[in] | output_terrain_radiometry | Output terrain radiometry |
[in] | rtc_area_mode | RTC area mode (AREA or AREA_FACTOR) |
[in] | geogrid_upsampling | Geogrid upsampling |
[in] | rtc_min_value_db | Minimum value for the RTC area normalization factor. Radar data with RTC area normalization factor below this limit will be set to NaN |
[out] | out_sigma | Output sigma surface area (rtc_area_mode = AREA) or area factor (rtc_area_mode = AREA_FACTOR) raster |
[in] | threshold | Azimuth time threshold for convergence (s) |
[in] | num_iter | Maximum number of Newton-Raphson iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
[in] | az_time_correction | Azimuth additive correction, in seconds, as a function of azimuth and range |
[in] | slant_range_correction | Slant range additive correction, in meters, as a function of azimuth and range |
DEMInterpolator isce3::geometry::DEMRasterToInterpolator | ( | isce3::io::Raster & | demRaster, |
const isce3::product::GeoGridParameters & | geoGrid, | ||
const int | demMarginInPixels = 50, | ||
const isce3::core::dataInterpMethod | demInterpMethod = isce3::core::BIQUINTIC_METHOD ) |
returns a DEM interpolator for a geocoded grid.
The geocoded grid and the input raster of the DEM can be in different or same projection systems
[in] | demRaster | a raster of the DEM |
[in] | geoGrid | parameters of the geocoded grid |
[in] | demMarginInPixels | DEM extra margin in pixels |
[in] | demInterpMethod | DEM interpolation method |
isce3::geometry::DEMInterpolator isce3::geometry::DEMRasterToInterpolator | ( | isce3::io::Raster & | demRaster, |
const isce3::product::GeoGridParameters & | geoGrid, | ||
const int | lineStart, | ||
const int | blockLength, | ||
const int | blockWidth, | ||
const int | demMarginInPixels = 50, | ||
const isce3::core::dataInterpMethod | demInterpMethod = isce3::core::BIQUINTIC_METHOD ) |
returns a DEM interpolator for a block of geocoded grid.
The geocoded grid and the inpit raster of the DEM can be in different or same projection systems
[in] | demRaster | a raster of the DEM |
[in] | geoGrid | parameters of the geocoded grid |
[in] | lineStart | start line of the block of interest in the geocoded grid |
[in] | blockLength | length of the block of interest in the geocoded grid |
[in] | blockWidth | width of the block of interest in the geocoded grid |
[in] | demMarginInPixels | DEM extra margin in pixels |
[in] | demInterpMethod | DEM interpolation method |
Vec3 isce3::geometry::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 geodetic location of spacecraft.
[in] | lon | : geodetic longitude in radians |
[in] | lat | : geodetic latitude in radians |
[in] | vel | : velocity vector or its unit vector in ECEF (x,y,z) |
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 |
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_bracket | ( | const isce3::core::Vec3 & | x, |
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | doppler, | ||
double & | aztime, | ||
double & | range, | ||
double | wavelength, | ||
isce3::core::LookSide | side, | ||
double | tolAzTime = isce3::geometry::detail::DEFAULT_TOL_AZ_TIME, | ||
std::optional< double > | timeStart = std::nullopt, | ||
std::optional< double > | timeEnd = std::nullopt ) |
Solve the inverse mapping problem using a derivative-free method.
[in] | x | Target position, XYZ in m |
[in] | orbit | Platform position vs time |
[in] | doppler | Geometric Doppler centroid as function of (time, range), Hz |
[out] | aztime | Time when Doppler centroid crosses target, s |
[out] | range | Distance to target at aztime, m |
[in] | wavelength | Radar wavelength used to scale Doppler, m |
[in] | side | Radar look direction, Left or Right |
[in] | tolAzTime | Azimuth convergence tolerance, s |
[in] | timeStart | Start of search interval, s Defaults to max of orbit and Doppler LUT start time |
[in] | timeEnd | End of search interval, s Defaults to min of orbit and Doppler LUT end time |
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 = detail::DEFAULT_TOL_HEIGHT, | ||
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 | Height threshold (m) for rdr2geo 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 = detail::DEFAULT_TOL_HEIGHT, | ||
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 | Margin to add to estimated bounding box in decimal degrees |
[in] | pointsPerEge | Number of points to use on each edge of radar grid |
[in] | threshold | Height threshold (m) for rdr2geo convergence |
[in] | height_threshold | Height threshold for convergence The output of this method is an OGREnvelope. |
void isce3::geometry::getGeolocationGrid | ( | isce3::io::Raster & | dem_raster, |
const isce3::product::RadarGridParameters & | radar_grid, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | native_doppler, | ||
const isce3::core::LUT2d< double > & | grid_doppler, | ||
const int | epsg, | ||
isce3::core::dataInterpMethod | dem_interp_method = isce3::core::dataInterpMethod::BIQUINTIC_METHOD, | ||
const isce3::geometry::detail::Rdr2GeoParams & | rdr2geo_params = {}, | ||
const isce3::geometry::detail::Geo2RdrParams & | geo2rdr_params = {}, | ||
isce3::io::Raster * | interpolated_dem_raster = nullptr, | ||
isce3::io::Raster * | coordinate_x_raster = nullptr, | ||
isce3::io::Raster * | coordinate_y_raster = nullptr, | ||
isce3::io::Raster * | incidence_angle_raster = nullptr, | ||
isce3::io::Raster * | los_unit_vector_x_raster = nullptr, | ||
isce3::io::Raster * | los_unit_vector_y_raster = nullptr, | ||
isce3::io::Raster * | along_track_unit_vector_x_raster = nullptr, | ||
isce3::io::Raster * | along_track_unit_vector_y_raster = nullptr, | ||
isce3::io::Raster * | elevation_angle_raster = nullptr, | ||
isce3::io::Raster * | ground_track_velocity_raster = nullptr ) |
Get geolocation grid from L1 products.
The target-to-sensor line-of-sight (LOS) and along-track unit vectors are referenced to ENU coordinates computed wrt targets.
[in] | dem_raster | DEM raster |
[in] | radar_grid | Radar grid |
[in] | orbit | Reference orbit |
[in] | native_doppler | Native image Doppler |
[in] | grid_doppler | Grid Doppler |
[in] | epsg | Output geolocation EPSG |
[in] | dem_interp_method | DEM interpolation method |
[in] | rdr2geo_params | Geo2rdr parameters |
[in] | geo2rdr_params | Rdr2geo parameters |
[in] | threshold_geo2rdr | Range threshold for geo2rdr |
[in] | numiter_geo2rdr | Geo2rdr maximum number of iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
[out] | interpolated_dem_raster | Interpolated DEM raster |
[out] | coordinate_x_raster | Coordinate-X raster |
[out] | coordinate_y_raster | Coordinate-Y raster |
[out] | incidence_angle_raster | Incidence angle (in degrees wrt ellipsoid normal at target) cube raster |
[out] | los_unit_vector_x_raster | LOS (target-to-sensor) unit vector X cube raster |
[out] | los_unit_vector_y_raster | LOS (target-to-sensor) unit vector Y cube raster |
[out] | along_track_unit_vector_x_raster | Along-track unit vector X cube raster |
[out] | along_track_unit_vector_y_raster | Along-track unit vector Y cube raster |
[out] | elevation_angle_raster | Elevation angle (in degrees wrt geodedic nadir) cube raster |
[out] | ground_track_velocity_raster | Ground-track velocity raster |
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 = detail::DEFAULT_TOL_HEIGHT ) |
Compute the perimeter of a radar grid in map coordinates.
[in] | radarGrid | RadarGridParameters object |
[in] | orbit | Orbit object |
[in] | proj | ProjectionBase object indicating desired map 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 | Height threshold (m) for rdr2geo convergence |
The output of this method is an OGRLinearRing.
The sequence of walking the perimeter is guaranteed to have counter-clockwise order on a map projection, which is required by the ISO 19125 specification for Simple Features used by OGR/GDAL. (Note that ESRI shapefiles use the opposite convention.)
A consequence of the consistently CCW map ordering is that the radar coordinate traversal order depends on the look direction. If we label the corners as follows:
isce3::geometry::RadarGridBoundingBox isce3::geometry::getRadarBoundingBox | ( | const isce3::product::GeoGridParameters & | geoGrid, |
const isce3::product::RadarGridParameters & | radarGrid, | ||
const isce3::core::Orbit & | orbit, | ||
const float | min_height = isce3::core::GLOBAL_MIN_HEIGHT, | ||
const float | max_height = isce3::core::GLOBAL_MAX_HEIGHT, | ||
const isce3::core::LUT2d< double > & | doppler = {}, | ||
const int | margin = 50, | ||
const isce3::geometry::detail::Geo2RdrParams & | g2r_params = {}, | ||
const int | geogrid_expansion_threshold = 100 ) |
Compute bounding box of a geocoded grid within radar grid.
The output of this function is a RadarGridBoundingBox object that defines the bounding box that is extended by the input margin in all directions. The output object contains index of: 1) the first line in azimuth 2) the last line in azimuth 3) first pixel in range 4) the last pixel in range An exception is raised if: 1) any corner fails to converge 2) computed bounding box min index >= max index along an axis 3) computed bounding box index is out of bounds
[in] | geoGrid | GeoGridParameters object |
[in] | radarGrid | RadarGridParameters object |
[in] | orbit | Orbit object |
[in] | min_height | minimum height value to be used with x and y corners for geo2rdr computations |
[in] | max_height | maximum height value to be used with x and y corners for geo2rdr computations |
[in] | doppler | LUT2d doppler model of the radar image grid |
[in] | margin | Margin to add to estimated bounding box in pixels (default 50) |
[in] | g2r_params | Stucture containing the following geo2rdr params: azimuth time threshold for convergence (default: 1e-8) maximum number of iterations for convergence (default: 50) step size used for computing derivative of doppler (default: 10.0) |
[in] | geogrid_expansion_threshold | Maximum number of iterations to outwardly expand each geogrid corner to search for geo2rdr convergence. An outward extension shifts a corner by (geogrid.dx, geogrid.dy) (default: 100) |
isce3::geometry::RadarGridBoundingBox isce3::geometry::getRadarBoundingBox | ( | const isce3::product::GeoGridParameters & | geoGrid, |
const isce3::product::RadarGridParameters & | radarGrid, | ||
const isce3::core::Orbit & | orbit, | ||
isce3::geometry::DEMInterpolator & | dem_interp, | ||
const isce3::core::LUT2d< double > & | doppler = {}, | ||
const int | margin = 50, | ||
const isce3::geometry::detail::Geo2RdrParams & | g2r_params = {}, | ||
const int | geogrid_expansion_threshold = 100 ) |
Compute bounding box of a geocoded grid within radar grid.
Use this for one off radar grid bounding box computations. Repeated use of this function with the same DEM wastefully recomputes min and max DEM height.
The output of this function is a RadarGridBoundingBox object that defines the bounding box that is extended by the input margin in all directions. The output object contains index of: 1) the first line in azimuth 2) the last line in azimuth 3) first pixel in range 4) the last pixel in range An exception is raised if: 1) any corner fails to converge 2) computed bounding box min index >= max index along an axis 3) computed bounding box index is out of bounds
[in] | geoGrid | GeoGridParameters object |
[in] | radarGrid | RadarGridParameters object |
[in] | orbit | Orbit object |
[in] | demInterp | DEMInterpolator object used to find min and max height values to used in geo2rdr computation. |
[in] | doppler | LUT2d doppler model of the radar image grid |
[in] | margin | Margin to add to estimated bounding box in pixels (default 50) |
[in] | g2r_params | Stucture containing the following geo2rdr params: azimuth time threshold for convergence (default: 1e-8) maximum number of iterations for convergence (default: 50) step size used for computing derivative of doppler (default: 10.0) |
[in] | geogrid_expansion_threshold | Maximum number of iterations to outwardly expand each geogrid corner to search for geo2rdr convergence. An outward extension shifts a corner by (geogrid.dx, geogrid.dy) (default: 100) |
double isce3::geometry::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.
[in] | lon | : geodetic longitude in radians |
[in] | lat | : geodetic latitude in radians |
[in] | vel | : velocity vector or its unit vector in ECEF (x,y,z) |
isce3::error::ErrorCode isce3::geometry::loadDemFromProj | ( | isce3::io::Raster & | dem_raster, |
const double | minX, | ||
const double | maxX, | ||
const double | minY, | ||
const double | maxY, | ||
isce3::geometry::DEMInterpolator * | dem_interp, | ||
isce3::core::ProjectionBase * | proj = nullptr, | ||
const int | dem_margin_x_in_pixels = 100, | ||
const int | dem_margin_y_in_pixels = 100, | ||
const int | dem_raster_band = 1 ) |
Load DEM raster into a DEMInterpolator object around a given bounding box in the same or different coordinate system as the DEM raster.
[in] | dem_raster | DEM raster |
[in] | x0 | Easting/longitude of western edge of bounding box, If the DEM is in geographic coordinates and the x0 coordinate is not from the polar stereo system EPSG 3031 or EPSG 3413, this point represents the minimum X coordinate value. In this case, the maximum longitude span that this function can handle is 180 degrees (when the DEM is in geographic coordinates and proj is in polar stereo) |
[in] | xf | Easting/longitude of eastern edge of bounding box If the DEM is in geographic coordinates and the xf coordinate is not from the polar stereo system EPSG 3031 or EPSG 3413, this point represents the maximum X coordinate value. In this case, the maximum longitude span that this function can handle is 180 degrees (when the DEM is in geographic coordinates and proj is in polar stereo) |
[in] | minY | Minimum Y/northing position |
[in] | maxY | Maximum Y/northing position |
[out] | dem_interp_block | DEM interpolation object |
[in] | proj | Projection object (nullptr to use same DEM projection) |
[in] | dem_margin_x_in_pixels | DEM X/easting margin in pixels |
[in] | dem_margin_y_in_pixels | DEM Y/northing margin in pixels |
[in] | dem_raster_band | DEM raster band (starting from 1) |
std::tuple< Eigen::ArrayXd, Eigen::ArrayXd > isce3::geometry::lookIncAngFromSlantRange | ( | const Eigen::Ref< const Eigen::ArrayXd > & | slant_range, |
const isce3::core::Orbit & | orbit, | ||
std::optional< double > | az_time = {}, | ||
const DEMInterpolator & | dem_interp = {}, | ||
const isce3::core::Ellipsoid & | ellips = {} ) |
Overloaded vectorized version of estimating look angle (off-nadir angle) and ellipsoidal incidence angle at a desired slant range from orbit (spacecraft/antenna statevector) and at a certain relative azimuth time.
Note that this is an approximate closed-form solution where an approximate sphere (like SCH coordinate) to Ellipsoid is formed at a certain spacecraft location and along its heading. Finally,the look angle is calculated in a closed-form expression via "Law of Cosines". See Law of Cosines The respective ellipsoidal incidence angle is formed in that approximate sphere by using a fixed height w.r.t its ellipsoid for all look angles via "Law of Sines". See Law of Sines.
[in] | slant_range | true slant range in meters from antenna phase center (or spacecraft position) to the ground. |
[in] | orbit | isce3 orbit object. |
[in] | az_time | (optional): relative azimuth time in seconds w.r.t reference epoch time of orbit object. If not speficied or set to {} or std::nullopt, the mid time of orbit will be used as azimuth time. |
[in] | dem_interp | (optional) : DEMInterpolator object wrt to reference ellipsoid. Default is global 0.0 (m) height. |
[in] | ellips | (optional) : Ellipsoid object. Default is WGS84 reference ellipsoid. |
std::tuple< double, double > isce3::geometry::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 orbit (spacecraft/antenna statevector) and at a certain relative azimuth time.
Note that this is an approximate closed-form solution where an approximate sphere (like SCH coordinate) to Ellipsoid is formed at a certain spacecraft location and along its heading. Finally,the look angle is calculated in a closed-form expression via "Law of Cosines". See Law of Cosines The respective ellipsoidal incidence angle is formed in that approximate sphere by using a fixed height w.r.t its ellipsoid for all look angles via "Law of Sines". See Law of Sines.
[in] | slant_range | true slant range in meters from antenna phase center (or spacecraft position) to the ground. |
[in] | orbit | isce3 orbit object. |
[in] | az_time | (optional): relative azimuth time in seconds w.r.t reference epoch time of orbit object. If not speficied or set to {} or std::nullopt, the mid time of orbit will be used as azimuth time. |
[in] | dem_interp | (optional) : DEMInterpolator object wrt to reference ellipsoid. Default is global 0.0 (m) height. |
[in] | ellips | (optional) : Ellipsoid object. Default is WGS84 reference ellipsoid. |
RuntimeError |
void isce3::geometry::makeGeolocationGridCubes | ( | const isce3::product::RadarGridParameters & | radar_grid, |
const std::vector< double > & | heights, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | native_doppler, | ||
const isce3::core::LUT2d< double > & | grid_doppler, | ||
const int | epsg, | ||
isce3::io::Raster * | coordinate_x_raster = nullptr, | ||
isce3::io::Raster * | coordinate_y_raster = nullptr, | ||
isce3::io::Raster * | incidence_angle_raster = nullptr, | ||
isce3::io::Raster * | los_unit_vector_x_raster = nullptr, | ||
isce3::io::Raster * | los_unit_vector_y_raster = nullptr, | ||
isce3::io::Raster * | along_track_unit_vector_x_raster = nullptr, | ||
isce3::io::Raster * | along_track_unit_vector_y_raster = nullptr, | ||
isce3::io::Raster * | elevation_angle_raster = nullptr, | ||
isce3::io::Raster * | ground_track_velocity_raster = nullptr, | ||
const double | threshold_geo2rdr = 1e-8, | ||
const int | numiter_geo2rdr = 100, | ||
const double | delta_range = 1e-8 ) |
Make metadata geolocation grid cubes.
Metadata geolocation grid cubes describe the radar geometry over a three-dimensional grid, defined by a reference radar grid and a vector of heights.
The representation as cubes, rather than two-dimensional rasters, is intended to reduce the amount of disk space required to store radar geometry values within NISAR L1 products.
This is possible because each cube contains slow-varying values in space, that can be described by a low-resolution three-dimensional grid with sufficient accuracy.
These values, however, are usually required at the terrain height, often characterized by a fast-varying surface representing the local topography. A high-resolution DEM can then be used to interpolate the metadata cubes and generate high-resolution maps of the corresponding radar geometry variable.
Each output layer is saved onto the first band of its associated raster file.
The line-of-sight (LOS) and along-track unit vectors are referenced to ENU coordinates computed wrt targets.
[in] | radar_grid | Cube radar grid |
[in] | heights | Cube heights |
[in] | orbit | Orbit |
[in] | native_doppler | Native image Doppler |
[in] | grid_doppler | Grid Doppler. |
[in] | epsg | Output geolocation EPSG |
[out] | coordinate_x_raster | Geolocation coordinate X raster |
[out] | coordinate_y_raster | Geolocation coordinage Y raster |
[out] | incidence_angle_raster | Incidence angle (in degrees wrt ellipsoid normal at target) cube raster |
[out] | los_unit_vector_x_raster | LOS (target-to-sensor) unit vector X cube raster |
[out] | los_unit_vector_y_raster | LOS (target-to-sensor) unit vector Y cube raster |
[out] | along_track_unit_vector_x_raster | Along-track unit vector X cube raster |
[out] | along_track_unit_vector_y_raster | Along-track unit vector Y cube raster |
[out] | elevation_angle_raster | Elevation angle (in degrees wrt geodedic nadir) cube raster |
[out] | ground_track_velocity_raster | Ground-track velocity raster |
[in] | threshold_geo2rdr | Azimuth time threshold for geo2rdr |
[in] | numiter_geo2rdr | Geo2rdr maximum number of iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
void isce3::geometry::makeRadarGridCubes | ( | const isce3::product::RadarGridParameters & | radar_grid, |
const isce3::product::GeoGridParameters & | geogrid, | ||
const std::vector< double > & | heights, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::LUT2d< double > & | native_doppler, | ||
const isce3::core::LUT2d< double > & | grid_doppler, | ||
isce3::io::Raster * | slant_range_raster = nullptr, | ||
isce3::io::Raster * | azimuth_time_raster = nullptr, | ||
isce3::io::Raster * | incidence_angle_raster = nullptr, | ||
isce3::io::Raster * | los_unit_vector_x_raster = nullptr, | ||
isce3::io::Raster * | los_unit_vector_y_raster = nullptr, | ||
isce3::io::Raster * | along_track_unit_vector_x_raster = nullptr, | ||
isce3::io::Raster * | along_track_unit_vector_y_raster = nullptr, | ||
isce3::io::Raster * | elevation_angle_raster = nullptr, | ||
isce3::io::Raster * | ground_track_velocity_raster = nullptr, | ||
const double | threshold_geo2rdr = 1e-8, | ||
const int | numiter_geo2rdr = 100, | ||
const double | delta_range = 1e-8, | ||
bool | flag_set_output_rasters_geolocation = false ) |
Make metadata radar grid cubes.
Metadata radar grid cubes describe the radar geometry over a three-dimensional grid, defined by a reference geogrid and a vector of heights.
The representation as cubes, rather than two-dimensional rasters, is intended to reduce the amount of disk space required to store radar geometry values within NISAR L2 products.
This is possible because each cube contains slow-varying values in space, that can be described by a low-resolution three-dimensional grid with sufficient accuracy.
These values, however, are usually required at the terrain height, often characterized by a fast-varying surface representing the local topography. A high-resolution DEM can then be used to interpolate the metadata cubes and generate high-resolution maps of the corresponding radar geometry variable.
Each output layer is saved onto the first band of its associated raster file.
The line-of-sight (LOS) and along-track unit vectors are referenced to ENU coordinates computed wrt targets.
[in] | radar_grid | Radar grid |
[in] | geogrid | Geogrid to be used as reference for output cubes |
[in] | heights | Cube heights [m] |
[in] | orbit | Reference orbit |
[in] | native_doppler | Native image Doppler |
[in] | grid_doppler | Grid Doppler |
[out] | slant_range_raster | Slant-range (in meters) cube raster |
[out] | azimuth_time_raster | Azimuth time (in seconds) cube raster |
[out] | incidence_angle_raster | Incidence angle (in degrees wrt ellipsoid normal at target) cube raster |
[out] | los_unit_vector_x_raster | LOS (target-to-sensor) unit vector X cube raster |
[out] | los_unit_vector_y_raster | LOS (target-to-sensor) unit vector Y cube raster |
[out] | along_track_unit_vector_x_raster | Along-track unit vector X cube raster |
[out] | along_track_unit_vector_y_raster | Along-track unit vector Y cube raster |
[out] | elevation_angle_raster | Elevation angle (in degrees wrt geodedic nadir) cube raster |
[out] | ground_track_velocity_raster | Ground-track velocity raster |
[in] | threshold_geo2rdr | Azimuth time threshold for geo2rdr |
[in] | numiter_geo2rdr | Geo2rdr maximum number of iterations |
[in] | delta_range | Step size used for computing derivative of doppler |
[in] | flag_set_output_rasters_geolocation | Set output rasters' geotransform and spatial reference |
Vec3 isce3::geometry::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 geodetic location of spacecraft.
[in] | lon | : geodetic longitude in radians |
[in] | lat | : geodetic latitude in radians |
[in] | vel | : velocity vector or its unit vector in ECEF (x,y,z) |
Vec3 isce3::geometry::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 geodetic location of spacecraft.
[in] | lon | : geodetic longitude in radians |
[in] | lat | : geodetic latitude in radians |
[in] | vel | : velocity vector or its unit vector in ECEF (x,y,z) |
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. |
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_bracket | ( | double | aztime, |
double | slantRange, | ||
double | doppler, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::geometry::DEMInterpolator & | dem, | ||
isce3::core::Vec3 & | targetXYZ, | ||
double | wavelength, | ||
isce3::core::LookSide | side, | ||
double | tolHeight = isce3::geometry::detail::DEFAULT_TOL_HEIGHT, | ||
double | lookMin = 0.0, | ||
double | lookMax = M_PI/2 ) |
Transform from radar coordinates (range, azimuth) to ECEF XYZ coordinates.
[in] | aztime | Target azimuth time w.r.t. orbit reference epoch (s) |
[in] | slantRange | Target slant range (m) |
[in] | doppler | Doppler centroid at target position (Hz) |
[in] | orbit | Platform orbit |
[in] | dem | DEM sampling interface |
[out] | targetXYZ | Output target ECEF XYZ position (m) |
[in] | wavelength | Radar wavelength (wrt requested Doppler) (m) |
[in] | side | Radar look side |
[in] | tolHeight | Allowable height error of solution (m) |
[in] | lookMin | Smallest possible pseudo-look angle (rad) |
[in] | lookMax | Largest possible pseudo-look angle (rad) |
Note: Usually the look angle is defined as the angle between the line of sight vector and the nadir vector. Here the pseudo-look angle is defined in a similar way, except it's with respect to the projection of the nadir vector into a plane perpendicular to the velocity. For simplicity, we use the geocentric nadir definition.
double isce3::geometry::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 look vector (poitning) in ECEF.
[in] | pos | : a non-zero x,y,z positions of antenna/platform in (m,m,m) |
[in] | lkvec | : looking/pointing unit vector in ECEF towards planet from Antenna/platform |
[in] | ellips | (optional) : Ellipsoid object. Default is WGS84 reference ellipsoid |
InvalidArgument,RuntimeError | See section 6.1 of reference [ReeTechDesDoc] |
std::pair< int, double > isce3::geometry::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 reference ellipsoid of planet for a look vector looking from a certain spacecraft position in ECEF towards the planet.
[out] | sr | : slant range (m) to the point on the planet at a certain height. |
[out] | tg_pos | : ECEF Position (x,y,z) of a point target on the planet at certain height |
[out] | llh | : geodetic Lon/lat/height Position of a point on the planet at certain height. |
[in] | sc_pos | : Spacecraft position in ECEF (x,y,z) all in (m) |
[in] | lkvec | : look unit vector in ECEF, looking from spacecraft towards the planet. |
[in] | dem_interp | (optional) : DEMInterpolator object wrt to reference ellipsoid. Default is global 0.0 (m) height. |
[in] | hgt_err | (optional) : Max error in height estimation (m) between desired input height and final output height. |
[in] | num_iter | (optional) : Max number of iterations in height estimation |
[in] | ellips | (optional) : Ellipsoid object. Default is WGS84 reference ellipsoid. |
[in] | initial_height | (optional): initial height wrt ellipsoid used in the iterative process. If not specified or set to {} or std::nullopt, stats of DEM raster is computed if not already. |
InvalidArgument,RuntimeError | See section 6.1 of reference [ReeTechDesDoc] |
void isce3::geometry::writeVectorDerivedCubes | ( | const int | array_pos_i, |
const int | array_pos_j, | ||
const double | native_azimuth_time, | ||
const isce3::core::Vec3 & | target_llh, | ||
const isce3::core::Orbit & | orbit, | ||
const isce3::core::Ellipsoid & | ellipsoid, | ||
isce3::io::Raster * | incidence_angle_raster, | ||
isce3::core::Matrix< float > & | incidence_angle_array, | ||
isce3::io::Raster * | los_unit_vector_x_raster, | ||
isce3::core::Matrix< float > & | los_unit_vector_x_array, | ||
isce3::io::Raster * | los_unit_vector_y_raster, | ||
isce3::core::Matrix< float > & | los_unit_vector_y_array, | ||
isce3::io::Raster * | along_track_unit_vector_x_raster, | ||
isce3::core::Matrix< float > & | along_track_unit_vector_x_array, | ||
isce3::io::Raster * | along_track_unit_vector_y_raster, | ||
isce3::core::Matrix< float > & | along_track_unit_vector_y_array, | ||
isce3::io::Raster * | elevation_angle_raster, | ||
isce3::core::Matrix< float > & | elevation_angle_array, | ||
isce3::io::Raster * | ground_track_velocity_raster, | ||
isce3::core::Matrix< double > & | ground_track_velocity_array, | ||
isce3::io::Raster * | local_incidence_angle_raster, | ||
isce3::core::Matrix< float > & | local_incidence_angle_array, | ||
isce3::io::Raster * | projection_angle_raster, | ||
isce3::core::Matrix< float > & | projection_angle_array, | ||
isce3::io::Raster * | simulated_radar_brightness_raster, | ||
isce3::core::Matrix< float > & | simulated_radar_brightness_array, | ||
isce3::core::Vec3 * | terrain_normal_unit_vec_enu = nullptr, | ||
isce3::core::LookSide * | lookside = nullptr ) |
Compute geometry vectors and metadata cube values for a given target point and write to the corresponding arrays.
Each output value is saved onto the first band of its associated raster file.
The target-to-sensor line-of-sight (LOS) and along-track unit vectors are referenced to ENU coordinates computed wrt targets.
The terrain normal unit vector is required for computing the local- incidence angle, projection angle, and simulated radar brightness. The lookside is required for computing the projection angle and simulated radar brightness.
[in] | array_pos_i | Output array line |
[in] | array_pos_j | Output array column |
[in] | native_azimuth_time | Native Doppler azimuth time |
[in] | target_llh | Target lat/lon/height vector |
[in] | orbit | Reference orbit |
[in] | ellipsoid | Reference ellipsoid |
[out] | incidence_angle_raster | Incidence angle cube raster |
[out] | incidence_angle_array | Incidence angle cube array |
[out] | los_unit_vector_x_raster | LOS (target-to-sensor) unit vector X cube raster |
[out] | los_unit_vector_x_array | LOS (target-to-sensor) unit vector X cube array |
[out] | los_unit_vector_y_raster | LOS (target-to-sensor) unit vector Y cube raster |
[out] | los_unit_vector_y_array | LOS (target-to-sensor) unit vector Y cube array |
[out] | along_track_unit_vector_x_raster | Along-track unit vector X raster |
[out] | along_track_unit_vector_x_array | Along-track unit vector X array |
[out] | along_track_unit_vector_y_raster | Along-track unit vector Y raster |
[out] | along_track_unit_vector_y_array | Along-track unit vector Y array |
[out] | elevation_angle_raster | Elevation cube raster |
[out] | elevation_angle_array | Elevation cube array |
[out] | ground_track_velocity_raster | Ground-track velocity raster |
[out] | ground_track_velocity_array | Ground-track velocity array |
[out] | local_incidence_angle_raster | Local-incidence angle raster |
[out] | local_incidence_angle_array | Local-incidence angle array |
[out] | projection_angle_raster | Projection angle raster |
[out] | projection_angle_array | Projection angle array |
[out] | simulated_radar_brightness_raster | Simulated radar brightness raster |
[out] | simulated_radar_brightness_array | Simulated radar brightness array |
[in] | terrain_normal_unit_vec_enu | Terrain normal vector (required to compute local-incidence, projection angles, and simulated radar brightness) |
[in] | lookside | Look side (required to compute the projection angle and simulated radar brightness) |