All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
interp.h
Go to the documentation of this file.
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 #ifndef O2SCL_INTERP_H
24 #define O2SCL_INTERP_H
25 
26 /** \file interp.h
27  \brief One-dimensional interpolation classes and interpolation types
28 */
29 
30 #include <iostream>
31 #include <string>
32 
33 #include <boost/numeric/ublas/vector.hpp>
34 #include <boost/numeric/ublas/vector_proxy.hpp>
35 
36 #include <o2scl/search_vec.h>
37 #include <o2scl/tridiag.h>
38 #include <o2scl/funct.h>
39 
40 #ifndef DOXYGEN_NO_O2NS
41 namespace o2scl {
42 #endif
43 
44  /** \brief Interpolation types
45  */
46  enum {
47  /// Linear
49  /// Cubic spline for natural boundary conditions
51  /// Cubic spline for periodic boundary conditions
53  /// Akima spline for natural boundary conditions
55  /// Akima spline for periodic boundary conditions
57  /// Monotonicity-preserving interpolation
59  };
60 
61  /** \brief Base low-level interpolation class [abstract base]
62 
63  See also the \ref intp_section section of the \o2 User's guide.
64 
65  To interpolate a set vectors \c x and \c y, call set() and then
66  the interpolation functions eval(), deriv(), deriv2() and
67  integ(). If the \c x and \c y vectors do not change, then you
68  may call the interpolation functions multiple times in
69  succession. These classes do not copy the user-specified vectors
70  but store pointers to them. Thus, if the vector is changed
71  without a new call to \ref interp_base::set(), the behavior of
72  the interpolation functions is undefined.
73 
74  \comment
75  AWS, 12/27/13: Copy constructors might be ill-advised for
76  this class since we store pointers. For now, we don't
77  allow the user to use them.
78  \endcomment
79  */
80  template<class vec_t, class vec2_t=vec_t> class interp_base {
81  //public funct {
82 
83 #ifdef O2SCL_NEVER_DEFINED
84  }{
85 #endif
86 
87 #ifndef DOXYGEN_INTERNAL
88 
89  protected:
90 
91  /** \brief To perform binary searches
92 
93  This pointer is set to zero in the constructor and should be
94  non-zero only if it has been allocated with \c new.
95  */
97 
98  /// Independent vector
99  const vec_t *px;
100 
101  /// Dependent vector
102  const vec2_t *py;
103 
104  /// Vector size
105  size_t sz;
106 
107  /** \brief An internal function to assist in computing the
108  integral for both the cspline and Akima types
109  */
110  double integ_eval(double ai, double bi, double ci, double di, double xi,
111  double a, double b) const {
112 
113  double r1=a-xi;
114  double r2=b-xi;
115  double r12=r1+r2;
116  double bterm=0.5*bi*r12;
117  double cterm=ci*(r1*r1+r2*r2+r1*r2)/3.0;
118  double dterm=0.25*di*r12*(r1*r1+r2*r2);
119 
120  return (b-a)*(ai+bterm+cterm+dterm);
121  }
122 
123 #endif
124 
125  public:
126 
127  interp_base() {
128  sz=0;
129  }
130 
131  virtual ~interp_base() {
132  }
133 
134  /** \brief The minimum size of the vectors to interpolate between
135 
136  This variable must be set in the constructor of the children
137  for access by the class user.
138  */
139  size_t min_size;
140 
141  /// Initialize interpolation routine
142  virtual void set(size_t size, const vec_t &x, const vec2_t &y)=0;
143 
144  /// Give the value of the function \f$ y(x=x_0) \f$ .
145  virtual double eval(double x0) const=0;
146 
147  /// Give the value of the function \f$ y(x=x_0) \f$ .
148  virtual double operator()(double x0) const {
149  return eval(x0);
150  }
151 
152  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
153  virtual double deriv(double x0) const=0;
154 
155  /** \brief Give the value of the second derivative
156  \f$ y^{\prime \prime}(x=x_0) \f$ .
157  */
158  virtual double deriv2(double x0) const=0;
159 
160  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
161  virtual double integ(double a, double b) const=0;
162 
163  /// Return the type
164  virtual const char *type() const=0;
165 
166 #ifndef DOXYGEN_INTERNAL
167 
168  private:
169 
172 
173 #endif
174 
175  };
176 
177  /** \brief Linear interpolation (GSL)
178 
179  See also the \ref intp_section section of the \o2 User's guide.
180 
181  Linear interpolation requires no calls to allocate() or free()
182  as there is no internal storage required.
183  */
184  template<class vec_t, class vec2_t=vec_t> class interp_linear :
185  public interp_base<vec_t,vec2_t> {
186 
187 #ifdef O2SCL_NEVER_DEFINED
188  }{
189 #endif
190 
191  public:
192 
193  interp_linear() {
194  this->min_size=2;
195  }
196 
197  virtual ~interp_linear() {}
198 
199  /// Initialize interpolation routine
200  virtual void set(size_t size, const vec_t &x, const vec2_t &y) {
201  if (size<this->min_size) {
202  O2SCL_ERR((((std::string)"Vector size, ")+szttos(size)+", is less"+
203  " than "+szttos(this->min_size)+" in interp_linear::"+
204  "set().").c_str(),exc_einval);
205  }
206  this->svx.set_vec(size,x);
207  this->px=&x;
208  this->py=&y;
209  this->sz=size;
210  return;
211  }
212 
213  /// Give the value of the function \f$ y(x=x_0) \f$ .
214  virtual double eval(double x0) const {
215 
216  size_t index=this->svx.find(x0);
217 
218  double x_lo=(*this->px)[index];
219  double x_hi=(*this->px)[index+1];
220  double y_lo=(*this->py)[index];
221  double y_hi=(*this->py)[index+1];
222  double dx=x_hi-x_lo;
223 
224  return y_lo+(x0-x_lo)/dx*(y_hi-y_lo);
225  }
226 
227  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
228  virtual double deriv(double x0) const {
229 
230  size_t index=this->svx.find(x0);
231 
232  double x_lo=(*this->px)[index];
233  double x_hi=(*this->px)[index+1];
234  double y_lo=(*this->py)[index];
235  double y_hi=(*this->py)[index+1];
236  double dx=x_hi-x_lo;
237  double dy=y_hi-y_lo;
238 
239  return dy/dx;
240  }
241 
242  /** \brief Give the value of the second derivative
243  \f$ y^{\prime \prime}(x=x_0) \f$ (always zero)
244  */
245  virtual double deriv2(double x0) const {
246  return 0.0;
247  }
248 
249  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
250  virtual double integ(double a, double b) const {
251 
252  size_t i, index_a, index_b;
253 
254  bool flip=false;
255  if (((*this->px)[0]<(*this->px)[this->sz-1] && a>b) ||
256  ((*this->px)[0]>(*this->px)[this->sz-1] && a<b)) {
257  double tmp=a;
258  a=b;
259  b=tmp;
260  flip=true;
261  }
262 
263  index_a=this->svx.find(a);
264  index_b=this->svx.find(b);
265 
266  double result=0.0;
267  for(i=index_a; i<=index_b; i++) {
268 
269  double x_lo=(*this->px)[i];
270  double x_hi=(*this->px)[i+1];
271  double y_lo=(*this->py)[i];
272  double y_hi=(*this->py)[i+1];
273  double dx=x_hi-x_lo;
274 
275  if(dx != 0.0) {
276 
277  if (i == index_a || i == index_b) {
278  double x1=(i == index_a) ? a : x_lo;
279  double x2=(i == index_b) ? b : x_hi;
280  double D=(y_hi-y_lo)/dx;
281  result += (x2-x1)*(y_lo+0.5*D*((x2-x_lo)+(x1-x_lo)));
282  } else {
283  result += 0.5*dx*(y_lo+y_hi);
284  }
285 
286  } else {
287  std::string str=((std::string)"Interval of length zero ")+
288  "between ("+o2scl::dtos(x_lo)+","+o2scl::dtos(y_lo)+
289  ") and ("+o2scl::dtos(x_hi)+","+o2scl::dtos(y_hi)+
290  " in interp_linear::integ().";
291  O2SCL_ERR(str.c_str(),exc_einval);
292  }
293  }
294 
295  if (flip) result=-result;
296  return result;
297  }
298 
299  /// Return the type, \c "interp_linear".
300  virtual const char *type() const { return "interp_linear"; }
301 
302 #ifndef DOXYGEN_INTERNAL
303 
304  private:
305 
308 
309 #endif
310 
311  };
312 
313  /** \brief Cubic spline interpolation (GSL)
314 
315  See also the \ref intp_section section of the \o2 User's guide.
316 
317  By default, this class uses natural boundary conditions, where
318  the second derivative vanishes at each end point. Extrapolation
319  effectively assumes that the second derivative is linear outside
320  of the endpoints.
321  */
322  template<class vec_t, class vec2_t=vec_t> class interp_cspline :
323  public interp_base<vec_t,vec2_t> {
324 
325 #ifdef O2SCL_NEVER_DEFINED
326  }{
327 #endif
328 
329  public:
330 
331  typedef boost::numeric::ublas::vector<double> ubvector;
332  typedef boost::numeric::ublas::vector_slice<ubvector> ubvector_slice;
333  typedef boost::numeric::ublas::vector_range<ubvector> ubvector_range;
334  typedef boost::numeric::ublas::slice slice;
335  typedef boost::numeric::ublas::range range;
336 
337 #ifndef DOXYGEN_INTERNAL
338 
339  protected:
340 
341  /// \name Storage for cubic spline interpolation
342  //@{
343  ubvector c;
344  ubvector g;
345  ubvector diag;
346  ubvector offdiag;
347  //@}
348 
349  /// Memory for the tridiagonalization
351 
352  /// Compute coefficients for cubic spline interpolation
353  void coeff_calc(const ubvector &c_array, double dy, double dx,
354  size_t index, double &b, double &c2, double &d) const {
355 
356  double c_i=c_array[index];
357  double c_ip1=c_array[index+1];
358  b=(dy/dx)-dx*(c_ip1+2.0*c_i)/3.0;
359  c2=c_i;
360  d=(c_ip1-c_i)/(3.0*dx);
361 
362  return;
363  }
364 
365 #endif
366 
367  public:
368 
369  /** \brief Create a base interpolation object with natural or
370  periodic boundary conditions
371  */
373  this->min_size=3;
374  }
375 
376  virtual ~interp_cspline() {
377  }
378 
379  /** \brief Initialize interpolation routine
380  */
381  virtual void set(size_t size, const vec_t &xa, const vec2_t &ya) {
382 
383  if (size<this->min_size) {
384  O2SCL_ERR((((std::string)"Vector size, ")+szttos(size)+", is less"+
385  " than "+szttos(this->min_size)+" in interp_cspline::"+
386  "set().").c_str(),exc_einval);
387  }
388 
389  if (size!=this->sz) {
390  c.resize(size);
391  g.resize(size);
392  diag.resize(size);
393  offdiag.resize(size);
394  p4m.resize(size);
395  }
396 
397  this->px=&xa;
398  this->py=&ya;
399  this->sz=size;
400 
401  this->svx.set_vec(size,xa);
402 
403  /// Natural boundary conditions
404 
405  size_t i;
406  size_t num_points=size;
407  size_t max_index=num_points-1;
408  size_t sys_size=max_index-1;
409 
410  c[0]=0.0;
411  c[max_index]=0.0;
412 
413  for (i=0; i < sys_size; i++) {
414  double h_i=xa[i+1]-xa[i];
415  double h_ip1=xa[i+2]-xa[i+1];
416  double ydiff_i=ya[i+1]-ya[i];
417  double ydiff_ip1=ya[i+2]-ya[i+1];
418  double g_i=(h_i != 0.0) ? 1.0/h_i : 0.0;
419  double g_ip1=(h_ip1 != 0.0) ? 1.0/h_ip1 : 0.0;
420  offdiag[i]=h_ip1;
421  diag[i]=2.0*(h_ip1+h_i);
422  g[i]=3.0*(ydiff_ip1*g_ip1-ydiff_i*g_i);
423  }
424 
425  if (sys_size == 1) {
426 
427  c[1]=g[0]/diag[0];
428 
429  return;
430  }
431 
432  ubvector_range cp1(c,range(1,c.size()));
433  o2scl_linalg::solve_tridiag_sym<ubvector,ubvector,ubvector,
434  ubvector_range,o2scl_linalg::ubvector_4_mem,ubvector>
435  (diag,offdiag,g,cp1,sys_size,p4m);
436 
437  return;
438  }
439 
440  /// Give the value of the function \f$ y(x=x_0) \f$ .
441  virtual double eval(double x0) const {
442 
443  size_t index=this->svx.find(x0);
444 
445  double x_lo=(*this->px)[index];
446  double x_hi=(*this->px)[index+1];
447  double dx=x_hi-x_lo;
448 
449  double y_lo=(*this->py)[index];
450  double y_hi=(*this->py)[index+1];
451  double dy=y_hi-y_lo;
452  double delx=x0-x_lo;
453  double b_i, c_i, d_i;
454 
455  coeff_calc(c,dy,dx,index,b_i,c_i,d_i);
456 
457  return y_lo+delx*(b_i+delx*(c_i+delx*d_i));
458  }
459 
460  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
461  virtual double deriv(double x0) const {
462 
463  size_t index=this->svx.find(x0);
464 
465  double x_lo=(*this->px)[index];
466  double x_hi=(*this->px)[index+1];
467  double dx=x_hi-x_lo;
468 
469  double y_lo=(*this->py)[index];
470  double y_hi=(*this->py)[index+1];
471  double dy=y_hi-y_lo;
472  double delx=x0-x_lo;
473  double b_i, c_i, d_i;
474 
475  coeff_calc(c,dy,dx,index,b_i,c_i,d_i);
476 
477  return b_i+delx*(2.0*c_i+3.0*d_i*delx);
478  }
479 
480  /** \brief Give the value of the second derivative
481  \f$ y^{\prime \prime}(x=x_0) \f$ .
482  */
483  virtual double deriv2(double x0) const {
484 
485  size_t index=this->svx.find(x0);
486 
487  double x_lo=(*this->px)[index];
488  double x_hi=(*this->px)[index+1];
489  double dx=x_hi-x_lo;
490 
491  double y_lo=(*this->py)[index];
492  double y_hi=(*this->py)[index+1];
493  double dy=y_hi-y_lo;
494  double delx=x0-x_lo;
495  double b_i, c_i, d_i;
496 
497  coeff_calc(c,dy,dx,index,b_i,c_i,d_i);
498 
499  return 2.0*c_i+6.0*d_i*delx;
500  }
501 
502  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
503  virtual double integ(double a, double b) const {
504 
505  size_t i, index_a, index_b;
506 
507  bool flip=false;
508  if (((*this->px)[0]<(*this->px)[this->sz-1] && a>b) ||
509  ((*this->px)[0]>(*this->px)[this->sz-1] && a<b)) {
510  double tmp=a;
511  a=b;
512  b=tmp;
513  flip=true;
514  }
515 
516  index_a=this->svx.find(a);
517  index_b=this->svx.find(b);
518 
519  double result=0.0;
520 
521  for(i=index_a; i<=index_b; i++) {
522 
523  double x_lo=(*this->px)[i];
524  double x_hi=(*this->px)[i+1];
525  double y_lo=(*this->py)[i];
526  double y_hi=(*this->py)[i+1];
527  double dx=x_hi-x_lo;
528  double dy=y_hi-y_lo;
529 
530  if(dx != 0.0) {
531  double b_i, c_i, d_i;
532  coeff_calc(c,dy,dx,i,b_i,c_i,d_i);
533  if (i == index_a || i == index_b) {
534  double x1=(i == index_a) ? a : x_lo;
535  double x2=(i == index_b) ? b : x_hi;
536  result += this->integ_eval(y_lo,b_i,c_i,d_i,x_lo,x1,x2);
537  } else {
538  result += dx*(y_lo+dx*(0.5*b_i+
539  dx*(c_i/3.0+0.25*d_i*dx)));
540  }
541  } else {
542  std::string str=((std::string)"Interval of length zero ")+
543  "between ("+o2scl::dtos(x_lo)+","+o2scl::dtos(y_lo)+
544  ") and ("+o2scl::dtos(x_hi)+","+o2scl::dtos(y_hi)+
545  " in interp_cspline::integ().";
546  O2SCL_ERR(str.c_str(),exc_einval);
547  }
548 
549  }
550 
551  if (flip) result*=-1.0;
552 
553  return result;
554  }
555 
556  /// Return the type, \c "interp_cspline".
557  virtual const char *type() const { return "interp_cspline"; }
558 
559 #ifndef DOXYGEN_INTERNAL
560 
561  private:
562 
565  (const interp_cspline<vec_t,vec2_t>&);
566 
567 #endif
568 
569  };
570 
571  /** \brief Cubic spline interpolation with periodic
572  boundary conditions (GSL)
573 
574  See also the \ref intp_section section of the \o2 User's guide.
575  */
576  template<class vec_t, class vec2_t=vec_t> class interp_cspline_peri :
577  public interp_cspline<vec_t,vec2_t> {
578 
579 #ifdef O2SCL_NEVER_DEFINED
580  }{
581 #endif
582 
583  protected:
584 
585  /// Memory for the tridiagonalization
587 
588  public:
589 
591  typedef boost::numeric::ublas::vector_slice<ubvector> ubvector_slice;
592  typedef boost::numeric::ublas::vector_range<ubvector> ubvector_range;
593  typedef boost::numeric::ublas::slice slice;
594  typedef boost::numeric::ublas::range range;
595 
596  interp_cspline_peri() : interp_cspline<vec_t,vec2_t>() {
597  this->min_size=2;
598  }
599 
600  virtual ~interp_cspline_peri() {
601  }
602 
603  /// Return the type, \c "interp_cspline_peri".
604  virtual const char *type() const { return "interp_cspline_peri"; }
605 
606  /** \brief Initialize interpolation routine
607  */
608  virtual void set(size_t size, const vec_t &xa, const vec2_t &ya) {
609 
610  if (size<this->min_size) {
611  O2SCL_ERR((((std::string)"Vector size, ")+szttos(size)+", is less"+
612  " than "+szttos(this->min_size)+" in interp_cspline"+
613  "_peri::set().").c_str(),exc_einval);
614  }
615 
616  if (size!=this->sz) {
617  this->c.resize(size);
618  this->g.resize(size);
619  this->diag.resize(size);
620  this->offdiag.resize(size);
621  p5m.resize(size);
622  }
623 
624  this->px=&xa;
625  this->py=&ya;
626  this->sz=size;
627 
628  this->svx.set_vec(size,xa);
629 
630  /// Periodic boundary conditions
631 
632  size_t i;
633  size_t num_points=size;
634  // Engeln-Mullges+Uhlig "n"
635  size_t max_index=num_points-1;
636  // linear system is sys_size x sys_size
637  size_t sys_size=max_index;
638 
639  if (sys_size == 2) {
640 
641  // solve 2x2 system
642 
643  double h0=xa[1]-xa[0];
644  double h1=xa[2]-xa[1];
645  double h2=xa[3]-xa[2];
646  double A=2.0*(h0+h1);
647  double B=h0+h1;
648  double gx[2];
649  double det;
650 
651  gx[0]=3.0*((ya[2]-ya[1])/h1-(ya[1]-ya[0])/h0);
652  gx[1]=3.0*((ya[1]-ya[2])/h2-(ya[2]-ya[1])/h1);
653 
654  det=3.0*(h0+h1)*(h0+h1);
655  this->c[1]=( A*gx[0]-B*gx[1])/det;
656  this->c[2]=(-B*gx[0]+A*gx[1])/det;
657  this->c[0]=this->c[2];
658 
659  return;
660 
661  } else {
662 
663  for (i=0; i < sys_size-1; i++) {
664  double h_i=xa[i+1]-xa[i];
665  double h_ip1=xa[i+2]-xa[i+1];
666  double ydiff_i=ya[i+1]-ya[i];
667  double ydiff_ip1=ya[i+2]-ya[i+1];
668  double g_i=(h_i != 0.0) ? 1.0/h_i : 0.0;
669  double g_ip1=(h_ip1 != 0.0) ? 1.0/h_ip1 : 0.0;
670  this->offdiag[i]=h_ip1;
671  this->diag[i]=2.0*(h_ip1+h_i);
672  this->g[i]=3.0*(ydiff_ip1*g_ip1-ydiff_i*g_i);
673  }
674 
675  i=sys_size-1;
676  {
677  double h_i=xa[i+1]-xa[i];
678  double h_ip1=xa[1]-xa[0];
679  double ydiff_i=ya[i+1]-ya[i];
680  double ydiff_ip1=ya[1]-ya[0];
681  double g_i=(h_i != 0.0) ? 1.0/h_i : 0.0;
682  double g_ip1=(h_ip1 != 0.0) ? 1.0/h_ip1 : 0.0;
683  this->offdiag[i]=h_ip1;
684  this->diag[i]=2.0*(h_ip1+h_i);
685  this->g[i]=3.0*(ydiff_ip1*g_ip1-ydiff_i*g_i);
686  }
687 
688  ubvector_range cp1(this->c,range(1,this->c.size()));
689  o2scl_linalg::solve_cyc_tridiag_sym<ubvector,ubvector,ubvector,
690  ubvector_range,o2scl_linalg::ubvector_5_mem,ubvector>
691  (this->diag,this->offdiag,this->g,cp1,sys_size,p5m);
692  this->c[0]=this->c[max_index];
693 
694  return;
695  }
696 
697  }
698 
699 #ifndef DOXYGEN_INTERNAL
700 
701  private:
702 
706  (const interp_cspline_peri<vec_t,vec2_t>&);
707 
708 #endif
709 
710  };
711 
712  /** \brief Akima spline interpolation (GSL)
713 
714  See also the \ref intp_section section of the \o2 User's guide.
715 
716  This class uses natural boundary conditions, where the second
717  derivative vanishes at each end point. Extrapolation effectively
718  assumes that the second derivative is linear outside of the
719  endpoints. Use \ref interp_akima_peri for perodic boundary
720  conditions.
721  */
722  template<class vec_t, class vec2_t=vec_t> class interp_akima :
723  public interp_base<vec_t,vec2_t> {
724 
725  public:
726 
728  typedef boost::numeric::ublas::vector_slice<ubvector> ubvector_slice;
729  typedef boost::numeric::ublas::vector_range<ubvector> ubvector_range;
730  typedef boost::numeric::ublas::slice slice;
731  typedef boost::numeric::ublas::range range;
732 
733 #ifndef DOXYGEN_INTERNAL
734 
735  protected:
736 
737  /// \name Storage for Akima spline interpolation
738  //@{
739  ubvector b;
740  ubvector c;
741  ubvector d;
742  ubvector um;
743  //@}
744 
745  /// For initializing the interpolation
746  void akima_calc(const vec_t &x_array, size_t size,
747  ubvector &umx) {
748 
749  for(size_t i=0;i<this->sz-1;i++) {
750 
751  double NE=fabs(umx[3+i]-umx[2+i])+fabs(umx[1+i]-umx[i]);
752 
753  if (NE == 0.0) {
754  b[i]=umx[2+i];
755  c[i]=0.0;
756  d[i]=0.0;
757  } else {
758  double h_i=(*this->px)[i+1]-(*this->px)[i];
759  double NE_next=fabs(umx[4+i]-umx[3+i])+
760  fabs(umx[2+i]-umx[1+i]);
761  double alpha_i=fabs(umx[1+i]-umx[i])/NE;
762  double alpha_ip1;
763  double tL_ip1;
764  if (NE_next == 0.0) {
765  tL_ip1=umx[2+i];
766  } else {
767  alpha_ip1=fabs(umx[2+i]-umx[1+i])/NE_next;
768  tL_ip1=(1.0-alpha_ip1)*umx[2+i]+alpha_ip1*umx[3+i];
769  }
770  b[i]=(1.0-alpha_i)*umx[1+i]+alpha_i*umx[2+i];
771  c[i]=(3.0*umx[2+i]-2.0*b[i]-tL_ip1)/h_i;
772  d[i]=(b[i]+tL_ip1-2.0*umx[2+i])/(h_i*h_i);
773  }
774  }
775  }
776 
777 #endif
778 
779  public:
780 
781  /** \brief Create a base interpolation object with or without
782  periodic boundary conditions
783  */
785  this->min_size=5;
786  }
787 
788  virtual ~interp_akima() {
789  }
790 
791  /** \brief Initialize interpolation routine
792  */
793  virtual void set(size_t size, const vec_t &xa, const vec2_t &ya) {
794 
795  if (size<this->min_size) {
796  O2SCL_ERR((((std::string)"Vector size, ")+szttos(size)+", is less"+
797  " than "+szttos(this->min_size)+" in interp_akima::"+
798  "set().").c_str(),exc_einval);
799  }
800 
801  if (size!=this->sz) {
802  b.resize(size);
803  c.resize(size);
804  d.resize(size);
805  um.resize(size+4);
806  }
807 
808  this->px=&xa;
809  this->py=&ya;
810  this->sz=size;
811 
812  this->svx.set_vec(size,xa);
813 
814  // Non-periodic boundary conditions
815 
816  ubvector_range m(um,range(2,um.size()));
817  size_t i;
818  for (i=0;i<=size-2;i++) {
819  m[i]=(ya[i+1]-ya[i])/(xa[i+1]-xa[i]);
820  }
821 
822  um[0]=3.0*m[0]-2.0*m[1];
823  um[1]=2.0*m[0]-m[1];
824  m[this->sz-1]=2.0*m[size-2]-m[size-3];
825  m[size]=3.0*m[size-2]-2.0*m[size-3];
826 
827  akima_calc(xa,size,um);
828 
829  return;
830  }
831 
832  /// Give the value of the function \f$ y(x=x_0) \f$ .
833  virtual double eval(double x0) const {
834 
835  size_t index=this->svx.find(x0);
836 
837  double x_lo=(*this->px)[index];
838  double delx=x0-x_lo;
839  double bb=b[index];
840  double cc=c[index];
841  double dd=d[index];
842 
843  return (*this->py)[index]+delx*(bb+delx*(cc+dd*delx));
844  }
845 
846  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
847  virtual double deriv(double x0) const {
848 
849  size_t index=this->svx.find(x0);
850 
851  double x_lo=(*this->px)[index];
852  double delx=x0-x_lo;
853  double bb=b[index];
854  double cc=c[index];
855  double dd=d[index];
856 
857  return bb+delx*(2.0*cc+3.0*dd*delx);
858  }
859 
860  /** \brief Give the value of the second derivative
861  \f$ y^{\prime \prime}(x=x_0) \f$ .
862  */
863  virtual double deriv2(double x0) const {
864 
865  size_t index=this->svx.find(x0);
866 
867  double x_lo=(*this->px)[index];
868  double delx=x0-x_lo;
869  double cc=c[index];
870  double dd=d[index];
871 
872  return 2.0*cc+6.0*dd*delx;
873  }
874 
875  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
876  virtual double integ(double aa, double bb) const {
877 
878  size_t i, index_a, index_b;
879 
880  bool flip=false;
881  if (((*this->px)[0]<(*this->px)[this->sz-1] && aa>bb) ||
882  ((*this->px)[0]>(*this->px)[this->sz-1] && aa<bb)) {
883  double tmp=aa;
884  aa=bb;
885  bb=tmp;
886  flip=true;
887  }
888 
889  index_a=this->svx.find(aa);
890  index_b=this->svx.find(bb);
891 
892  double result=0.0;
893 
894  for(i=index_a; i<=index_b; i++) {
895 
896  double x_lo=(*this->px)[i];
897  double x_hi=(*this->px)[i+1];
898  double y_lo=(*this->py)[i];
899  double dx=x_hi-x_lo;
900 
901  if (dx != 0.0) {
902 
903  if (i==index_a || i==index_b) {
904  double x1=(i==index_a) ? aa : x_lo;
905  double x2=(i==index_b) ? bb : x_hi;
906  result += this->integ_eval(y_lo,b[i],c[i],d[i],x_lo,x1,x2);
907  } else {
908  result+=dx*(y_lo+dx*(0.5*b[i]+dx*(c[i]/3.0+0.25*d[i]*dx)));
909  }
910 
911  } else {
912  double y_hi=(*this->py)[i+1];
913  std::string str=((std::string)"Interval of length zero ")+
914  "between ("+o2scl::dtos(x_lo)+","+o2scl::dtos(y_lo)+
915  ") and ("+o2scl::dtos(x_hi)+","+o2scl::dtos(y_hi)+
916  " in interp_akima::integ().";
917  O2SCL_ERR(str.c_str(),exc_einval);
918  }
919  }
920 
921  if (flip) result*=-1.0;
922 
923  return result;
924  }
925 
926  /// Return the type, \c "interp_akima".
927  virtual const char *type() const { return "interp_akima"; }
928 
929 #ifndef DOXYGEN_INTERNAL
930 
931  private:
932 
935 
936 #endif
937 
938  };
939 
940  /** \brief Akima spline interpolation with periodic
941  boundary conditions (GSL)
942 
943  See also the \ref intp_section section of the \o2 User's guide.
944  */
945  template<class vec_t, class vec2_t=vec_t> class interp_akima_peri :
946  public interp_akima<vec_t,vec2_t> {
947 
948  public:
949 
951  typedef boost::numeric::ublas::vector_slice<ubvector> ubvector_slice;
952  typedef boost::numeric::ublas::vector_range<ubvector> ubvector_range;
953  typedef boost::numeric::ublas::slice slice;
954  typedef boost::numeric::ublas::range range;
955 
956  public:
957 
959  }
960 
961  virtual ~interp_akima_peri() {
962  }
963 
964  /// Return the type, \c "interp_akima_peri".
965  virtual const char *type() const { return "interp_akima_peri"; }
966 
967  /// Initialize interpolation routine
968  virtual void set(size_t size, const vec_t &xa, const vec2_t &ya) {
969 
970  if (size<this->min_size) {
971  O2SCL_ERR((((std::string)"Vector size, ")+szttos(size)+", is less"+
972  " than "+szttos(this->min_size)+" in interp_akima"+
973  "_peri::set().").c_str(),exc_einval);
974  }
975 
976  if (size!=this->sz) {
977  this->b.resize(size);
978  this->c.resize(size);
979  this->d.resize(size);
980  this->um.resize(size+4);
981  }
982 
983  this->px=&xa;
984  this->py=&ya;
985  this->sz=size;
986 
987  this->svx.set_vec(size,xa);
988 
989  // Periodic boundary conditions
990 
991  ubvector_range m(this->um,range(2,this->um.size()));
992 
993  // Form the required set of divided differences
994  for (size_t i=0;i<=size-2;i++) {
995  m[i]=(ya[i+1]-ya[i])/(xa[i+1]-xa[i]);
996  }
997 
998  this->um[0]=m[this->sz-3];
999  this->um[1]=m[this->sz-2];
1000  m[this->sz-1]=m[0];
1001  m[this->sz]=m[1];
1002 
1003  this->akima_calc(xa,size,this->um);
1004 
1005  return;
1006  }
1007 
1008 #ifndef DOXYGEN_INTERNAL
1009 
1010  private:
1011 
1015 
1016 #endif
1017 
1018  };
1019 
1020  /** \brief Monotonicity-preserving interpolation
1021 
1022  \warning This class is experimental. Integrals don't work yet.
1023 
1024  This class uses a method based on cubic Hermite interpolation,
1025  modifying the slopes to guarantee monotonicity. In the
1026  notation of \ref Fritsch80, if
1027  \f[
1028  \alpha_i^2+\beta_i^2 \geq 9
1029  \f]
1030  then \f$ \alpha \f$ and \f$ \beta \f$ are decreased by
1031  the factor \f$ \tau \f$ as described at the end of
1032  section 4.
1033 
1034  \note The results of the interpolation will only be monotonic in
1035  the regions where the original data set is also monotonic. Also,
1036  this interpolation routine does not enforce strict monotonicity,
1037  and the results of the interpolation will be flat where the data
1038  is also flat.
1039 
1040  Based on \ref Fritsch80 .
1041 
1042  \future Convert into fewer loops over the data
1043  */
1044  template<class vec_t, class vec2_t=vec_t> class interp_monotonic :
1045  public interp_base<vec_t,vec2_t> {
1046 
1047 #ifdef O2SCL_NEVER_DEFINED
1048  }{
1049 #endif
1050 
1051  public:
1052 
1053  typedef boost::numeric::ublas::vector<double> ubvector;
1054 
1055  protected:
1056 
1057  /// Slopes
1058  ubvector m;
1059  /// Finite differences
1060  ubvector Delta;
1061  /// Ratio
1062  ubvector alpha;
1063  /// Staggered ratio
1064  ubvector beta;
1065 
1066  public:
1067 
1068  interp_monotonic() {
1069  this->min_size=2;
1070  }
1071 
1072  virtual ~interp_monotonic() {
1073  }
1074 
1075  /// Initialize interpolation routine
1076  virtual void set(size_t size, const vec_t &x, const vec2_t &y) {
1077 
1078  // Verify size
1079  if (size<this->min_size) {
1080  O2SCL_ERR((((std::string)"Vector size, ")+szttos(size)+", is less"+
1081  " than "+szttos(this->min_size)+" in interp_monotonic::"+
1082  "set().").c_str(),exc_einval);
1083  }
1084 
1085  // Setup search_vec object
1086  this->svx.set_vec(size,x);
1087 
1088  // Resize internal vectors
1089  if (this->sz!=size) {
1090  m.resize(size);
1091  Delta.resize(size-1);
1092  alpha.resize(size-1);
1093  beta.resize(size-1);
1094  }
1095 
1096  // Copy pointers
1097  this->px=&x;
1098  this->py=&y;
1099  this->sz=size;
1100 
1101  // Create the interpolation arrays in this sub-interval
1102 
1103  // Compute Delta and m
1104  for(size_t i=0;i<size-1;i++) {
1105  Delta[i]=(y[i+1]-y[i])/(x[i+1]-x[i]);
1106  if (i>0) {
1107  m[i]=(Delta[i]+Delta[i-1])/2.0;
1108  }
1109  }
1110  m[0]=Delta[0];
1111  m[size-1]=Delta[size-2];
1112 
1113  // Check to see if the data is flat anywhere
1114  for(size_t i=0;i<size-1;i++) {
1115  if (y[i]==y[i+1]) {
1116  m[i]=0.0;
1117  m[i+1]=0.0;
1118  }
1119  }
1120 
1121  // Compute alpha and beta
1122  for(size_t i=0;i<size-1;i++) {
1123  alpha[i]=m[i]/Delta[i];
1124  beta[i]=m[i+1]/Delta[i];
1125  }
1126 
1127  // Constrain m to ensure monotonicity
1128  for(size_t i=0;i<size-1;i++) {
1129  double norm2=alpha[i]*alpha[i]+beta[i]*beta[i];
1130  if (norm2>9.0) {
1131  double tau=3.0/sqrt(norm2);
1132  m[i]=tau*alpha[i]*Delta[i];
1133  m[i+1]=tau*beta[i]*Delta[i];
1134  }
1135  }
1136 
1137  return;
1138  }
1139 
1140  /// Give the value of the function \f$ y(x=x_0) \f$ .
1141  virtual double eval(double x0) const {
1142 
1143  size_t index=this->svx.find(x0);
1144 
1145  double x_lo=(*this->px)[index];
1146  double x_hi=(*this->px)[index+1];
1147  double y_lo=(*this->py)[index];
1148  double y_hi=(*this->py)[index+1];
1149  double h=x_hi-x_lo;
1150  double t=(x0-x_lo)/h;
1151  double t2=t*t, t3=t2*t;
1152 
1153  double h00=2.0*t3-3.0*t2+1.0;
1154  double h10=t3-2.0*t2+t;
1155  double h01=-2.0*t3+3.0*t2;
1156  double h11=t3-t2;
1157  double interp=y_lo*h00+h*m[index]*h10+y_hi*h01+h*m[index+1]*h11;
1158 
1159  return interp;
1160  }
1161 
1162  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
1163  virtual double deriv(double x0) const {
1164 
1165  size_t index=this->svx.find(x0);
1166 
1167  double x_lo=(*this->px)[index];
1168  double x_hi=(*this->px)[index+1];
1169  double y_lo=(*this->py)[index];
1170  double y_hi=(*this->py)[index+1];
1171  double h=x_hi-x_lo;
1172  double t=(x0-x_lo)/h;
1173  double t2=t*t;
1174 
1175  double dh00=6.0*t2-6.0*t;
1176  double dh10=3.0*t2-4.0*t+1.0;
1177  double dh01=-6.0*t2+6.0*t;
1178  double dh11=3.0*t2-2.0*t;
1179  double deriv=(y_lo*dh00+h*m[index]*dh10+y_hi*dh01+
1180  h*m[index+1]*dh11)/h;
1181 
1182  return deriv;
1183  }
1184 
1185  /** \brief Give the value of the second derivative
1186  \f$ y^{\prime \prime}(x=x_0) \f$
1187  */
1188  virtual double deriv2(double x0) const {
1189 
1190  size_t index=this->svx.find(x0);
1191 
1192  double x_lo=(*this->px)[index];
1193  double x_hi=(*this->px)[index+1];
1194  double y_lo=(*this->py)[index];
1195  double y_hi=(*this->py)[index+1];
1196  double h=x_hi-x_lo;
1197  double t=(x0-x_lo)/h;
1198 
1199  double ddh00=12.0*t-6.0;
1200  double ddh10=6.0*t-4.0;
1201  double ddh01=-12.0*t+6.0;
1202  double ddh11=6.0*t-2.0;
1203  double deriv2=(y_lo*ddh00+h*m[index]*ddh10+y_hi*ddh01+
1204  h*m[index+1]*ddh11)/h/h;
1205 
1206  return deriv2;
1207  }
1208 
1209  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
1210  virtual double integ(double a, double b) const {
1211 
1212  size_t i, index_a, index_b;
1213 
1214  bool flip=false;
1215  if (((*this->px)[0]<(*this->px)[this->sz-1] && a>b) ||
1216  ((*this->px)[0]>(*this->px)[this->sz-1] && a<b)) {
1217  double tmp=a;
1218  a=b;
1219  b=tmp;
1220  flip=true;
1221  }
1222 
1223  index_a=this->svx.find(a);
1224  index_b=this->svx.find(b);
1225 
1226  double result=0.0;
1227 
1228  for(i=index_a; i<=index_b; i++) {
1229 
1230  double x_hi=(*this->px)[i+1];
1231  double x_lo=(*this->px)[i];
1232  double y_lo=(*this->py)[i];
1233  double y_hi=(*this->py)[i+1];
1234  double h=x_hi-x_lo;
1235 
1236  if (h != 0.0) {
1237 
1238  if (i == index_a) {
1239  x_lo=a;
1240  }
1241  if (i == index_b) {
1242  x_hi=b;
1243  }
1244 
1245  double t=(x_hi-x_lo)/h;
1246  double t2=t*t, t3=t2*t, t4=t3*t;
1247 
1248  double ih00=t4/2.0-t3+t;
1249  double ih10=t4/4.0-2.0*t3/3.0+t2/2.0;
1250  double ih01=-t4/2.0+t3;
1251  double ih11=t4/4.0-t3/3.0;
1252  double intres=h*(y_lo*ih00+h*m[i]*ih10+y_hi*ih01+
1253  h*m[i+1]*ih11);
1254  result+=intres;
1255 
1256  } else {
1257  std::string str=((std::string)"Interval of length zero ")+
1258  "between ("+o2scl::dtos(x_lo)+","+o2scl::dtos(y_lo)+
1259  ") and ("+o2scl::dtos(x_hi)+","+o2scl::dtos(y_hi)+
1260  " in interp_monotonic::integ().";
1261  O2SCL_ERR(str.c_str(),exc_einval);
1262  }
1263 
1264  }
1265 
1266  if (flip) result*=-1.0;
1267 
1268  return result;
1269  }
1270 
1271  /// Return the type, \c "interp_monotonic".
1272  virtual const char *type() const { return "interp_monotonic"; }
1273 
1274 #ifndef DOXYGEN_INTERNAL
1275 
1276  private:
1277 
1280  (const interp_monotonic<vec_t,vec2_t>&);
1281 
1282 #endif
1283 
1284  };
1285 
1286  /** \brief Interpolation class for general vectors
1287 
1288  See also the \ref intp_section section of the \o2 User's guide.
1289 
1290  Interpolation of ublas vector like objects is performed with the
1291  default template parameters, and \ref interp_array is provided
1292  for simple interpolation on C-style arrays.
1293 
1294  The type of interpolation to be performed can be specified using
1295  the set_type() function or in the constructor. The default is
1296  cubic splines with natural boundary conditions.
1297  */
1298  template<class vec_t=boost::numeric::ublas::vector<double>,
1299  class vec2_t=vec_t> class interp {
1300 
1301 #ifndef DOXYGEN_INTERNAL
1302 
1303  protected:
1304 
1305  /// Pointer to base interpolation object
1307 
1308 #endif
1309 
1310  public:
1311 
1312  /// Create with base interpolation object \c it
1313  interp(size_t interp_type=itp_cspline) {
1314  if (interp_type==itp_linear) {
1316  } else if (interp_type==itp_cspline) {
1318  } else if (interp_type==itp_cspline_peri) {
1320  } else if (interp_type==itp_akima) {
1322  } else if (interp_type==itp_akima_peri) {
1324  } else if (interp_type==itp_monotonic) {
1326  } else {
1327  O2SCL_ERR((((std::string)"Invalid interpolation type, ")+
1328  o2scl::szttos(interp_type)+", in "+
1329  "interp::interp().").c_str(),exc_einval);
1330  }
1331  }
1332 
1333  virtual ~interp() {
1334  delete itp;
1335  }
1336 
1337  /// Give the value of the function \f$ y(x=x_0) \f$ .
1338  virtual double eval(const double x0, size_t n, const vec_t &x,
1339  const vec2_t &y) {
1340  itp->set(n,x,y);
1341  return itp->eval(x0);
1342  }
1343 
1344  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
1345  virtual double deriv(const double x0, size_t n, const vec_t &x,
1346  const vec2_t &y) {
1347  itp->set(n,x,y);
1348  return itp->deriv(x0);
1349  }
1350 
1351  /** \brief Give the value of the second derivative
1352  \f$ y^{\prime \prime}(x=x_0) \f$ .
1353  */
1354  virtual double deriv2(const double x0, size_t n, const vec_t &x,
1355  const vec2_t &y) {
1356  itp->set(n,x,y);
1357  return itp->deriv2(x0);
1358  }
1359 
1360  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
1361  virtual double integ(const double x1, const double x2, size_t n,
1362  const vec_t &x, const vec2_t &y) {
1363  itp->set(n,x,y);
1364  return itp->integ(x1,x2);
1365  }
1366 
1367  /// Set base interpolation type
1368  void set_type(size_t interp_type) {
1369  delete itp;
1370  if (interp_type==itp_linear) {
1372  } else if (interp_type==itp_cspline) {
1374  } else if (interp_type==itp_cspline_peri) {
1376  } else if (interp_type==itp_akima) {
1378  } else if (interp_type==itp_akima_peri) {
1380  } else {
1381  O2SCL_ERR((((std::string)"Invalid interpolation type, ")+
1382  o2scl::szttos(interp_type)+", in "+
1383  "interp::set().").c_str(),exc_einval);
1384  }
1385  return;
1386  }
1387 
1388 #ifndef DOXYGEN_INTERNAL
1389 
1390  private:
1391 
1393  interp<vec_t,vec2_t>& operator=(const interp<vec_t,vec2_t>&);
1394 
1395 #endif
1396 
1397  };
1398 
1399  /** \brief Interpolation class for pre-specified vector
1400 
1401  See also the \ref intp_section section of the \o2 User's guide.
1402 
1403  This interpolation class is intended to be sufficiently general
1404  to handle any vector type. Interpolation of ublas vector like
1405  objects is performed with the default template parameters.
1406 
1407  This class does not double check the vector to ensure that
1408  all of the intervals for the abcissa are all increasing or
1409  all decreasing.
1410 
1411  The type of interpolation to be performed can be specified
1412  using the set_type() function. The default is cubic splines
1413  with natural boundary conditions.
1414 
1415  \future Make version which copies vectors
1416  rather than storing pointers might be better and then
1417  has copy constructors.
1418  */
1419  template<class vec_t=boost::numeric::ublas::vector<double>,
1420  class vec2_t=vec_t> class interp_vec :
1421  public interp_base<vec_t,vec2_t> {
1422 
1423 #ifndef DOXYGEN_INTERNAL
1424 
1425  protected:
1426 
1427  /// Base interpolation object
1429 
1430  /// Interpolation type
1431  size_t itype;
1432 
1433 #endif
1434 
1435  public:
1436 
1437  /// Blank interpolator
1439  itp=0;
1441  }
1442 
1443  /// Create with base interpolation object \c it
1444  interp_vec(size_t n, const vec_t &x,
1445  const vec2_t &y, size_t interp_type=itp_cspline) {
1446 
1447  if (x[0]==x[n-1]) {
1448  O2SCL_ERR((((std::string)"Vector endpoints equal (value=")+
1449  o2scl::dtos(x[0])+") in interp_vec()::"+
1450  "interp_vec().").c_str(),exc_einval);
1451  }
1452 
1453  if (interp_type==itp_linear) {
1455  } else if (interp_type==itp_cspline) {
1457  } else if (interp_type==itp_cspline_peri) {
1459  } else if (interp_type==itp_akima) {
1461  } else if (interp_type==itp_akima_peri) {
1463  } else if (interp_type==itp_monotonic) {
1465  } else {
1466  O2SCL_ERR((((std::string)"Invalid interpolation type, ")+
1467  o2scl::szttos(interp_type)+", in "+
1468  "interp_vec::interp_vec().").c_str(),exc_einval);
1469  }
1470  itype=interp_type;
1471 
1472  itp->set(n,x,y);
1473  }
1474 
1475  virtual ~interp_vec() {
1476  if (itp!=0) {
1477  delete itp;
1478  }
1479  }
1480 
1481  /** \brief Set a new vector to interpolate
1482  */
1483  void set(size_t n, const vec_t &x, const vec2_t &y) {
1484  set(n,x,y,itype);
1485  return;
1486  }
1487 
1488  /** \brief Set a new vector to interpolate
1489  */
1490  void set(size_t n, const vec_t &x,
1491  const vec2_t &y, size_t interp_type) {
1492 
1493  if (x[0]==x[n-1]) {
1494  O2SCL_ERR((((std::string)"Vector endpoints equal (value=")+
1495  o2scl::dtos(x[0])+") in interp_vec()::"+
1496  "interp_vec().").c_str(),exc_einval);
1497  }
1498 
1499  delete itp;
1500  if (interp_type==itp_linear) {
1502  } else if (interp_type==itp_cspline) {
1504  } else if (interp_type==itp_cspline_peri) {
1506  } else if (interp_type==itp_akima) {
1508  } else if (interp_type==itp_akima_peri) {
1510  } else if (interp_type==itp_monotonic) {
1512  } else {
1513  O2SCL_ERR((((std::string)"Invalid interpolation type, ")+
1514  o2scl::szttos(interp_type)+", in "+
1515  "interp_vec::set().").c_str(),exc_einval);
1516  }
1517  itype=interp_type;
1518 
1519  itp->set(n,x,y);
1520  }
1521 
1522  /** \brief Manually clear the pointer to the user-specified vector
1523  */
1524  void clear() {
1525  if (itp!=0) {
1526  delete itp;
1527  itp=0;
1528  }
1529  return;
1530  }
1531 
1532  /// Give the value of the function \f$ y(x=x_0) \f$ .
1533  virtual double eval(const double x0) const {
1534  if (itp==0) {
1535  O2SCL_ERR("No vector set in interp_vec::eval().",
1536  exc_einval);
1537  }
1538  return itp->eval(x0);
1539  }
1540 
1541  /// Give the value of the function \f$ y(x=x_0) \f$ .
1542  virtual double operator()(double x0) const {
1543  if (itp==0) {
1544  O2SCL_ERR("No vector set in interp_vec::operator().",
1545  exc_einval);
1546  }
1547  return itp->eval(x0);
1548  }
1549 
1550  /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
1551  virtual double deriv(const double x0) const {
1552  if (itp==0) {
1553  O2SCL_ERR("No vector set in interp_vec::deriv().",
1554  exc_einval);
1555  }
1556  return itp->deriv(x0);
1557  }
1558 
1559  /** \brief Give the value of the second derivative
1560  \f$ y^{\prime \prime}(x=x_0) \f$ .
1561  */
1562  virtual double deriv2(const double x0) const {
1563  if (itp==0) {
1564  O2SCL_ERR("No vector set in interp_vec::deriv2().",
1565  exc_einval);
1566  }
1567  return itp->deriv2(x0);
1568  }
1569 
1570  /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
1571  virtual double integ(const double x1, const double x2) const {
1572  if (itp==0) {
1573  O2SCL_ERR("No vector set in interp_vec::integ().",
1574  exc_einval);
1575  }
1576  return itp->integ(x1,x2);
1577  }
1578 
1579  /// Return the type, "interp_vec"
1580  virtual const char *type() const {
1581  return "interp_vec";
1582  }
1583 
1584 #ifndef DOXYGEN_INTERNAL
1585 
1586  private:
1587 
1589  interp_vec<vec_t,vec2_t> &operator=
1590  (const interp_vec<vec_t,vec2_t> &it);
1591 
1592 #endif
1593 
1594  };
1595 
1596  /** \brief A specialization of interp for C-style double arrays
1597 
1598  See also the \ref intp_section section of the \o2 User's guide.
1599  */
1600  template<size_t n> class interp_array :
1601  public interp<double[n]> {
1602 
1603  public:
1604 
1605  /// Create with base interpolation objects \c it and \c rit
1606  interp_array(size_t interp_type)
1607  : interp<double[n]>(interp_type) {}
1608 
1609  /** \brief Create an interpolator using the default base
1610  interpolation objects
1611  */
1612  interp_array() : interp<double[n]>() {}
1613 
1614  };
1615 
1616  /** \brief A specialization of \ref o2scl::interp_vec for C-style arrays
1617 
1618  See also the \ref intp_section section of the \o2 User's guide.
1619  */
1620  template<class arr_t> class interp_array_vec :
1621  public interp_vec<arr_t> {
1622 
1623  public:
1624 
1625  /// Create with base interpolation object \c it
1626  interp_array_vec(size_t nv, const arr_t &x, const arr_t &y,
1627  size_t interp_type) :
1628  interp_vec<arr_t>(nv,x,y,interp_type) {}
1629  };
1630 
1631  /** \brief Count level crossings
1632 
1633  This function counts the number of times the function \f$ y(x) =
1634  \mathrm{level} \f$ where the function is defined by the vectors
1635  \c x and \c y.
1636 
1637  The value returned is exactly the same as the size of the
1638  \c locs vector computed by \ref vector_find_level().
1639 
1640  This function will call the error handler if \c n is
1641  less than two.
1642 
1643  If one of the entries in \c y is either larger or smaller than
1644  it's neighbors (i.e. if the function \f$ y(x) \f$ has an
1645  extremum), and if the value of \c level is exactly equal to the
1646  extremum, then this is counted as 1 level crossing. The same
1647  applies if either the first or last entry in \c y is exactly
1648  equal to \c level .
1649  */
1650  template<class vec_t, class vec2_t> size_t vector_level_count
1651  (double level, size_t n, vec_t &x, vec2_t &y) {
1652 
1653  if (n<=1) {
1654  O2SCL_ERR2("Need at least two data points in ",
1655  "vector_find_count().",exc_einval);
1656  }
1657 
1658  size_t count=0;
1659 
1660  // Handle left boundary
1661  if (y[0]==level) count++;
1662 
1663  // Find intersections
1664  for(size_t i=0;i<n-1;i++) {
1665 
1666  if ((y[i]<level && y[i+1]>=level) ||
1667  (y[i]>level && y[i+1]<=level)) {
1668  count++;
1669  }
1670  }
1671 
1672  return count;
1673  }
1674 
1675  /** \brief Perform inverse linear interpolation
1676 
1677  This function performs inverse linear interpolation of the data
1678  defined by \c x and \c y, finding all points in \c x which have
1679  the property \f$ \mathrm{level} = y(x) \f$. All points for which
1680  this relation holds are put into the vector \c locs. The
1681  previous information contained in vector \c locs before the
1682  function call is destroyed.
1683 
1684  This is the 1-dimensional analog of finding contour levels as
1685  the \ref contour class does for 2 dimensions.
1686 
1687  This function will call the error handler if \c n is
1688  less than two.
1689 
1690  This function is inclusive of the endpoints, so that if either
1691  <tt>y[0]</tt> and/or <tt>y[n-1]</tt> is equal to level, then
1692  <tt>x[0]</tt> and/or <tt>x[n-1]</tt> are added to \c locs,
1693  respectively.
1694  */
1695  template<class vec_t, class vec2_t> void vector_find_level
1696  (double level, size_t n, vec_t &x, vec2_t &y, std::vector<double> &locs) {
1697 
1698  if (n<=1) {
1699  O2SCL_ERR2("Need at least two data points in ",
1700  "vector_find_level().",exc_einval);
1701  }
1702 
1703  // Ensure that the location vector is empty
1704  locs.clear();
1705 
1706  // Handle left boundary
1707  if (y[0]==level) {
1708  locs.push_back(x[0]);
1709  }
1710 
1711  // Find intersections
1712  for(size_t i=0;i<n-1;i++) {
1713 
1714  if ((y[i]<level && y[i+1]>level) ||
1715  (y[i]>level && y[i+1]<level)) {
1716  // For each intersection, add the location using linear
1717  // interpolation
1718  double x0=x[i]+(x[i+1]-x[i])*(level-y[i])/(y[i+1]-y[i]);
1719  locs.push_back(x0);
1720  } else if (y[i+1]==level) {
1721  locs.push_back(x[i+1]);
1722  }
1723  }
1724 
1725  return;
1726  }
1727 
1728  /** \brief Compute the integral over <tt>y(x)</tt> using linear
1729  interpolation
1730 
1731  Assuming vectors \c y and \c x define a function \f$ y(x) \f$
1732  then this computes
1733  \f[
1734  \int_{x_0}^{x^{n-1}} y(x)~dx
1735  \f]
1736  using the trapezoidal rule implied by linear interpolation.
1737 
1738  \future It might be faster to compute the sum directly rather
1739  than going through an \ref o2scl::interp object .
1740  */
1741  template<class vec_t, class vec2_t>
1742  double vector_integ_linear(size_t n, vec_t &x, vec2_t &y) {
1743 
1744  // Interpolation object
1746 
1747  // Compute full integral
1748  double total=si.integ(x[0],x[n-1],n,x,y);
1749 
1750  return total;
1751  }
1752 
1753  /** \brief Compute the endpoints which enclose the regions whose
1754  integral is equal to \c sum
1755 
1756  Defining a new function, \f$ g(y_0) \f$ which takes as input any
1757  y-value, \f$ y_0 \f$ from the function \f$ y(x) \f$ (specified
1758  with the parameters \c x and \c y) and outputs the integral of
1759  the function \f$ y(x) \f$ over all regions where \f$ y(x) > y_0
1760  \f$. This function inverts \f$ g(y) \f$, taking the value
1761  of an integral as input, and returns the corresponding y-value
1762  in the variable \c lev.
1763 
1764  In order to make sure that the intepretation of the integral is
1765  unambiguous, this function requires that the first and last values
1766  of \c y are equal, i.e. <tt>y[0]==y[n-1]</tt>. This function
1767  also requires at least two data points, i.e. \c n cannot be
1768  0 or 1.
1769 
1770  This function is particularly useful, for example, in computing
1771  the region which defines 68\% around a peak of data, thus
1772  providing \f$ 1~\sigma \f$ confidence limits.
1773 
1774  Linear interpolation is used to describe the function \f$ g \f$,
1775  and the precision of this function is limited by this assumption.
1776  This function may also sometimes fail if \c sum is very close to
1777  the minimum or maximum value of the function \f$ g \f$.
1778 
1779  \todo This function may also require that all of the
1780  y-vector values have the same sign or are all positive.
1781  Check this.
1782 
1783  \comment
1784  Note that the two vector types for x and y must be the
1785  same in order to use o2scl_interp.
1786  \endcomment
1787  */
1788  template<class vec_t, class vec2_t> void vector_invert_enclosed_sum
1789  (double sum, size_t n, vec_t &x, vec2_t &y, double &lev,
1790  int verbose=0) {
1791 
1792  if (n<=1) {
1793  O2SCL_ERR2("Need at least two data points in ",
1794  "vector_invert_enclosed_sum().",exc_einval);
1795  }
1796 
1797  if (y[0]!=y[n-1]) {
1798  O2SCL_ERR2("The first and last y-values must be equal in ",
1799  "vector_invert_enclosed_sum().",exc_einval);
1800  }
1801 
1802  // Construct a sorted list of function values
1803  typedef boost::numeric::ublas::vector<double> ubvector;
1804  ubvector ysort(n);
1805  vector_copy(n,y,ysort);
1806  vector_sort<ubvector,double>(ysort.size(),ysort);
1807 
1808  // Create list of y-values to perform y-value and integral
1809  // interpolation. We choose values in between the grid points to
1810  // ensure that we don't accidentally land precisely on a peak or
1811  // valley.
1812  std::vector<double> ylist;
1813  for(size_t i=0;i<ysort.size()-1;i++) {
1814  if (ysort[i]!=ysort[i+1]) {
1815  ylist.push_back((ysort[i+1]+3.0*ysort[i])/4.0);
1816  ylist.push_back((ysort[i+1]*3.0+ysort[i])/4.0);
1817  }
1818  }
1819 
1820  // Interpolation objects. We need two, one for the
1821  // user-specified vector type, and one for std::vector<double>
1823  interp<std::vector<double>,std::vector<double> > itp2(itp_linear);
1824 
1825  // Construct vectors for interpolation
1826  std::vector<double> xi, yi;
1827 
1828  size_t nfail=0;
1829 
1830  for(size_t k=0;k<ylist.size();k++) {
1831  double lev_tmp=ylist[k];
1832  std::vector<double> locs;
1833  vector_find_level(lev_tmp,n,x,y,locs);
1834  if ((locs.size()%2)!=0) {
1835  nfail++;
1836  } else {
1837  double sum_temp=0.0;
1838  for(size_t i=0;i<locs.size()/2;i++) {
1839  double x0=locs[2*i];
1840  double x1=locs[2*i+1];
1841  sum_temp+=itp.integ(x0,x1,n,x,y);
1842  }
1843  xi.push_back(sum_temp);
1844  yi.push_back(lev_tmp);
1845  }
1846  }
1847  if (nfail>10) {
1848  O2SCL_ERR2("More than 10 failures in ",
1849  "vector_invert_enclosed_sum().",exc_einval);
1850  }
1851 
1852  if (verbose>1) {
1853  std::cout << "i, xi, yi: " << std::endl;
1854  for(size_t i=0;i<xi.size();i++) {
1855  std::cout << i << " " << xi[i] << " " << yi[i] << std::endl;
1856  }
1857  }
1858 
1859  lev=itp2.eval(sum,xi.size(),xi,yi);
1860 
1861  return;
1862  }
1863 
1864  /** \brief Find the region enclosing a partial integral
1865  */
1866  template<class vec_t, class vec2_t> void vector_region_parint
1867  (size_t n, vec_t &x, vec2_t &y, double frac, std::vector<double> &locs,
1868  int verbose=0) {
1869 
1870  if (frac<0.0 || frac>1.0) {
1871  O2SCL_ERR2("Fraction must be between 0 and 1 (exclusive) in ",
1872  "vector_region_parint().",exc_efailed);
1873  }
1874 
1875  // Total integral
1876  double total=vector_integ_linear(n,x,y);
1877  if (verbose>0) {
1878  std::cout << "Total integral: " << total << std::endl;
1879  }
1880  // Specified partial integral
1881  double partial=frac*total;
1882  if (verbose>0) {
1883  std::cout << "Partial integral: " << partial << std::endl;
1884  }
1885  // Find correct level
1886  double lev;
1887  vector_invert_enclosed_sum(partial,n,x,y,lev,verbose);
1888  if (verbose>0) {
1889  std::cout << "Level from vector_invert: " << lev << std::endl;
1890  }
1891  // Inverse interpolate to find locations corresponding to level
1892  vector_find_level(lev,n,x,y,locs);
1893  if (verbose>0) {
1894  std::cout << "Locations from vector_find_level: ";
1895  for(size_t i=0;i<locs.size();i++) {
1896  std::cout << locs[i];
1897  if (i!=locs.size()-1) {
1898  std::cout << " ";
1899  }
1900  }
1901  std::cout << std::endl;
1902  }
1903  return;
1904  }
1905 
1906  /** \brief Find the boundaries of the region enclosing a partial integral
1907  */
1908  template<class vec_t, class vec2_t> void vector_bound_parint
1909  (size_t n, vec_t &x, vec2_t &y, double frac, double &low, double &high) {
1910 
1911  std::vector<double> locs;
1912  vector_region_parint(n,x,y,frac,locs);
1913  if (locs.size()==0) {
1914  O2SCL_ERR("Zero level crossings in vector_bound_sigma().",
1915  exc_efailed);
1916  }
1917  // Return min and max location
1918  size_t ix;
1919  vector_min(locs.size(),locs,ix,low);
1920  vector_max(locs.size(),locs,ix,high);
1921  return;
1922  }
1923 
1924 #ifndef DOXYGEN_NO_O2NS
1925 }
1926 #endif
1927 
1928 #endif
void vector_region_parint(size_t n, vec_t &x, vec2_t &y, double frac, std::vector< double > &locs, int verbose=0)
Find the region enclosing a partial integral.
Definition: interp.h:1867
Interpolation class for pre-specified vector.
Definition: interp.h:1420
virtual double operator()(double x0) const
Give the value of the function .
Definition: interp.h:148
o2scl_linalg::ubvector_5_mem p5m
Memory for the tridiagonalization.
Definition: interp.h:586
virtual double integ(double a, double b) const
Give the value of the integral .
Definition: interp.h:250
virtual void set(size_t size, const vec_t &xa, const vec2_t &ya)
Initialize interpolation routine.
Definition: interp.h:381
virtual double eval(const double x0, size_t n, const vec_t &x, const vec2_t &y)
Give the value of the function .
Definition: interp.h:1338
virtual const char * type() const
Return the type, "interp_akima_peri".
Definition: interp.h:965
interp_vec(size_t n, const vec_t &x, const vec2_t &y, size_t interp_type=itp_cspline)
Create with base interpolation object it.
Definition: interp.h:1444
virtual double integ(double aa, double bb) const
Give the value of the integral .
Definition: interp.h:876
double integ_eval(double ai, double bi, double ci, double di, double xi, double a, double b) const
An internal function to assist in computing the integral for both the cspline and Akima types...
Definition: interp.h:110
interp(size_t interp_type=itp_cspline)
Create with base interpolation object it.
Definition: interp.h:1313
virtual void set(size_t size, const vec_t &xa, const vec2_t &ya)
Initialize interpolation routine.
Definition: interp.h:608
void solve_tridiag_sym(const vec_t &diag, const vec2_t &offdiag, const vec3_t &b, vec4_t &x, size_t N, mem_t &m)
Solve a symmetric tridiagonal linear system with user-specified memory.
Definition: tridiag_base.h:118
void akima_calc(const vec_t &x_array, size_t size, ubvector &umx)
For initializing the interpolation.
Definition: interp.h:746
size_t min_size
The minimum size of the vectors to interpolate between.
Definition: interp.h:139
A specialization of interp for C-style double arrays.
Definition: interp.h:1600
virtual double integ(const double x1, const double x2) const
Give the value of the integral .
Definition: interp.h:1571
interp_base< vec_t, vec2_t > * itp
Pointer to base interpolation object.
Definition: interp.h:1306
Akima spline interpolation with periodic boundary conditions (GSL)
Definition: interp.h:945
void clear()
Manually clear the pointer to the user-specified vector.
Definition: interp.h:1524
virtual double deriv(double x0) const
Give the value of the derivative .
Definition: interp.h:228
virtual const char * type() const
Return the type, "interp_cspline".
Definition: interp.h:557
size_t sz
Vector size.
Definition: interp.h:105
virtual double deriv2(double x0) const
Give the value of the second derivative (always zero)
Definition: interp.h:245
virtual double deriv2(const double x0) const
Give the value of the second derivative .
Definition: interp.h:1562
virtual double deriv(double x0) const
Give the value of the derivative .
Definition: interp.h:1163
virtual double eval(double x0) const
Give the value of the function .
Definition: interp.h:441
invalid argument supplied by user
Definition: err_hnd.h:59
const vec_t * px
Independent vector.
Definition: interp.h:99
virtual double deriv2(double x0) const
Give the value of the second derivative .
Definition: interp.h:863
Monotonicity-preserving interpolation.
Definition: interp.h:1044
Base low-level interpolation class [abstract base].
Definition: interp.h:80
interp_array_vec(size_t nv, const arr_t &x, const arr_t &y, size_t interp_type)
Create with base interpolation object it.
Definition: interp.h:1626
virtual const char * type() const
Return the type, "interp_cspline_peri".
Definition: interp.h:604
Akima spline for periodic boundary conditions.
Definition: interp.h:56
virtual double integ(double a, double b) const
Give the value of the integral .
Definition: interp.h:1210
void vector_min(size_t n, const vec_t &data, size_t &index, data_t &val)
Compute the minimum of the first n elements of a vector.
Definition: vector.h:740
interp_vec()
Blank interpolator.
Definition: interp.h:1438
virtual double deriv(double x0) const
Give the value of the derivative .
Definition: interp.h:461
generic failure
Definition: err_hnd.h:61
interp_akima()
Create a base interpolation object with or without periodic boundary conditions.
Definition: interp.h:784
virtual void set(size_t size, const vec_t &xa, const vec2_t &ya)
Initialize interpolation routine.
Definition: interp.h:968
Cubic spline interpolation with periodic boundary conditions (GSL)
Definition: interp.h:576
virtual double integ(const double x1, const double x2, size_t n, const vec_t &x, const vec2_t &y)
Give the value of the integral .
Definition: interp.h:1361
virtual double eval(double x0) const
Give the value of the function .
Definition: interp.h:214
Cubic spline for natural boundary conditions.
Definition: interp.h:50
virtual double integ(double a, double b) const
Give the value of the integral .
Definition: interp.h:503
void vector_bound_parint(size_t n, vec_t &x, vec2_t &y, double frac, double &low, double &high)
Find the boundaries of the region enclosing a partial integral.
Definition: interp.h:1909
virtual double deriv2(double x0) const
Give the value of the second derivative .
Definition: interp.h:1188
ubvector beta
Staggered ratio.
Definition: interp.h:1064
Allocation object for 5 arrays of equal size.
Definition: tridiag_base.h:86
virtual double eval(double x0) const
Give the value of the function .
Definition: interp.h:833
void set(size_t n, const vec_t &x, const vec2_t &y, size_t interp_type)
Set a new vector to interpolate.
Definition: interp.h:1490
interp_array()
Create an interpolator using the default base interpolation objects.
Definition: interp.h:1612
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
Definition: err_hnd.h:281
size_t itype
Interpolation type.
Definition: interp.h:1431
std::string dtos(double x, int prec=6, bool auto_prec=false)
Convert a double to a string.
virtual const char * type() const
Return the type, "interp_monotonic".
Definition: interp.h:1272
const vec2_t * py
Dependent vector.
Definition: interp.h:102
o2scl_linalg::ubvector_4_mem p4m
Memory for the tridiagonalization.
Definition: interp.h:350
virtual double deriv(double x0) const
Give the value of the derivative .
Definition: interp.h:847
size_t find(const double x0) const
Search an increasing or decreasing vector for the interval containing x0
Definition: search_vec.h:136
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
Definition: err_hnd.h:273
Allocation object for 4 arrays of equal size.
Definition: tridiag_base.h:70
void set_type(size_t interp_type)
Set base interpolation type.
Definition: interp.h:1368
size_t vector_level_count(double level, size_t n, vec_t &x, vec2_t &y)
Count level crossings.
Definition: interp.h:1651
Akima spline for natural boundary conditions.
Definition: interp.h:54
double vector_integ_linear(size_t n, vec_t &x, vec2_t &y)
Compute the integral over y(x) using linear interpolation.
Definition: interp.h:1742
Akima spline interpolation (GSL)
Definition: interp.h:722
void coeff_calc(const ubvector &c_array, double dy, double dx, size_t index, double &b, double &c2, double &d) const
Compute coefficients for cubic spline interpolation.
Definition: interp.h:353
void vector_invert_enclosed_sum(double sum, size_t n, vec_t &x, vec2_t &y, double &lev, int verbose=0)
Compute the endpoints which enclose the regions whose integral is equal to sum.
Definition: interp.h:1789
void vector_max(size_t n, const vec_t &data, size_t &index, data_t &val)
Compute the maximum of the first n elements of a vector.
Definition: vector.h:683
interp_cspline()
Create a base interpolation object with natural or periodic boundary conditions.
Definition: interp.h:372
Cubic spline for periodic boundary conditions.
Definition: interp.h:52
virtual double deriv(const double x0) const
Give the value of the derivative .
Definition: interp.h:1551
Monotonicity-preserving interpolation.
Definition: interp.h:58
virtual double deriv2(double x0) const
Give the value of the second derivative .
Definition: interp.h:483
virtual double deriv2(const double x0, size_t n, const vec_t &x, const vec2_t &y)
Give the value of the second derivative .
Definition: interp.h:1354
ubvector Delta
Finite differences.
Definition: interp.h:1060
Cubic spline interpolation (GSL)
Definition: interp.h:322
virtual double eval(double x0) const
Give the value of the function .
Definition: interp.h:1141
interp_array(size_t interp_type)
Create with base interpolation objects it and rit.
Definition: interp.h:1606
virtual const char * type() const
Return the type, "interp_akima".
Definition: interp.h:927
void vector_copy(vec_t &src, vec2_t &dest)
Simple generic vector copy.
Definition: vector.h:73
void solve_cyc_tridiag_sym(const vec_t &diag, const vec2_t &offdiag, const vec3_t &b, vec4_t &x, size_t N, mem_t &m)
Solve a symmetric cyclic tridiagonal linear system with user specified memory.
Definition: tridiag_base.h:241
virtual void set(size_t size, const vec_t &x, const vec2_t &y)
Initialize interpolation routine.
Definition: interp.h:200
virtual double eval(const double x0) const
Give the value of the function .
Definition: interp.h:1533
virtual const char * type() const
Return the type, "interp_linear".
Definition: interp.h:300
virtual double operator()(double x0) const
Give the value of the function .
Definition: interp.h:1542
void set_vec(size_t nn, const vec_t &x)
Set the vector to be searched.
Definition: search_vec.h:119
search_vec< const vec_t > svx
To perform binary searches.
Definition: interp.h:96
ubvector alpha
Ratio.
Definition: interp.h:1062
static const double x2[5]
Definition: inte_qng_gsl.h:66
static const double x1[5]
Definition: inte_qng_gsl.h:48
Linear interpolation (GSL)
Definition: interp.h:184
std::string szttos(size_t x)
Convert a size_t to a string.
virtual const char * type() const
Return the type, "interp_vec".
Definition: interp.h:1580
void vector_find_level(double level, size_t n, vec_t &x, vec2_t &y, std::vector< double > &locs)
Perform inverse linear interpolation.
Definition: interp.h:1696
Linear.
Definition: interp.h:48
ubvector m
Slopes.
Definition: interp.h:1058
A specialization of o2scl::interp_vec for C-style arrays.
Definition: interp.h:1620
interp_base< vec_t, vec2_t > * itp
Base interpolation object.
Definition: interp.h:1428
virtual void set(size_t size, const vec_t &xa, const vec2_t &ya)
Initialize interpolation routine.
Definition: interp.h:793
Interpolation class for general vectors.
Definition: interp.h:1299
virtual double deriv(const double x0, size_t n, const vec_t &x, const vec2_t &y)
Give the value of the derivative .
Definition: interp.h:1345
void set(size_t n, const vec_t &x, const vec2_t &y)
Set a new vector to interpolate.
Definition: interp.h:1483
virtual void set(size_t size, const vec_t &x, const vec2_t &y)
Initialize interpolation routine.
Definition: interp.h:1076

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..