LSST Applications 27.0.0,g0265f82a02+469cd937ee,g02d81e74bb+21ad69e7e1,g1470d8bcf6+cbe83ee85a,g2079a07aa2+e67c6346a6,g212a7c68fe+04a9158687,g2305ad1205+94392ce272,g295015adf3+81dd352a9d,g2bbee38e9b+469cd937ee,g337abbeb29+469cd937ee,g3939d97d7f+72a9f7b576,g487adcacf7+71499e7cba,g50ff169b8f+5929b3527e,g52b1c1532d+a6fc98d2e7,g591dd9f2cf+df404f777f,g5a732f18d5+be83d3ecdb,g64a986408d+21ad69e7e1,g858d7b2824+21ad69e7e1,g8a8a8dda67+a6fc98d2e7,g99cad8db69+f62e5b0af5,g9ddcbc5298+d4bad12328,ga1e77700b3+9c366c4306,ga8c6da7877+71e4819109,gb0e22166c9+25ba2f69a1,gb6a65358fc+469cd937ee,gbb8dafda3b+69d3c0e320,gc07e1c2157+a98bf949bb,gc120e1dc64+615ec43309,gc28159a63d+469cd937ee,gcf0d15dbbd+72a9f7b576,gdaeeff99f8+a38ce5ea23,ge6526c86ff+3a7c1ac5f1,ge79ae78c31+469cd937ee,gee10cc3b42+a6fc98d2e7,gf1cff7945b+21ad69e7e1,gfbcc870c63+9a11dc8c8f
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Attributes | Friends | List of all members
lsst::afw::math::InterpolateGsl Class Reference
Inheritance diagram for lsst::afw::math::InterpolateGsl:
lsst::afw::math::Interpolate

Public Types

enum  Style {
  UNKNOWN = -1 , CONSTANT = 0 , LINEAR = 1 , NATURAL_SPLINE = 2 ,
  CUBIC_SPLINE = 3 , CUBIC_SPLINE_PERIODIC = 4 , AKIMA_SPLINE = 5 , AKIMA_SPLINE_PERIODIC = 6 ,
  NUM_STYLES
}
 

Public Member Functions

 ~InterpolateGsl () override
 
double interpolate (double const x) const override
 
std::vector< double > interpolate (std::vector< double > const &x) const
 
ndarray::Array< double, 1 > interpolate (ndarray::Array< double const, 1 > const &x) const
 

Protected Attributes

std::vector< double > const _x
 
std::vector< double > const _y
 
Interpolate::Style const _style
 

Friends

std::shared_ptr< InterpolatemakeInterpolate (std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style)
 A factory function to make Interpolate objects.
 

Detailed Description

Definition at line 170 of file Interpolate.cc.

Member Enumeration Documentation

◆ Style

Enumerator
UNKNOWN 
CONSTANT 
LINEAR 
NATURAL_SPLINE 
CUBIC_SPLINE 
CUBIC_SPLINE_PERIODIC 
AKIMA_SPLINE 
AKIMA_SPLINE_PERIODIC 
NUM_STYLES 

Definition at line 38 of file Interpolate.h.

Constructor & Destructor Documentation

◆ ~InterpolateGsl()

lsst::afw::math::InterpolateGsl::~InterpolateGsl ( )
override

Definition at line 219 of file Interpolate.cc.

219 {
220 ::gsl_interp_free(_interp);
221 ::gsl_interp_accel_free(_acc);
222}

Member Function Documentation

◆ interpolate() [1/3]

double lsst::afw::math::InterpolateGsl::interpolate ( double const x) const
overridevirtual

Implements lsst::afw::math::Interpolate.

Definition at line 224 of file Interpolate.cc.

