LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
Interpolate.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
30 #include <limits>
31 #include <algorithm>
32 #include <map>
33 #include "boost/format.hpp"
34 #include "boost/shared_ptr.hpp"
35 #include "gsl/gsl_errno.h"
36 #include "gsl/gsl_interp.h"
37 #include "gsl/gsl_spline.h"
38 #include "lsst/pex/exceptions.h"
40 
41 namespace lsst {
42 namespace afw {
43 namespace math {
44 
45 /************************************************************************************************************/
46 
47 namespace {
48  std::pair<std::vector<double>, std::vector<double> >
49  recenter(std::vector<double> const &x,
50  std::vector<double> const &y)
51  {
52  if (x.size() != y.size()) {
53  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
54  str(boost::format("Dimensions of x and y must match; %ul != %ul")
55  % x.size() % y.size()));
56  }
57  std::size_t const len = x.size();
58  if (len == 0) {
59  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
60  "You must provide at least 1 point");
61  } else if (len == 1) {
62  return std::make_pair(x, y);
63  }
64 
65  std::vector<double> recentered_x(len + 1);
66  std::vector<double> recentered_y(len + 1);
67 
68  recentered_x[0] = 0.5*(3*x[0] - x[1]);
69  recentered_y[0] = y[0];
70 
71  for (std::size_t i = 0, j = 1; i < len - 1; ++i, ++j) {
72  recentered_x[j] = 0.5*(x[i] + x[i + 1]);
73  recentered_y[j] = 0.5*(y[i] + y[i + 1]);
74  }
75  recentered_x[len] = 0.5*(3*x[len - 1] - x[len - 2]);
76  recentered_y[len] = y[len - 1];
77 
78  return std::make_pair(recentered_x, recentered_y);
79  }
80 }
81 
83  friend PTR(Interpolate) makeInterpolate(std::vector<double> const &x, std::vector<double> const &y,
84  Interpolate::Style const style);
85 public:
86  virtual ~InterpolateConstant() {}
87  virtual double interpolate(double const x) const;
88 private:
89  InterpolateConstant(std::vector<double> const &x,
90  std::vector<double> const &y,
91  Interpolate::Style const style
92  ) :
93  Interpolate(recenter(x, y)), _old(_x.begin()) {}
94  mutable std::vector<double>::const_iterator _old; // last position we found xInterp at
95 };
96 
97 
99 double InterpolateConstant::interpolate(double const xInterp // the value we want to interpolate to
100  ) const
101 {
102  //
103  // Look for the interval wherein lies xInterp. We could naively use std::upper_bound, but that requires a
104  // logarithmic time lookup so we'll cache the previous answer in _old -- this is a good idea if people
105  // usually call this routine repeatedly for a range of x
106  //
107  // We start by searching up from _old
108  //
109  if (xInterp < *_old) { // We're to the left of the cache
110  if (_old == _x.begin()) { // ... actually off the array
111  return _y[0];
112  }
113  _old = _x.begin(); // reset the cached point to the start of the array
114  } else { // see if we're still in the same interval
115  if (_old < _x.end() - 1 and xInterp < *(_old + 1)) { // we are, so we're done
116  return _y[_old - _x.begin()];
117  }
118  }
119  // We're to the right of the cached point and not in the same inverval, so search up from _old
120  std::vector<double>::const_iterator low = std::upper_bound(_old, _x.end(), xInterp);
121  //
122  // Did that work?
123  if (low == _old && _old != _x.begin()) {
124  // No. Sigh. Search the entire range.
125  low = std::upper_bound(_x.begin(), low + 1, xInterp);
126  }
127  //
128  // OK, we've found the right interval. Return the desired value, being careful at the ends
129  //
130  if (low == _x.end()) {
131  return _y[_y.size() - 1];
132  } else if (low == _x.begin()) {
133  return _y[0];
134  } else {
135  --low;
136  _old = low;
137  return _y[low - _x.begin()];
138  }
139 }
140 
141 /************************************************************************************************************/
142 namespace {
143 /*
144  * Conversion function to switch an Interpolate::Style to a gsl_interp_type.
145  */
146 ::gsl_interp_type const *
147 styleToGslInterpType(Interpolate::Style const style)
148 {
149  switch (style) {
151  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError, "CONSTANT interpolation not supported.");
152  case Interpolate::LINEAR:
153  return ::gsl_interp_linear;
155  return ::gsl_interp_cspline;
157  return ::gsl_interp_cspline;
159  return ::gsl_interp_cspline_periodic;
161  return ::gsl_interp_akima;
163  return ::gsl_interp_akima_periodic;
165  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
166  "I am unable to make an interpolator of type UNKNOWN");
168  throw LSST_EXCEPT(pex::exceptions::LogicError,
169  str(boost::format("You can't get here: style == %") % style));
170  }
171 }
172 }
173 
174 class InterpolateGsl : public Interpolate {
175  friend PTR(Interpolate) makeInterpolate(std::vector<double> const &x, std::vector<double> const &y,
176  Interpolate::Style const style);
177 public:
178  virtual ~InterpolateGsl();
179  virtual double interpolate(double const x) const;
180 private:
181  InterpolateGsl(std::vector<double> const &x, std::vector<double> const &y, Interpolate::Style const style);
182 
183  ::gsl_interp_type const *_interpType;
184  ::gsl_interp_accel *_acc;
185  ::gsl_interp *_interp;
186 };
187 
188 InterpolateGsl::InterpolateGsl(std::vector<double> const &x,
189  std::vector<double> const &y,
190  Interpolate::Style const style
191  ) :
192  Interpolate(x, y), _interpType(styleToGslInterpType(style))
193 {
194  _acc = ::gsl_interp_accel_alloc();
195  if (!_acc) {
196  throw LSST_EXCEPT(pex::exceptions::MemoryError, "gsl_interp_accel_alloc failed");
197  }
198 
199  _interp = ::gsl_interp_alloc(_interpType, _y.size());
200  if (!_interp) {
201  throw LSST_EXCEPT(pex::exceptions::OutOfRangeError,
202  str(boost::format("Failed to initialise spline for type %s, length %d")
203  % _interpType->name % _y.size()));
204 
205  }
206  // Note, "x" and "y" are vector<double>; gsl_inter_init requires double[].
207  // The &(x[0]) here is valid because std::vector guarantees that the values are
208  // stored contiguously in memory (for types other than bool); C++0X 23.3.6.1 for
209  // those of you reading along.
210  int const status = ::gsl_interp_init(_interp, &x[0], &y[0], _y.size());
211  if (status != 0) {
212  throw LSST_EXCEPT(pex::exceptions::RuntimeError,
213  str(boost::format("gsl_interp_init failed: %s [%d]")
214  % ::gsl_strerror(status) % status));
215  }
216 }
217 
219  ::gsl_interp_free(_interp);
220  ::gsl_interp_accel_free(_acc);
221 }
222 
223 double InterpolateGsl::interpolate(double const xInterp) const
224 {
225  // New GSL versions refuse to extrapolate.
226  // gsl_interp_init() requires x to be ordered, so can just check
227  // the array endpoints for out-of-bounds.
228  if ((xInterp < _x.front() || (xInterp > _x.back()))) {
229  // do our own quadratic extrapolation.
230  // (GSL only provides first and second derivative functions)
231  /* could also just fail via:
232  throw LSST_EXCEPT(
233  pex::exceptions::InvalidParameterError,
234  (boost::format("Interpolation point %f outside range [%f, %f]")
235  % x % _x.front() % _x.back()).str()
236  );
237  */
238  double x0, y0;
239  if (xInterp < _x.front()) {
240  x0 = _x.front();
241  y0 = _y.front();
242  } else {
243  x0 = _x.back();
244  y0 = _y.back();
245  }
246  // first derivative at endpoint
247  double d = ::gsl_interp_eval_deriv(_interp, &_x[0], &_y[0], x0, _acc);
248  // second derivative at endpoint
249  double d2 = ::gsl_interp_eval_deriv2(_interp, &_x[0], &_y[0], x0, _acc);
250  return y0 + (xInterp - x0)*d + 0.5*(xInterp - x0)*(xInterp - x0)*d2;
251  }
252  assert(xInterp >= _x.front());
253  assert(xInterp <= _x.back());
254  return ::gsl_interp_eval(_interp, &_x[0], &_y[0], xInterp, _acc);
255 }
256 
257 /************************************************************************************************************/
262 Interpolate::Style stringToInterpStyle(std::string const &style
263  )
264 {
265  static std::map<std::string, Interpolate::Style> gslInterpTypeStrings;
266  if (gslInterpTypeStrings.empty()) {
267  gslInterpTypeStrings["CONSTANT"] = Interpolate::CONSTANT;
268  gslInterpTypeStrings["LINEAR"] = Interpolate::LINEAR;
269  gslInterpTypeStrings["CUBIC_SPLINE"] = Interpolate::CUBIC_SPLINE;
270  gslInterpTypeStrings["NATURAL_SPLINE"] = Interpolate::NATURAL_SPLINE;
271  gslInterpTypeStrings["CUBIC_SPLINE_PERIODIC"] = Interpolate::CUBIC_SPLINE_PERIODIC;
272  gslInterpTypeStrings["AKIMA_SPLINE"] = Interpolate::AKIMA_SPLINE;
273  gslInterpTypeStrings["AKIMA_SPLINE_PERIODIC"] = Interpolate::AKIMA_SPLINE_PERIODIC;
274  }
275 
276  if ( gslInterpTypeStrings.find(style) == gslInterpTypeStrings.end()) {
277  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError, "Interp style not found: "+style);
278  }
279  return gslInterpTypeStrings[style];
280 }
281 
286  ) {
287  if (n < 1) {
288  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError, "n must be greater than 0");
289  } else if (n > 4) {
291  } else {
292  static std::vector<Interpolate::Style> styles;
293  if (styles.empty()) {
294  styles.resize(5);
295 
296  styles[0] = Interpolate::UNKNOWN; // impossible to reach as we check for n < 1
297  styles[1] = Interpolate::CONSTANT;
298  styles[2] = Interpolate::LINEAR;
299  styles[3] = Interpolate::CUBIC_SPLINE;
300  styles[4] = Interpolate::CUBIC_SPLINE;
301  }
302  return styles[n];
303  }
304 }
305 
310  ) {
311  static std::vector<int> minPoints;
312  if (minPoints.empty()) {
313  minPoints.resize(Interpolate::NUM_STYLES);
314  minPoints[Interpolate::CONSTANT] = 1;
315  minPoints[Interpolate::LINEAR] = 2;
316  minPoints[Interpolate::NATURAL_SPLINE] = 3;
317  minPoints[Interpolate::CUBIC_SPLINE] = 3;
318  minPoints[Interpolate::CUBIC_SPLINE_PERIODIC] = 3;
319  minPoints[Interpolate::AKIMA_SPLINE] = 5;
320  minPoints[Interpolate::AKIMA_SPLINE_PERIODIC] = 5;
321  }
322 
323  if (style >= 0 && style < Interpolate::NUM_STYLES) {
324  return minPoints[style];
325  } else {
326  throw LSST_EXCEPT(pex::exceptions::OutOfRangeError,
327  str(boost::format("Style %d is out of range 0..%d")
328  % style % (Interpolate::NUM_STYLES - 1)));
329  }
330 }
331 
332 /************************************************************************************************************/
342  std::pair<std::vector<double>, std::vector<double> > const xy,
343  Interpolate::Style const style
345  ) : _x(xy.first), _y(xy.second), _style(style)
346 {
347  ;
348 }
349 
353 PTR(Interpolate) makeInterpolate(std::vector<double> const &x,
354  std::vector<double> const &y,
355  Interpolate::Style const style
356  )
357 {
358  switch (style) {
360  return PTR(Interpolate)(new InterpolateConstant(x, y, style));
361  default: // use GSL
362  return PTR(Interpolate)(new InterpolateGsl(x, y, style));
363  }
364 }
365 
366 }}}
int y
Interpolate::Style stringToInterpStyle(std::string const &style)
Conversion function to switch a string to an Interpolate::Style.
Definition: Interpolate.cc:262
Interpolate values for a set of x,y vector&lt;&gt;s.
Definition: Interpolate.h:37
std::vector< double >::const_iterator _old
Definition: Interpolate.cc:94
virtual double interpolate(double const x) const
Interpolate a constant to the point xInterp.
Definition: Interpolate.cc:99
friend boost::shared_ptr< Interpolate > makeInterpolate(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style)
Include files required for standard LSST Exception handling.
boost::shared_ptr< Interpolate > makeInterpolate(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style=Interpolate::AKIMA_SPLINE)
Definition: Interpolate.cc:353
virtual double interpolate(double const x) const
Definition: Interpolate.cc:223
#define PTR(...)
Definition: base.h:41
int const x0
Definition: saturated.cc:45
::gsl_interp_type const * _interpType
Definition: Interpolate.cc:183
::gsl_interp_accel * _acc
Definition: Interpolate.cc:184
double x
std::vector< double > const _x
Definition: Interpolate.h:67
InterpolateConstant(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style)
Definition: Interpolate.cc:89
InterpolateGsl(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style)
Definition: Interpolate.cc:188
Interpolate::Style lookupMaxInterpStyle(int const n)
Get the highest order Interpolation::Style available for &#39;n&#39; points.
Definition: Interpolate.cc:285
int lookupMinInterpPoints(Interpolate::Style const style)
Get the minimum number of points needed to use the requested interpolation style. ...
Definition: Interpolate.cc:309
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
Interpolate(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style=UNKNOWN)
Definition: Interpolate.h:60
friend boost::shared_ptr< Interpolate > makeInterpolate(std::vector< double > const &x, std::vector< double > const &y, Interpolate::Style const style)
int status
std::vector< double > const _y
Definition: Interpolate.h:68
int const y0
Definition: saturated.cc:45