isce3  0.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
Classes | Typedefs | Enumerations | Functions | Variables
isce3::core Namespace Reference

The isce3::core namespace. More...

Classes

class  Attitude
 Base class for attitude data representation. More...
 
class  Baseline
 Data structure for computing interferometric baselines. More...
 
class  Basis
 Simple class to store three-dimensional basis vectors. More...
 
struct  FixedString
 Struct with fixed-length string for serialization. More...
 
struct  double_promote
 Precision-promotion to double/complex<double> More...
 
struct  double_promote< float >
 Template specialization for float. More...
 
struct  double_promote< double >
 Template specialization for double. More...
 
struct  double_promote< std::complex< float > >
 Template specialization for complex<float> More...
 
struct  double_promote< std::complex< double > >
 Template specialization for complex<double> More...
 
class  Cube
 
class  DateTime
 Data structure to store date time to nano-sec precision. More...
 
class  DenseMatrix
 
class  Ellipsoid
 Data structure to store Ellipsoid information. More...
 
class  EulerAngles
 Data structure for Euler Angle representation of attitude information. More...
 
class  Vector
 
class  Linspace
 A uniformly-spaced sequence of values over some interval. More...
 
class  LUT1d
 Data structure to hold a 1D Lookup table. More...
 
class  LUT2d
 Data structure to store 2D Lookup table. More...
 
class  Matrix
 Data structure for a 2D row-major matrix. More...
 
class  Interpolator
 Definition of parent Interpolator. More...
 
class  BilinearInterpolator
 Definition of BilinearInterpolator. More...
 
class  BicubicInterpolator
 Definition of BicubicInterpolator. More...
 
class  NearestNeighborInterpolator
 Definition of NearestNeighborInterpolator. More...
 
class  Spline2dInterpolator
 Definition of Spline2dInterpolator. More...
 
class  Sinc2dInterpolator
 Definition of Sinc2dInterpolator. More...
 
class  Kernel
 Abstract base class for all kernels. More...
 
class  BartlettKernel
 Bartlett kernel (triangle function). More...
 
class  KnabKernel
 Kernel based on the paper by Knab for interpolating band-limited signals. More...
 
class  LinearKernel
 Linear kernel, which is just a special case of Bartlett. More...
 
class  NFFTKernel
 NFFT time-domain kernel. More...
 
class  TabulatedKernel
 Tabulated Kernel. More...
 
class  ChebyKernel
 Polynomial Kernel. More...
 
class  Metadata
 Data structure for storing basic radar geometry image metadata. More...
 
class  Orbit
 Sequence of platform ephemeris samples (state vectors) with uniform temporal spacing, supporting efficient lookup and interpolation. More...
 
class  Peg
 Data structure to store a peg point. More...
 
struct  Pegtrans
 Data structure to assist with Peg point transformations. More...
 
class  Pixel
 Helper datastructure to handle slant range information for a pixel. More...
 
class  Poly1d
 Data structure for representing 1D polynomials. More...
 
class  Poly2d
 Data structure for representing 1D polynomials. More...
 
class  ProjectionBase
 Abstract base class for individual projections. More...
 
class  LonLat
 Standard WGS84 Lon/Lat Projection extension of ProjBase - EPSG:4326. More...
 
class  Geocent
 Standard WGS84 ECEF coordinates extension of ProjBase - EPSG:4978. More...
 
class  UTM
 UTM coordinate extension of ProjBase. More...
 
class  PolarStereo
 Polar stereographic extension of ProjBase. More...
 
class  CEA
 Equal Area Projection extension of ProjBase. More...
 
class  Quaternion
 Quaternion representation of attitude information. More...
 
struct  StateVector
 
class  TimeDelta
 Data structure to store TimeDelta to double precision seconds. More...
 

Typedefs

template<typename T , int rows = Eigen::Dynamic, int cols = Eigen::Dynamic>
using EMatrix2D = typename Eigen::Matrix< T, rows, cols, Eigen::RowMajor >
 
template<typename T , int rows = Eigen::Dynamic, int cols = Eigen::Dynamic>
using EArray2D = typename Eigen::Array< T, rows, cols, Eigen::RowMajor >
 
using Mat3 = DenseMatrix< 3 >
 
using Vec3 = Vector< 3 >
 
using cartmat_t = Mat3
 
using cartesian_t = Vec3
 

Enumerations

enum  dataInterpMethod {
  SINC_METHOD = 0, BILINEAR_METHOD = 1, BICUBIC_METHOD = 2, NEAREST_METHOD = 3,
  BIQUINTIC_METHOD = 4
}
 Enumeration type to indicate interpolation method.
 
