A collection of functions for EL Null analysis in range direction.
More...
#include <complex>
#include <optional>
#include <vector>
#include <Eigen/Dense>
#include <isce3/core/EMatrix.h>
#include <isce3/core/Linspace.h>
#include "detail/BinarySearchFunc.h"
Go to the source code of this file.
|
namespace | isce3 |
| base interpolator is an abstract base class
|
|
|
using | isce3::antenna::RowMatrixXcf = isce3::core::EMatrix2D<std::complex<float>> |
|
using | isce3::antenna::Linspace_t = isce3::core::Linspace<double> |
|
using | isce3::antenna::tuple_ant = std::tuple<Eigen::ArrayXcd, Eigen::ArrayXcd, Eigen::ArrayXd> |
|
using | isce3::antenna::tuple_echo = std::tuple<Eigen::ArrayXd, Eigen::ArrayXd, detail::ArrayXui> |
|
|
Eigen::ArrayXcd | isce3::antenna::linearInterpComplex1d (const Eigen::Ref< const Eigen::ArrayXd > &x0, const Linspace_t &x, const Eigen::Ref< const Eigen::ArrayXcd > &y) |
| 1-D linear interpolation of a complex array of uniformly-sampled y(x) at new array of points x0.
|
|
tuple_ant | isce3::antenna::genAntennaPairCoefs (const Eigen::Ref< const Eigen::ArrayXcd > &el_cut_left, const Eigen::Ref< const Eigen::ArrayXcd > &el_cut_right, double el_ang_start, double el_ang_step, std::optional< double > el_res_max={}) |
| Generate weighting coefficients used in digital beamforming (DBF) and null formation for a pair of adjacent, partially overlapping beams.
|
|
std::tuple< double, Eigen::Index, double, Eigen::ArrayXd > | isce3::antenna::locateAntennaNull (const Eigen::Ref< const Eigen::ArrayXcd > &coef_left, const Eigen::Ref< const Eigen::ArrayXcd > &coef_right, const Eigen::Ref< const Eigen::ArrayXd > &el_ang_vec) |
| Locate antenna null by forming antenna null and get its min location in EL direction from a pair of EL antenna pattern or its equivalent complex conjugate coeffs.
|
|
tuple_echo | isce3::antenna::formEchoNull (const std::vector< std::complex< float > > &chirp_ref, const Eigen::Ref< const RowMatrixXcf > &echo_left, const Eigen::Ref< const RowMatrixXcf > &echo_right, double sr_start, double sr_spacing, const Eigen::Ref< const Eigen::ArrayXcd > &coef_left, const Eigen::Ref< const Eigen::ArrayXcd > &coef_right, const Eigen::Ref< const Eigen::ArrayXd > &sr_coef) |
| Form x-track echo null power (linear) averaged over range lines and formed from a pair of echoes and a pair of weighting coefs as a function of both EL angles (rad) and slant ranges (m).
|
|
A collection of functions for EL Null analysis in range direction.
◆ formEchoNull()
tuple_echo isce3::antenna::formEchoNull |
( |
const std::vector< std::complex< float > > & | chirp_ref, |
|
|
const Eigen::Ref< const RowMatrixXcf > & | echo_left, |
|
|
const Eigen::Ref< const RowMatrixXcf > & | echo_right, |
|
|
double | sr_start, |
|
|
double | sr_spacing, |
|
|
const Eigen::Ref< const Eigen::ArrayXcd > & | coef_left, |
|
|
const Eigen::Ref< const Eigen::ArrayXcd > & | coef_right, |
|
|
const Eigen::Ref< const Eigen::ArrayXd > & | sr_coef ) |
Form x-track echo null power (linear) averaged over range lines and formed from a pair of echoes and a pair of weighting coefs as a function of both EL angles (rad) and slant ranges (m).
- Parameters
-
[in] | chirp_ref | is complex chirp ref samples used in range compression. |
[in] | echo_left | is complex 2-D array of raw echo samples (pulse by range) for the left RX channel corresponding to the left beam. |
[in] | echo_right | is complex 2-D array of raw echo samples (pulse by range) for the right RX channel corresponding to the right beam. |
[in] | sr_start | is start slant range (m) for both uniformly-sampled echoes in range. |
[in] | sr_spacing | is slant range spacing (m) for both uniformly-sampled echoes in range. |
[in] | coef_left | is complex array of coef for the left beam. |
[in] | coef_right | is complex array of coef for the right beam. |
[in] | sr_coef | array of slant ranges (m) for both left/right coeffs |
- Returns
- array of echo null power pattern in (linear).
-
array of slant range values related to the null power pattern.
-
array of indices used for mapping null slant ranges to antenna EL angles.
- Exceptions
-
◆ genAntennaPairCoefs()
tuple_ant isce3::antenna::genAntennaPairCoefs |
( |
const Eigen::Ref< const Eigen::ArrayXcd > & | el_cut_left, |
|
|
const Eigen::Ref< const Eigen::ArrayXcd > & | el_cut_right, |
|
|
double | el_ang_start, |
|
|
double | el_ang_step, |
|
|
std::optional< double > | el_res_max = {} ) |
Generate weighting coefficients used in digital beamforming (DBF) and null formation for a pair of adjacent, partially overlapping beams.
Computes DBF coefficients by forming the complex conjugate of each input antenna elevation pattern within the inner peak-to-peak region between the two beams. The resulting coefficients may be resampled to finer resolution than the input spacing, if desired. The function returns a vector of coefficients for each beam, as well as a vector of the corresponding elevation angle positions.
- Parameters
-
[in] | el_cut_left | is complex array of uniformly-sampled relative or absolute EL-cut antenna pattern on the left side. |
[in] | el_cut_right | is complex array of uniformly-sampled relative or absolute EL-cut antenna pattern on the right side. It must have the same size as left one! |
[in] | el_ang_start | elevation angle of the first sample in the left/right patterns, in (rad) |
[in] | el_ang_step | elevation angular spacing for left/right patterns in (rad) |
[in] | el_res_max | (optional) max EL angle resolution in (rad). If this value is smaller than the input spacing, el_ang_step , then the outputs will be resampled via linear interpolation to either this finer resolution or a value close to that depending on whether ratio el_ang_step/el_res_max is an integer or not. |
- Returns
- array of complex coeff for the left beam
-
array of complex coeff for the right beam
-
array of uniformly-sampled EL angles in (rad)
- Exceptions
-
InvalidArgument,RuntimeError | |
◆ linearInterpComplex1d()
Eigen::ArrayXcd isce3::antenna::linearInterpComplex1d |
( |
const Eigen::Ref< const Eigen::ArrayXd > & | x0, |
|
|
const Linspace_t & | x, |
|
|
const Eigen::Ref< const Eigen::ArrayXcd > & | y ) |
1-D linear interpolation of a complex array of uniformly-sampled y(x) at new array of points x0.
- Parameters
-
[in] | x0 | is array of values to be interpolated |
[in] | x | is isce3 Linspace object for uniformly sampled input "x" |
[in] | y | is array of function value "y(x)". |
- Returns
- array of interpolated complex values at samples "x0".
- Exceptions
-
◆ locateAntennaNull()
std::tuple< double, Eigen::Index, double, Eigen::ArrayXd > isce3::antenna::locateAntennaNull |
( |
const Eigen::Ref< const Eigen::ArrayXcd > & | coef_left, |
|
|
const Eigen::Ref< const Eigen::ArrayXcd > & | coef_right, |
|
|
const Eigen::Ref< const Eigen::ArrayXd > & | el_ang_vec ) |
Locate antenna null by forming antenna null and get its min location in EL direction from a pair of EL antenna pattern or its equivalent complex conjugate coeffs.
- Parameters
-
[in] | coef_left | is complex array of coef for the left beam |
[in] | coef_right | is complex array of coef for the right beam |
[in] | el_ang_vec | is array of EL angles in (rad) |
- Returns
- el angle of the null location in (rad)
-
index of null location
-
peak-normalized null magnitude in (linear)
-
1-D peak-normalized antenna null power pattern (linear) with the same size as
el_ang_vec
.
- Exceptions
-
1.13.2.