11 #include "Constants.h"
13 #include "Utilities.h"
24 LUT2d(
double xstart,
double ystart,
double dx,
double dy,
27 bool boundsError =
true);
28 LUT2d(
const std::valarray<double> & xcoord,
29 const std::valarray<double> & ycoord,
32 bool boundsError =
true);
35 inline LUT2d(
const LUT2d<T> & lut);
38 inline LUT2d & operator=(
const LUT2d<T> & lut);
41 void setFromData(
const std::valarray<double> & xcoord,
42 const std::valarray<double> & ycoord,
52 inline double xStart()
const {
return _xstart; }
54 inline double yStart()
const {
return _ystart; }
56 inline double xSpacing()
const {
return _dx; }
58 inline double ySpacing()
const {
return _dy; }
60 inline size_t length()
const {
return _data.
length(); }
62 inline size_t width()
const {
return _data.
width(); }
64 inline T refValue()
const {
return _refValue; }
66 inline bool haveData()
const {
return _haveData; }
68 inline bool boundsError()
const {
return _boundsError; }
73 inline void boundsError(
bool flag) { _boundsError = flag; }
76 T
eval(
double y,
double x)
const;
81 const auto i = (y - _ystart) / _dy;
82 const auto j = (x - _xstart) / _dx;
83 return (i >= 0.0) && (i < _data.
length()) &&
84 (j >= 0.0) && (j < _data.
width());
89 bool _haveData, _boundsError;
92 double _xstart, _ystart, _dx, _dy;
109 template <typename U = T, std::enable_if_t<!std::is_compound<U>::value,
int> = 0>
112 bool equal = _data.
width() == other.width();
113 equal *= _data.
length() == other.length();
126 for (
size_t i = 0; i < _data.
length(); ++i) {
127 for (
size_t j = 0; j < _data.
width(); ++j) {
135 template <typename U = T, std::enable_if_t<std::is_compound<U>::value,
int> = 0>
138 bool equal = _data.
width() == other.width();
139 equal *= _data.
length() == other.length();
152 for (
size_t i = 0; i < _data.
length(); ++i) {
153 for (
size_t j = 0; j < _data.
width(); ++j) {
162 template <
typename T>
164 LUT2d() : _haveData(false), _boundsError(true), _refValue(0.0) {
165 _setInterpolator(isce3::core::BILINEAR_METHOD);
169 template <
typename T>
173 _setInterpolator(method);
177 template <
typename T>
180 _boundsError(lut.boundsError()),
181 _refValue(lut.refValue()),
182 _xstart(lut.xStart()), _ystart(lut.yStart()),
183 _dx(lut.xSpacing()), _dy(lut.ySpacing()),
185 _setInterpolator(lut.interpMethod());
189 template <
typename T>
193 _refValue = lut.refValue();
194 _xstart = lut.xStart();
195 _ystart = lut.yStart();
196 _dx = lut.xSpacing();
197 _dy = lut.ySpacing();
199 _haveData = lut.haveData();
200 _boundsError = lut.boundsError();
201 _setInterpolator(lut.interpMethod());
bool contains(double y, double x) const
Check if point resides in domain of LUT.
Definition: LUT2d.h:79
bool compareComplex(std::complex< T > first, std::complex< T > second)
Tolerance test for real and imaginary components for scalar value.
Definition: Utilities.h:237
dataInterpMethod
Enumeration type to indicate interpolation method.
Definition: Constants.h:23
bool compareFloatingPoint(T first, T second)
Combined absolute and relative tolerance test.
Definition: Utilities.h:226
Data structure to store 2D Lookup table.
Definition: forward.h:30
void setFromData(const std::valarray< double > &xcoord, const std::valarray< double > &ycoord, const isce3::core::Matrix< T > &data)
Definition: LUT2d.cpp:52
size_t length() const
Get matrix length.
Definition: Matrix.icc:357
T eval(double y, double x) const
Definition: LUT2d.cpp:108
size_t width() const
Get matrix width.
Definition: Matrix.icc:349