enum  LookSide { LookSide::Left = 1, LookSide::Right = -1 }
 Side that radar looks at, Left or Right. More...
 
enum  OrbitInterpMethod { OrbitInterpMethod::Hermite, OrbitInterpMethod::Legendre }
 Orbit interpolation method. More...
 
enum  OrbitInterpBorderMode { OrbitInterpBorderMode::Error, OrbitInterpBorderMode::Extrapolate, OrbitInterpBorderMode::FillNaN }
 Mode determining how interpolation outside of orbit domain is handled. More...
 

Functions

dataInterpMethod parseDataInterpMethod (const std::string &method)
 Convert string to dataInterpMethod.
 
double decimaldeg2meters (double deg)
 Convert decimal degrees to meters approximately.
 
std::ostream & operator<< (std::ostream &os, const DateTime &datetime)
 
bool _is_leap (int)
 
int _days_in_month (int, int)
 
int _days_before_year (int)
 
int _days_before_month (int, int)
 
int _ymd_to_ord (int, int, int)
 
void _ord_to_ymd (int, int &, int &, int &)
 
template<typename TK , typename TD >
TD interp1d (const Kernel< TK > &kernel, const TD *x, size_t length, size_t stride, double t, bool periodic=false)
 Interpolate sequence x at point t. More...
 
template<typename TK , typename TD >
TD interp1d (const Kernel< TK > &kernel, const std::valarray< TD > &x, double t, bool periodic=false)
 Interpolate sequence x at point t. More...
 
template<typename U >
Interpolator< U > * createInterpolator (dataInterpMethod method, size_t order=6, int sincLen=SINC_LEN, int sincSub=SINC_SUB)
 Utility function to create interpolator pointer given an interpolator enum type.
 
template<class U , class V >
sinc_eval (const Matrix< U > &, const Matrix< V > &, int, int, double, int)
 1-D sinc evaluation
 
template<typename T >
_sampling_window (T t, T halfwidth, T bandwidth)
 
template<typename T , typename U >
CUDA_HOSTDEV constexpr bool operator== (const Linspace< T > &, const Linspace< U > &)
 
template<typename T , typename U >
CUDA_HOSTDEV constexpr bool operator!= (const Linspace< T > &, const Linspace< U > &)
 
LookSide parseLookSide (const std::string &str)
 Parse string (e.g., "left" or "right") to enum LookSide. More...
 
std::string to_string (LookSide d)
 Convert enum LookSide to string ("left" or "right"). More...
 
std::ostream & operator<< (std::ostream &out, const LookSide d)
 
bool operator== (const Orbit &lhs, const Orbit &rhs)
 
bool operator!= (const Orbit &lhs, const Orbit &rhs)
 
constexpr int minStateVecs (OrbitInterpMethod method)
 Get minimum number of orbit state vectors required for interpolation with specified method.
 
ProjectionBasecreateProj (int epsg)
 
std::unique_ptr< ProjectionBasemakeProjection (int epsg)
 
int projTransform (ProjectionBase *in, ProjectionBase *out, const Vec3 &inpts, Vec3 &outpts)
 
template<typename T >
void load_archive (std::string metadata, char *objectTag, T *object)
 
template<typename T >
void load_archive_reference (std::string metadata, char *objectTag, T &object)
 
template<class Archive >
void save (Archive &archive, const Ellipsoid &ellps)
 
template<class Archive >
void load (Archive &archive, Ellipsoid &ellps)
 
void loadFromH5 (isce3::io::IGroup &group, Ellipsoid &ellps)
 Load Ellipsoid parameters from HDF5. More...
 
template<class Archive >
void save (Archive &archive, const Orbit &orbit)
 
template<class Archive >
void load (Archive &archive, Orbit &orbit)
 
void loadFromH5 (isce3::io::IGroup &group, Orbit &orbit)
 Load orbit data from HDF5 product. More...
 
void saveToH5 (isce3::io::IGroup &group, const Orbit &orbit)
 Save orbit data to HDF5 product. More...
 
void loadFromH5 (isce3::io::IGroup &group, EulerAngles &euler)
 Load Euler angle data from HDF5 product. More...
 
void saveToH5 (isce3::io::IGroup &group, const EulerAngles &euler)
 Save Euler angle data to HDF5 product. More...
 
void loadFromH5 (isce3::io::IGroup &group, Quaternion &qs)
 
void saveToH5 (isce3::io::IGroup &group, const Quaternion &qs)
 
template<class Archive >
void save (Archive &archive, const Metadata &meta)
 
template<class Archive >
void load (Archive &archive, Metadata &meta)
 