224 {
225 // NaNs interpolate as NaN.
226 if (std::isnan(xInterp)) {
227 return xInterp;
228 }
229 // New GSL versions refuse to extrapolate.
230 // gsl_interp_init() requires x to be ordered, so can just check
231 // the array endpoints for out-of-bounds.
232 if ((xInterp < _x.front() || (xInterp > _x.back()))) {
233 // do our own quadratic extrapolation.
234 // (GSL only provides first and second derivative functions)
235 /* could also just fail via:
236 throw LSST_EXCEPT(
237 pex::exceptions::InvalidParameterError,
238 (boost::format("Interpolation point %f outside range [%f, %f]")
239 % x % _x.front() % _x.back()).str()
240 );
241 */
242 double x0, y0;
243 if (xInterp < _x.front()) {
244 x0 = _x.front();
245 y0 = _y.front();
246 } else {
247 x0 = _x.back();
248 y0 = _y.back();
249 }
250 // first derivative at endpoint
251 double d = ::gsl_interp_eval_deriv(_interp, &_x[0], &_y[0], x0, _acc);
252 // second derivative at endpoint
253 double d2 = ::gsl_interp_eval_deriv2(_interp, &_x[0], &_y[0], x0, _acc);
254 return y0 + (xInterp - x0) * d + 0.5 * (xInterp - x0) * (xInterp - x0) * d2;
255 }
256 assert(xInterp >= _x.front());
257 assert(xInterp <= _x.back());
258 return ::gsl_interp_eval(_interp, &_x[0], &_y[0], xInterp, _acc);
259}
T back(T... args)
std::vector< double > const _x
Definition Interpolate.h:83
std::vector< double > const _y
Definition Interpolate.h:84
T front(T... args)
T isnan(T... args)

◆ interpolate() [2/3]

ndarray::Array< double, 1 > lsst::afw::math::Interpolate::interpolate ( ndarray::Array< double const, 1 > const & x) const
inherited

Definition at line 308 of file Interpolate.cc.

308 {
309 int const num = x.getShape()[0];
310 ndarray::Array<double, 1> out = ndarray::allocate(ndarray::makeVector(num));
311 for (int i = 0; i < num; ++i) {
312 std::cout << "Interpolating " << x[i] << std::endl;
313 out[i] = interpolate(x[i]);
314 }
315 return out;
316}
virtual double interpolate(double const x) const =0
T endl(T... args)

◆ interpolate() [3/3]

std::vector< double > lsst::afw::math::Interpolate::interpolate ( std::vector< double > const & x) const
inherited

Definition at line 299 of file Interpolate.cc.

299 {
300 size_t const num = x.size();
301 std::vector<double> out(num);
302 for (size_t i = 0; i < num; ++i) {
303 out[i] = interpolate(x[i]);
304 }
305 return out;
306}

Friends And Related Symbol Documentation

◆ makeInterpolate

std::shared_ptr< Interpolate > makeInterpolate ( std::vector< double > const & x,
std::vector< double > const & y,
Interpolate::Style const style = Interpolate::AKIMA_SPLINE )
friend

A factory function to make Interpolate objects.

Parameters
xthe x-values of points
ythe values at x[]
styledesired interpolator

Definition at line 347 of file Interpolate.cc.

348 {
349 switch (style) {
351 return std::shared_ptr<Interpolate>(new InterpolateConstant(x, y, style));
352 default: // use GSL
353 return std::shared_ptr<Interpolate>(new InterpolateGsl(x, y, style));
354 }
355}
int y
Definition SpanSet.cc:48

Member Data Documentation

◆ _style

Interpolate::Style const lsst::afw::math::Interpolate::_style
protectedinherited

Definition at line 85 of file Interpolate.h.

◆ _x

std::vector<double> const lsst::afw::math::Interpolate::_x
protectedinherited

Definition at line 83 of file Interpolate.h.

◆ _y

std::vector<double> const lsst::afw::math::Interpolate::_y
protectedinherited

Definition at line 84 of file Interpolate.h.


The documentation for this class was generated from the following file: