A class for estimating one-way or two-way elevation (EL) power pattern from 2-D raw echo data over quasi-homogenous scene and provide approximate look (off-nadir or EL) angle and ellipsoidal incidence angle (no local slope).
More...
#include <ElPatternEst.h>
|
using | RowMatrixXcf = isce3::core::EMatrix2D<std::complex<float>> |
|
|
| ElPatternEst (double sr_start, const isce3::core::Orbit &orbit, int polyfit_deg=6, const isce3::geometry::DEMInterpolator &dem_interp={}, double win_ped=0.0, const isce3::core::Ellipsoid &ellips={}, bool center_scale_pf=false) |
| A constructor with full input arguments.
|
|
| ElPatternEst (double sr_start, const isce3::core::Orbit &orbit, const isce3::geometry::DEMInterpolator &dem_interp) |
| A more concise constructor of key inputs.
|
|
tuple5_t | powerPattern2way (const Eigen::Ref< const RowMatrixXcf > &echo_mat, double sr_spacing, double chirp_rate, double chirp_dur, std::optional< double > az_time={}, int size_avg=8, bool inc_corr=true) const |
| Estimated averaged two-way time-series power pattern in Elevation from 2-D raw echo data uniformly sampled in slant range direction.
|
|
tuple5_t | powerPattern1way (const Eigen::Ref< const RowMatrixXcf > &echo_mat, double sr_spacing, double chirp_rate, double chirp_dur, std::optional< double > az_time={}, int size_avg=8, bool inc_corr=true) const |
| Estimated averaged one-way time-series power pattern in Elevation from 2-D raw echo data uniformly sampled in slant range direction.
|
|
double | winPed () const |
| Get raised-cosine window pedestal set at the constructor.
|
|
int | polyfitDeg () const |
| Get degree of polyfit set at the constructor.
|
|
bool | isCenterScalePolyfit () const |
| Check whether the polyfit process is centeralized and scaled.
|
|
A class for estimating one-way or two-way elevation (EL) power pattern from 2-D raw echo data over quasi-homogenous scene and provide approximate look (off-nadir or EL) angle and ellipsoidal incidence angle (no local slope).
The final power in dB scale is fit into N-order polynomials as a function of look angle in radians. Required relative radiometric calibration from \(\beta^{0}\) to \(\gamma^{0}\) is done approximately by using slant ranges and ellipsoidal incidence angle.
◆ tuple4_t
using isce3::antenna::ElPatternEst::tuple4_t |
|
protected |
Initial value: std::tuple<Eigen::ArrayXd, isce3::core::Linspace<double>,
Eigen::ArrayXd, Eigen::ArrayXd>
◆ tuple5_t
using isce3::antenna::ElPatternEst::tuple5_t |
|
protected |
Initial value: std::tuple<Eigen::ArrayXd, isce3::core::Linspace<double>,
Data structure for representing 1D polynomials.
Definition Poly1d.h:23
◆ ElPatternEst() [1/2]
A constructor with full input arguments.
- Parameters
-
[in] | sr_start | start slant range in (m) |
[in] | orbit | isce3 Orbit object |
[in] | polyfit_deg | (optional) polyfit degree of polynomial for polyfitting estimated power patterns. The value must be > 1. Default is 6. @pram[in] dem_interp (optional) isce3 DEMInterpolator object. Default is global 0.0 (m) height (reference ellipsoid). Note that simply averaged DEM will be used. No local slope is taken into account! |
[in] | win_ped | (optional) Raised-cosine window pedestal. A value within [0, 1]. Default is Hann window.https://en.wikipedia.org/wiki/Hann_function See raised cosine window. |
[in] | ellips | (optional) isce3 Ellipsoid object. Default is WGS84. |
[in] | center_scale_pf | (optional) whether or not use center and scaling in polyfit. Default is false. |
- Exceptions
-
◆ ElPatternEst() [2/2]
A more concise constructor of key inputs.
- Parameters
-
- Exceptions
-
◆ isCenterScalePolyfit()
bool isce3::antenna::ElPatternEst::isCenterScalePolyfit |
( |
| ) |
const |
|
inline |
Check whether the polyfit process is centeralized and scaled.
- Returns
- bool
◆ polyfitDeg()
int isce3::antenna::ElPatternEst::polyfitDeg |
( |
| ) |
const |
|
inline |
Get degree of polyfit set at the constructor.
- Returns
- degree of polyfit used for 1-way or 2-way power pattern (dB) as a function of look angles (rad).
◆ powerPattern1way()
ElPatternEst::tuple5_t isce3::antenna::ElPatternEst::powerPattern1way |
( |
const Eigen::Ref< const RowMatrixXcf > & | echo_mat, |
|
|
double | sr_spacing, |
|
|
double | chirp_rate, |
|
|
double | chirp_dur, |
|
|
std::optional< double > | az_time = {}, |
|
|
int | size_avg = 8, |
|
|
bool | inc_corr = true ) const |
Estimated averaged one-way time-series power pattern in Elevation from 2-D raw echo data uniformly sampled in slant range direction.
Uniform sampling in azimuth direction is not required! Note that it is assumed the data either has no TX gap (bad values) or its TX gap (bad values) have been already replaced by meaningful data.
- Parameters
-
[in] | echo_mat | raw echo matrix, a row-major Eigen matrix of type complex float. The rows represent range lines. The matrix shape is pulses (azimuth bins) by range bins. |
[in] | sr_spacing | slant range spacing in (m). |
[in] | chirp_rate | transmit chirp rate in (Hz/sec). |
[in] | chirp_dur | transmit chirp duration in (sec). |
[in] | az_time | (optional) relative azimuth time in seconds w.r.t reference epoch time of orbit object. Default is the mid orbit time if not specified or if set to {} or std::nullopt. |
[in] | size_avg | (optional) the block size for averaging in slant range direction. Default is 8. |
[in] | inc_corr | (optional) whether or not apply correction for incidence angles calculated at a certain mean DEM height wrt reference ellipsoid w/o taking into account any local slope! Default is true. |
- Returns
- eigen vector of times-series range-averaged peak-normalized 1-way power pattern in (dB) which is uniformly sampled in slant range with new range spacing defined by "sr_spacing * size_avg".
-
slant ranges in meters in the form of isce3 Linspace object
-
look angles vector in radians.
-
ellipsoidal incidence angles vector in radians.
-
isce3 Poly1d object mapping power in dB as a function look angle in radians
- Exceptions
-
InvalidArgument,RuntimeError | |
- See also
- powerPattern2way()
◆ powerPattern2way()
ElPatternEst::tuple5_t isce3::antenna::ElPatternEst::powerPattern2way |
( |
const Eigen::Ref< const RowMatrixXcf > & | echo_mat, |
|
|
double | sr_spacing, |
|
|
double | chirp_rate, |
|
|
double | chirp_dur, |
|
|
std::optional< double > | az_time = {}, |
|
|
int | size_avg = 8, |
|
|
bool | inc_corr = true ) const |
Estimated averaged two-way time-series power pattern in Elevation from 2-D raw echo data uniformly sampled in slant range direction.
Uniform sampling in azimuth direction is not required! Note that it is assumed the data either has no TX gap (bad values) or its TX gap (bad values) have been already replaced by meaningful data.
- Parameters
-
[in] | echo_mat | raw echo matrix, a row-major Eigen matrix of type complex float. The rows represent range lines. The matrix shape is pulses (azimuth bins) by range bins. |
[in] | sr_spacing | slant range spacing in (m). |
[in] | chirp_rate | transmit chirp rate in (Hz/sec). |
[in] | chirp_dur | transmit chirp duration in (sec). |
[in] | az_time | (optional) relative azimuth time in seconds w.r.t reference epoch time of orbit object. Default is the mid orbit time if not specified or if set to {} or std::nullopt. |
[in] | size_avg | (optional) the block size for averaging in slant range direction. Default is 8. |
[in] | inc_corr | (optional) whether or not apply correction for incidence angles calculated at a certain mean DEM height wrt reference ellipsoid w/o taking into account any local slope! Default is true. |
- Returns
- eigen vector of times-series range-averaged peak-normalized 2-way power pattern in (dB) which is uniformly sampled in slant range with new range spacing defined by "sr_spacing * size_avg".
-
slant ranges in meters in the form of isce3 Linspace object
-
look angles vector in radians.
-
ellipsoidal incidence angles vector in radians.
-
isce3 Poly1d object mapping power in dB as a function look angle in radians
- Exceptions
-
InvalidArgument,RuntimeError | |
- See also
- powerPattern1way()
◆ winPed()
double isce3::antenna::ElPatternEst::winPed |
( |
| ) |
const |
|
inline |
Get raised-cosine window pedestal set at the constructor.
- Returns
- window pedestal used for weighting ref chirp in range compression.
The documentation for this class was generated from the following files:
- /home/runner/work/isce3/isce3/cxx/isce3/antenna/ElPatternEst.h
- /home/runner/work/isce3/isce3/cxx/isce3/antenna/ElPatternEst.cpp
1.13.2.