template<class Archive >
void serialize (Archive &archive, Poly2d &poly)
 
void loadFromH5 (isce3::io::IGroup &group, Poly2d &poly, std::string name)
 Load polynomial coefficients from HDF5 product. More...
 
template<typename T >
void loadCalGrid (isce3::io::IGroup &group, const std::string &dsetName, isce3::core::LUT2d< T > &lut)
 Load LUT2d data from HDF5 product. More...
 
template<typename T >
void saveCalGrid (isce3::io::IGroup &group, const std::string &dsetName, const isce3::core::LUT2d< T > &lut, const isce3::core::DateTime &refEpoch, const std::string &units="")
 Save LUT2d data to HDF5 product. More...
 
template<class Archive , typename T >
void save (Archive &archive, LUT1d< T > const &lut)
 
template<class Archive , typename T >
void load (Archive &archive, LUT1d< T > &lut)
 
template<typename T >
void loadFromH5 (isce3::io::IGroup &group, LUT1d< T > &lut, std::string name_coords, std::string name_values)
 Load polynomial coefficients from HDF5 product. More...
 
template<class Archive >
void save (Archive &archive, const StateVector &sv)
 
template<class Archive >
void load (Archive &archive, StateVector &sv)
 
bool operator== (const StateVector &lhs, const StateVector &rhs)
 
bool operator!= (const StateVector &lhs, const StateVector &rhs)
 
TimeDelta operator* (double lhs, const TimeDelta &rhs)
 
template<typename T >
void checkVecLenDebug (const std::vector< T > &vec, size_t len, const char *vec_name, const char *parent_func)
 Inline function for input checking on vector lengths (primarily to check to see if 3D vector has the correct number of inputs, but is generalized to any length). More...
 
template<typename T >
void check2dVecLenDebug (const std::vector< std::vector< T >> &vec, size_t len, size_t width, const char *vec_name, const char *parent_func)
 
template<typename T >
std::vector< T > arange (T low, T high, T increment)
 Function to return a Python style arange vector. More...
 
template<typename T >
std::vector< T > linspace (T low, T high, std::size_t number)
 Function to return a Matlab/Python style linspace vector. More...
 
template<typename T >
void linspace (T low, T high, std::vector< T > &data)
 Function to fill a Matlab/Python style linspace vector. More...
 
template<typename T >
void linspace (T low, T high, std::valarray< T > &data)
 Function to fill a Matlab/Python style linspace valarray. More...
 
template<typename T >
std::vector< T > logspace (T first, T last, std::size_t number, T base=10)
 Function to return a Matlab/Python style logspace vector. More...
 
template<typename T >
std::vector< T > stringToVector (const std::string &str)
 Parse a space-separated string into a vector. More...
 
template<typename T >
std::string vectorToString (const std::vector< T > &vec)
 Generate a string of space-separated values from a vector. More...
 
std::string stringFromVRT (const char *filename, int bandNum=1)
 Read metadata from a VRT file and return an std::string. More...
 
template<typename T >
bool compareFloatingPoint (T first, T second)
 Combined absolute and relative tolerance test. More...
 
template<typename T >
bool compareComplex (std::complex< T > first, std::complex< T > second)
 Tolerance test for real and imaginary components for scalar value. More...
 
template<typename T >
std::valarray< bool > makeMask (const std::valarray< std::complex< T >> &v, std::complex< T > noDataValue)
 making a boolean mask for an array based on noDataValue More...
 
template<typename T >
std::valarray< bool > makeMask (const std::valarray< T > &v, T noDataValue)
 making a boolean mask for an array based on noDataValue More...
 
void insertionSort (std::valarray< double > &a, std::valarray< double > &b, std::valarray< double > &c)
 Sort arrays a, b, c by the values in array a.
 
void insertionSort (std::valarray< double > &a, std::valarray< double > &b)
 Sort arrays a and b by the values in array a.
 
int binarySearch (const std::valarray< double > &array, double value)
 Searches array for index closest to provided value.
 
template<class T >
const T & clamp (const T &x, const T &lower, const T &upper)
 Clip a number between an upper and lower range (implements std::clamp for older GCC)
 
CUDA_HOSTDEV Vec3 normalPlane (const Vec3 &p1, const Vec3 &p2, const Vec3 &p3)
 

Variables

const int SINC_HALF = 4
 Default sinc parameters.
 
const int SINC_LEN = 8
 
const int SINC_ONE = 9
 
const int SINC_SUB = 8192
 
const double EarthSemiMajorAxis = 6378137.0
 Semi-major axis for WGS84.
 
const double EarthEccentricitySquared = 0.0066943799901
 Eccentricity^2 for WGS84.
 
