12 #include "Constants.h"
21 using Map =
typename Eigen::Map<const EArray2D<U>>;
30 virtual U
interp_impl(
double x,
double y,
const Map& map)
const = 0;
50 const Map z {&z_data[0],
51 static_cast<Eigen::Index
>(z_data.size() / width),
52 static_cast<Eigen::Index>(width)};
60 const Map z {&z_data[0],
61 static_cast<Eigen::Index
>(z_data.size() / width),
62 static_cast<Eigen::Index>(width)};
79 using super_t = Interpolator<U>;
80 using typename super_t::Map;
83 U interp_impl(
double x,
double y,
const Map& z)
const override;
97 using super_t = Interpolator<U>;
98 using typename super_t::Map;
101 U interp_impl(
double x,
double y,
const Map& z)
const override;
116 using super_t = Interpolator<U>;
117 using typename super_t::Map;
120 U interp_impl(
double x,
double y,
const Map& z)
const override;
134 using super_t = Interpolator<U>;
135 using typename super_t::Map;
142 U
interp_impl(
double x,
double y,
const Map& z)
const override;
153 void _initSpline(
const std::valarray<U>&,
int, std::valarray<U>&,
154 std::valarray<U>&)
const;
156 U _spline(
double,
const std::valarray<U>&,
int,
157 const std::valarray<U>&)
const;
164 using super_t = Interpolator<U>;
165 using typename super_t::Map;
168 U interp_impl(
double x,
double y,
const Map& z)
const override;
179 void _sinc_coef(
double beta,
double relfiltlen,
int decfactor,
180 double pedestal,
int weight,
181 std::valarray<double>& filter)
const;
184 U _sinc_eval_2d(
const Map& z,
int intpx,
int intpy,
double frpx,
188 Matrix<double> _kernel;
189 int _kernelLength, _kernelWidth, _sincHalf;
193 namespace isce3 {
namespace core {
198 inline Interpolator<U>*
200 int sincLen = SINC_LEN,
int sincSub = SINC_SUB)
202 if (method == BILINEAR_METHOD) {
204 }
else if (method == BICUBIC_METHOD) {
206 }
else if (method == BIQUINTIC_METHOD) {
208 }
else if (method == NEAREST_METHOD) {
210 }
else if (method == SINC_METHOD) {
218 template<
class U,
class V>
219 U
sinc_eval(
const Matrix<U>&,
const Matrix<V>&,
int,
int,
double,
int);
dataInterpMethod
Enumeration type to indicate interpolation method.
Definition: Constants.h:23
virtual U interp_impl(double x, double y, const Map &map) const =0
Base implementation for all types.
Definition of BicubicInterpolator.
Definition: forward.h:35
U sinc_eval(const Matrix< U > &, const Matrix< V > &, int, int, double, int)
1-D sinc evaluation
Definition: Interpolator.cpp:11
virtual ~Interpolator()
Virtual destructor (allow destruction of base Interpolator pointer)
Definition: Interpolator.h:26
U interpolate(double x, double y, const Matrix< U > &z) const
Interpolate at a given coordinate for an input isce3::core::Matrix.
Definition: Interpolator.h:41
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.
Definition: Interpolator.h:199
NearestNeighborInterpolator()
Default constructor.
Definition: Interpolator.h:124
U interpolate(double x, double y, std::valarray< U > &z_data, size_t width) const
Interpolate at a given coordinate for data passed as a valarray.
Definition: Interpolator.h:47
Sinc2dInterpolator(int sincLen, int sincSub)
Default constructor.
Definition: Sinc2dInterpolator.cpp:14
U interp_impl(double x, double y, const Map &z) const override
Interpolate at a given coordinate.
Definition: Spline2dInterpolator.cpp:33
Data structure for a 2D row-major matrix.
Definition: forward.h:31
BicubicInterpolator()
Default constructor.
Definition: Interpolator.h:105
Spline2dInterpolator(size_t order)
Default constructor.
Definition: Spline2dInterpolator.cpp:12
U interpolate(double x, double y, const Map &map) const
Interpolate at a given coordinate for an input Eigen::Map.
Definition: Interpolator.h:35
BilinearInterpolator()
Default constructor.
Definition: Interpolator.h:87
Definition of NearestNeighborInterpolator.
Definition: forward.h:36
dataInterpMethod method() const
Return interpolation method.
Definition: Interpolator.h:67
Definition of Sinc2dInterpolator.
Definition: forward.h:38
Definition of BilinearInterpolator.
Definition: forward.h:34
Definition of Spline2dInterpolator.
Definition: forward.h:37
U interpolate(double x, double y, std::vector< U > &z_data, size_t width) const
Interpolate at a given coordinate for data passed as a vector.
Definition: Interpolator.h:57
Definition of parent Interpolator.
Definition: forward.h:33