const double GLOBAL_MIN_HEIGHT = -500.0
 Global minimum height.
 
const double GLOBAL_MAX_HEIGHT = 9000.0
 Global maximum height.
 
const short SHADOW_VALUE = 1
 Layover and shadow values.
 
const short LAYOVER_VALUE = 2
 
const int AREA_PROJECTION_RADAR_GRID_MARGIN = 100
 Area projection algorithm radar grid margin.
 
const float AREA_PROJECTION_MIN_VALID_SAMPLES_RATIO = 0.75
 
const DateTime MIN_DATE_TIME = DateTime(1970, 1, 1)
 
const std::string UNINITIALIZED_STRING = "uninitialized"
 
const int MIN_TO_SEC =60
 
const int HOUR_TO_SEC =3600
 
const int HOUR_TO_MIN =60
 
const int DAY_TO_SEC =86400
 
const int DAY_TO_MIN =1440
 
const int DAY_TO_HOUR =24
 

Detailed Description

The isce3::core namespace.

Enumeration Type Documentation

enum isce3::core::LookSide
strong

Side that radar looks at, Left or Right.

Enumerator
Left 

Radar points to left/port side of vehicle.

Right 

Radar points to right/starboard side of vehicle.

Mode determining how interpolation outside of orbit domain is handled.

Enumerator
Error 

Raise error if interpolation attempted outside orbit domain.

Extrapolate 

Allow extrapolation to points outside orbit domain.

FillNaN 

Output NaN for interpolation points outside orbit domain.

Orbit interpolation method.

See the documentation for a description of each method.

Enumerator
Hermite 

Third-order Hermite polynomial interpolation.

Legendre 

Eighth-order Legendre polynomial interpolation.

Function Documentation

template<typename T >
std::vector<T> isce3::core::arange ( low,
high,
increment 
)
inline

Function to return a Python style arange vector.

Parameters
[in]lowStarting value of vector
[in]highEnding value of vector (non-inclusive)
[in]incrementSpacing between vector elements
template<typename T >
void isce3::core::checkVecLenDebug ( const std::vector< T > &  vec,
size_t  len,
const char *  vec_name,
const char *  parent_func 
)
inline

Inline function for input checking on vector lengths (primarily to check to see if 3D vector has the correct number of inputs, but is generalized to any length).

'vec_name' is passed in by the wrapper macro (stringified vector name), and 'parent_func' is passed in the same way through PRETTY_FUNCTION

template<typename T >
bool isce3::core::compareComplex ( std::complex< T >  first,
std::complex< T >  second 
)
inline

Tolerance test for real and imaginary components for scalar value.

Parameters
[in]firstfirst complex value
[in]secondsecond complex value
[out]returnscomparison result (true or false)
template<typename T >
bool isce3::core::compareFloatingPoint ( first,
second 
)
inline

Combined absolute and relative tolerance test.

Parameters
[in]firstfirst value
[in]secondsecond value
[out]returnscomparison result (true or false)
template<typename TK , typename TD >
TD isce3::core::interp1d ( const Kernel< TK > &  kernel,
const TD *  x,
size_t  length,
size_t  stride,
double  t,
bool  periodic = false 
)

Interpolate sequence x at point t.

Parameters
[in]kernelKernel function to use for interpolation.
[in]xSequence to interpolate.
[in]lengthLength of sequence.
[in]strideStride between elements of sequence.
[in]tDesired time sample (0 <= t <= x.size()-1).
[in]periodicUse periodic boundary condition. Default = false.
Returns
Interpolated value or 0 if kernel would run off array.

Sequence x will be addressed as x[i*stride] for 0 <= i < length. Template parameter TK is kernel element type, TD is data element type.

template<typename TK , typename TD >
TD isce3::core::interp1d ( const Kernel< TK > &  kernel,
const std::valarray< TD > &  x,
double  t,
bool  periodic = false 
)

Interpolate sequence x at point t.

Parameters
[in]kernelKernel function to use for interpolation.
[in]xSequence to interpolate.
[in]tDesired time sample (0 <= t <= x.size()-1).
[in]periodicUse periodic boundary condition. Default = false.
Returns
Interpolated value or 0 if kernel would run off array.

Template parameter TK is kernel element type, TD is data element type.

template<typename T >
std::vector<T> isce3::core::linspace ( low,
high,
std::size_t  number 
)
inline

Function to return a Matlab/Python style linspace vector.

Parameters
[in]lowStarting value of vector
[in]highEnding value of vector (inclusive)
[in]numberNumber of vector elements
template<typename T >
void isce3::core::linspace ( low,
high,
std::vector< T > &  data 
)
inline

Function to fill a Matlab/Python style linspace vector.

Parameters
[in]lowStarting value of vector
[in]highEnding value of vector (inclusive)
[out]dataVector to fill in
template<typename T >
void isce3::core::linspace ( low,
high,
std::valarray< T > &  data 
)
inline

Function to fill a Matlab/Python style linspace valarray.

Parameters
[in]lowStarting value of vector
[in]highEnding value of vector (inclusive)
[in]dataVector to fill in
template<typename T >
void isce3::core::loadCalGrid ( isce3::io::IGroup group,
const std::string &  dsetName,
isce3::core::LUT2d< T > &  lut 
)
inline

Load LUT2d data from HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]dsetNameDataset name within group
[in]lutLUT2d to be configured.
void isce3::core::loadFromH5 ( isce3::io::IGroup group,
Ellipsoid ellps 
)
inline

Load Ellipsoid parameters from HDF5.

Parameters
[in]groupHDF5 group object.
[in]ellpsEllipsoid object to be configured.
void isce3::core::loadFromH5 ( isce3::io::IGroup group,
Orbit orbit 
)
inline

Load orbit data from HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]orbitOrbit object to be configured.
void isce3::core::loadFromH5 ( isce3::io::IGroup group,
EulerAngles euler 
)
inline

Load Euler angle data from HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]eulerEulerAngles object to be configured.
void isce3::core::loadFromH5 ( isce3::io::IGroup group,
Poly2d poly,
std::string  name 
)
inline

Load polynomial coefficients from HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]polyPoly2d to be configured.
[in]nameDataset name within group.
template<typename T >
void isce3::core::loadFromH5 ( isce3::io::IGroup group,
LUT1d< T > &  lut,
std::string  name_coords,
std::string  name_values 
)
inline

Load polynomial coefficients from HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]polyPoly2d to be configured.
[in]nameDataset name within group.
template<typename T >
std::vector<T> isce3::core::logspace ( first,
last,
std::size_t  number,
base = 10 
)
inline

Function to return a Matlab/Python style logspace vector.

Parameters
[in]firstbase^first is the starting value of the vector
[in]lastbase^last is the ending value of the vector
[in]Numberof vector elements
[in]baseBase of the log space
template<typename T >
std::valarray<bool> isce3::core::makeMask ( const std::valarray< std::complex< T >> &  v,
std::complex< T >  noDataValue 
)
inline

making a boolean mask for an array based on noDataValue

Parameters
[in]vinout array
[in]noDataValuethe value used to create the mask
[out]returnsa boolean mask
template<typename T >
std::valarray<bool> isce3::core::makeMask ( const std::valarray< T > &  v,
noDataValue 
)
inline

making a boolean mask for an array based on noDataValue

Parameters
[in]vinout array
[in]noDataValuethe value used to create the mask
[out]returnsa boolean mask
LookSide isce3::core::parseLookSide ( const std::string &  str)

Parse string (e.g., "left" or "right") to enum LookSide.

template<typename T >
void isce3::core::saveCalGrid ( isce3::io::IGroup group,
const std::string &  dsetName,
const isce3::core::LUT2d< T > &  lut,
const isce3::core::DateTime refEpoch,
const std::string &  units = "" 
)
inline

Save LUT2d data to HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]dsetNameDataset name within group
[in]lutLUT2d to be saved.
[in]unitsUnits of LUT2d data.
void isce3::core::saveToH5 ( isce3::io::IGroup group,
const Orbit orbit 
)
inline

Save orbit data to HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]orbitOrbit object to be saved.
void isce3::core::saveToH5 ( isce3::io::IGroup group,
const EulerAngles euler 
)
inline

Save Euler angle data to HDF5 product.

Parameters
[in]groupHDF5 group object.
[in]eulerEulerAngles object to be save.
std::string isce3::core::stringFromVRT ( const char *  filename,
int  bandNum = 1 
)

Read metadata from a VRT file and return an std::string.

Parameters
[in]filenameVRT product filename
[in]bandNumBand number to retrieve metadata from
[out]metastd::string containing metadata
template<typename T >
std::vector<T> isce3::core::stringToVector ( const std::string &  str)
inline

Parse a space-separated string into a vector.

Parameters
[in]strString of values
[out]vecVector of values
std::string isce3::core::to_string ( LookSide  d)

Convert enum LookSide to string ("left" or "right").

template<typename T >
std::string isce3::core::vectorToString ( const std::vector< T > &  vec)
inline

Generate a string of space-separated values from a vector.

Parameters
[in]vecVector of values
[out]strString of values

Generated for ISCE3.0 by doxygen 1.8.5.