interp.h

00001 /*
00002   -------------------------------------------------------------------
00003   
00004   Copyright (C) 2006, 2007, Andrew W. Steiner
00005   
00006   This file is part of O2scl.
00007   
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012   
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017   
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 #ifndef O2SCL_INTERP_H
00024 #define O2SCL_INTERP_H
00025 
00026 #include <iostream>
00027 #include <string>
00028 
00029 // for gsl_linalg_solve_symm_tridiag for cspline_interp
00030 #include <gsl/gsl_linalg.h>
00031 
00032 #include <o2scl/err_hnd.h>
00033 #include <o2scl/collection.h>
00034 #include <o2scl/ovector_tlate.h>
00035 #include <o2scl/search_vec.h>
00036 
00037 #ifndef DOXYGENP
00038 namespace o2scl {
00039 #endif
00040   
00041   /** 
00042       \brief Base low-level interpolation class
00043 
00044       The descendants of this class are intended to be fast
00045       interpolation routines for increasing functions, leaving the
00046       some error handling, user-friendliness, and other more 
00047       complicated improvements for other classes.
00048 
00049       For any pair of vectors x and y into which you would like to
00050       interpolate, you need to call alloc() and init() first, and then 
00051       the interpolation functions, and then free(). If the next pair
00052       of vectors has the same size, then you need only to call init()
00053       before the next call to an interpolation function. If the vectors
00054       do not change, then you may call the interpolation functions in
00055       succession.
00056 
00057       All of the descendants are based on the GSL interpolation 
00058       routines and give identical results.
00059 
00060       \future These might work for decreasing functions by just
00061       replacing calls to search_vec::bsearch_inc() with
00062       search_vec::bsearch_dec(). If this is the case, then this should
00063       be rewritten accordingly.
00064   */
00065   template<class vec_t> class base_interp {
00066 
00067 #ifndef DOXYGEN_INTERNAL
00068 
00069   protected:
00070 
00071     /// The binary search object
00072     search_vec<vec_t> sv;
00073 
00074 #endif
00075 
00076   public:
00077 
00078     virtual ~base_interp() {}
00079 
00080     /** 
00081         \brief The minimum size of the vectors to interpolate between
00082 
00083         This needs to be set in the constructor of the children for
00084         access by the class user
00085     */
00086     size_t min_size;
00087 
00088     /// Allocate memory, assuming x and y have size \c size
00089     virtual int alloc(size_t size) { return 0; }
00090 
00091     /// Free allocated memory
00092     virtual int free() { return 0; }
00093 
00094     /// Initialize interpolation routine
00095     virtual int init(const vec_t &x, const vec_t &y, size_t size) {
00096       return 0;
00097     }
00098     
00099     /// Give the value of the function \f$ y(x=x_0) \f$ .
00100     virtual int interp(const vec_t &x, const vec_t &y, size_t size, 
00101                        double x0, double &y0) {
00102       return 0;
00103     }
00104     
00105     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
00106     virtual int deriv(const vec_t &x, const vec_t &y, size_t size,
00107                       double x0, double &dydx) {
00108       return 0;
00109     }
00110     
00111     /** \brief Give the value of the second derivative 
00112         \f$ y^{\prime \prime}(x=x_0) \f$ .
00113     */
00114     virtual int deriv2(const vec_t &x, const vec_t &y, size_t size, 
00115                        double x0, double &d2ydx2) {
00116       return 0;
00117     }
00118     
00119     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
00120     virtual int integ(const vec_t &x, const vec_t &y, size_t size, 
00121                       double a, double b, double &result) {
00122       return 0;
00123     }
00124   };
00125   
00126   /** 
00127       \brief Linear interpolation (GSL)
00128 
00129       Linear interpolation requires no calls to alloc(), free() or
00130       init(), as there is no internal storage required.
00131   */
00132   template<class vec_t> class linear_interp : public base_interp<vec_t> {
00133     
00134   public:
00135 
00136     linear_interp() {
00137       this->min_size=2;
00138     }
00139 
00140     /// Give the value of the function \f$ y(x=x_0) \f$ .
00141     virtual int interp(const vec_t &x_array, const vec_t &y_array, 
00142                        size_t size, double x, double &y) {
00143       double x_lo, x_hi;
00144       double y_lo, y_hi;
00145       double dx;
00146       size_t index;
00147 
00148       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00149       
00150       /* evaluate */
00151       x_lo = x_array[index];
00152       x_hi = x_array[index + 1];
00153       y_lo = y_array[index];
00154       y_hi = y_array[index + 1];
00155       dx = x_hi - x_lo;
00156       if (dx > 0.0)
00157         {
00158           y = y_lo + (x - x_lo) / dx * (y_hi - y_lo);
00159           return gsl_success;
00160         }
00161       else
00162         {
00163           y = 0.0;
00164           return GSL_EINVAL;
00165         }
00166     }
00167     
00168     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
00169     virtual int deriv(const vec_t &x_array, const vec_t &y_array, 
00170                       size_t size,
00171                       double x, double &dydx) {
00172       
00173       double x_lo, x_hi;
00174       double y_lo, y_hi;
00175       double dx;
00176       double dy;
00177       size_t index;
00178       
00179       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00180       
00181       /* evaluate */
00182       x_lo = x_array[index];
00183       x_hi = x_array[index + 1];
00184       y_lo = y_array[index];
00185       y_hi = y_array[index + 1];
00186       dx = x_hi - x_lo;
00187       dy = y_hi - y_lo;
00188       if (dx > 0.0)
00189         {
00190           dydx = dy / dx;;
00191           return gsl_success;
00192         }
00193       else
00194         {
00195           dydx = 0.0;
00196           return GSL_EINVAL;
00197         }
00198     }
00199 
00200     /** \brief Give the value of the second derivative  
00201         \f$ y^{\prime \prime}(x=x_0) \f$ .
00202     */
00203     virtual int deriv2(const vec_t &x, const vec_t &y, size_t size, 
00204                        double x0, double &d2ydx2) {
00205       d2ydx2=0.0;
00206       return gsl_success;
00207     }
00208 
00209 
00210     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
00211     virtual int integ(const vec_t &x_array, const vec_t &y_array, 
00212                       size_t size, double a, double b, 
00213                       double &result) {
00214       size_t i, index_a, index_b;
00215       
00216       index_a=this->sv.bsearch_inc(a,x_array,0,size-1);
00217       index_b=this->sv.bsearch_inc(b,x_array,0,size-1);
00218       
00219       /* endpoints span more than one interval */
00220       
00221       result = 0.0;
00222       
00223       /* interior intervals */
00224       for(i=index_a; i<=index_b; i++) {
00225         const double x_hi = x_array[i + 1];
00226         const double x_lo = x_array[i];
00227         const double y_lo = y_array[i];
00228         const double y_hi = y_array[i + 1];
00229         const double dx = x_hi - x_lo;
00230         
00231         if(dx != 0.0) {
00232           if (i == index_a || i == index_b)
00233             {
00234               double x1 = (i == index_a) ? a : x_lo;
00235               double x2 = (i == index_b) ? b : x_hi;
00236               const double D = (y_hi-y_lo)/dx;
00237               result += (x2-x1) * (y_lo + 0.5*D*((x2-x_lo)+(x1-x_lo)));
00238             }
00239           else
00240             {
00241               result += 0.5 * dx * (y_lo + y_hi);
00242             }
00243         }
00244       }
00245       
00246       return gsl_success;
00247     }
00248   };
00249   
00250   /** 
00251       \brief Cubic spline interpolation (GSL)
00252 
00253   */
00254   template<class vec_t> class cspline_interp : public base_interp<vec_t> {
00255     
00256 #ifndef DOXYGEN_INTERNAL
00257     
00258   protected:
00259 
00260     /// \name Storage for cubic spline interpolation
00261     //@{
00262     double *c;
00263     double *g;
00264     double *diag;
00265     double *offdiag;
00266     //@}
00267 
00268     /// True for periodic boundary conditions
00269     bool peri;
00270 
00271     /// Compute coefficients for cubic spline interpolation
00272     void coeff_calc(const double c_array[], double dy, double dx, 
00273                     size_t index,  double * b, double * c2, double * d) {
00274       const double c_i = c_array[index];
00275       const double c_ip1 = c_array[index + 1];
00276       *b = (dy / dx) - dx * (c_ip1 + 2.0 * c_i) / 3.0;
00277       *c2 = c_i;
00278       *d = (c_ip1 - c_i) / (3.0 * dx);
00279     }
00280 
00281 #endif
00282     
00283   public:
00284 
00285     /** \brief Create a base interpolation object with natural or
00286         periodic boundary conditions
00287     */
00288     cspline_interp(bool periodic=false) {
00289       
00290       peri=periodic;
00291       if (peri) this->min_size=2;
00292       else this->min_size=3;
00293       
00294     }
00295 
00296     virtual ~cspline_interp() {}
00297     
00298     /// Allocate memory, assuming x and y have size \c size
00299     virtual int alloc(size_t size) {
00300       
00301       c = (double *) malloc (size * sizeof (double));
00302       
00303       if ( c == NULL)
00304         {
00305           GSL_ERROR_NULL("failed to allocate space for c", GSL_ENOMEM);
00306         }
00307       
00308       g = (double *) malloc (size * sizeof (double));
00309       
00310       if ( g == NULL)
00311         {
00312           std::free(c);
00313           GSL_ERROR_NULL("failed to allocate space for g", GSL_ENOMEM);
00314         }
00315       
00316       diag = (double *) malloc (size * sizeof (double));
00317       
00318       if ( diag == NULL)
00319         {
00320           std::free(g);
00321           std::free(c);
00322           GSL_ERROR_NULL("failed to allocate space for diag", GSL_ENOMEM);
00323         }
00324       
00325       offdiag = (double *) malloc (size * sizeof (double));
00326       
00327       if ( offdiag == NULL)
00328         {
00329           std::free(diag);
00330           std::free(g);
00331           std::free(c);
00332           GSL_ERROR_NULL("failed to allocate space for offdiag", GSL_ENOMEM);
00333         }
00334       
00335       return 0;
00336     }
00337     
00338     /// Initialize interpolation routine
00339     virtual int init(const vec_t &xa, const vec_t &ya, size_t size) {
00340 
00341       if (peri) {
00342 
00343         /// Periodic boundary conditions
00344          
00345         size_t i;
00346         size_t num_points = size;
00347         /* Engeln-Mullges + Uhlig "n" */
00348         size_t max_index = num_points - 1;  
00349         /* linear system is sys_size x sys_size */
00350         size_t sys_size = max_index;    
00351         
00352         if (sys_size == 2) {
00353           /* solve 2x2 system */
00354           
00355           const double h0 = xa[1] - xa[0];
00356           const double h1 = xa[2] - xa[1];
00357           const double h2 = xa[3] - xa[2];
00358           const double A = 2.0*(h0 + h1);
00359           const double B = h0 + h1;
00360           double gx[2];
00361           double det;
00362           
00363           gx[0] = 3.0 * ((ya[2] - ya[1]) / h1 - (ya[1] - ya[0]) / h0);
00364           gx[1] = 3.0 * ((ya[1] - ya[2]) / h2 - (ya[2] - ya[1]) / h1);
00365           
00366           det = 3.0 * (h0 + h1) * (h0 + h1);
00367           c[1] = ( A * gx[0] - B * gx[1])/det;
00368           c[2] = (-B * gx[0] + A * gx[1])/det;
00369           c[0] =  c[2];
00370           
00371           return gsl_success;
00372 
00373         } else {
00374           
00375           for (i = 0; i < sys_size-1; i++) {
00376             const double h_i       = xa[i + 1] - xa[i];
00377             const double h_ip1     = xa[i + 2] - xa[i + 1];
00378             const double ydiff_i   = ya[i + 1] - ya[i];
00379             const double ydiff_ip1 = ya[i + 2] - ya[i + 1];
00380             const double g_i = (h_i != 0.0) ? 1.0 / h_i : 0.0;
00381             const double g_ip1 = (h_ip1 != 0.0) ? 1.0 / h_ip1 : 0.0;
00382             offdiag[i] = h_ip1;
00383             diag[i] = 2.0 * (h_ip1 + h_i);
00384             g[i] = 3.0 * (ydiff_ip1 * g_ip1 - ydiff_i * g_i);
00385           }
00386           
00387           i = sys_size - 1;
00388           {
00389             const double h_i       = xa[i + 1] - xa[i];
00390             const double h_ip1     = xa[1] - xa[0];
00391             const double ydiff_i   = ya[i + 1] - ya[i];
00392             const double ydiff_ip1 = ya[1] - ya[0];
00393             const double g_i = (h_i != 0.0) ? 1.0 / h_i : 0.0;
00394             const double g_ip1 = (h_ip1 != 0.0) ? 1.0 / h_ip1 : 0.0;
00395             offdiag[i] = h_ip1;
00396             diag[i] = 2.0 * (h_ip1 + h_i);
00397             g[i] = 3.0 * (ydiff_ip1 * g_ip1 - ydiff_i * g_i);
00398           }
00399           
00400           {
00401             gsl_vector_view g_vec = gsl_vector_view_array( g, sys_size);
00402             gsl_vector_view diag_vec = 
00403               gsl_vector_view_array( diag,sys_size);
00404             gsl_vector_view offdiag_vec = 
00405               gsl_vector_view_array( offdiag, sys_size);
00406             gsl_vector_view solution_vec = 
00407               gsl_vector_view_array (( c) + 1, sys_size);
00408             int status = 
00409               gsl_linalg_solve_symm_cyc_tridiag(&diag_vec.vector,
00410                                                 &offdiag_vec.vector,
00411                                                 &g_vec.vector,
00412                                                 &solution_vec.vector);
00413             c[0] =  c[max_index];
00414             
00415             return status;
00416           }
00417         }
00418 
00419       } else {
00420 
00421         /// Natural boundary conditions
00422 
00423         size_t i;
00424         size_t num_points = size;
00425         size_t max_index = num_points - 1;
00426         size_t sys_size = max_index - 1;
00427 
00428         c[0] = 0.0;
00429         c[max_index] = 0.0;
00430 
00431         for (i = 0; i < sys_size; i++)
00432           {
00433             const double h_i   = xa[i + 1] - xa[i];
00434             const double h_ip1 = xa[i + 2] - xa[i + 1];
00435             const double ydiff_i   = ya[i + 1] - ya[i];
00436             const double ydiff_ip1 = ya[i + 2] - ya[i + 1];
00437             const double g_i = (h_i != 0.0) ? 1.0 / h_i : 0.0;
00438             const double g_ip1 = (h_ip1 != 0.0) ? 1.0 / h_ip1 : 0.0;
00439             offdiag[i] = h_ip1;
00440             diag[i] = 2.0 * (h_ip1 + h_i);
00441             g[i] = 3.0 * (ydiff_ip1 * g_ip1 -  ydiff_i * g_i);
00442           }
00443 
00444         if (sys_size == 1)
00445           {
00446             c[1] =  g[0] /  diag[0];
00447             return gsl_success;
00448           }
00449         else
00450           {
00451             gsl_vector_view g_vec = gsl_vector_view_array( g, sys_size);
00452             gsl_vector_view diag_vec = gsl_vector_view_array( diag, sys_size);
00453             gsl_vector_view offdiag_vec = 
00454               gsl_vector_view_array( offdiag, sys_size - 1);
00455             gsl_vector_view solution_vec = 
00456               gsl_vector_view_array (( c) + 1, sys_size);
00457             
00458             int status = gsl_linalg_solve_symm_tridiag(&diag_vec.vector, 
00459                                                        &offdiag_vec.vector, 
00460                                                        &g_vec.vector, 
00461                                                        &solution_vec.vector);
00462             return status;
00463           }
00464       }
00465     }
00466       
00467     
00468     /// Free allocated memory
00469     virtual int free() {
00470       std::free(c);
00471       std::free(g);
00472       std::free(diag);
00473       std::free(offdiag);
00474 
00475       return 0;
00476     }
00477 
00478     /// Give the value of the function \f$ y(x=x_0) \f$ .
00479     virtual int interp(const vec_t &x_array, const vec_t &y_array, 
00480                        size_t size, double x, double &y) {
00481       
00482       double x_lo, x_hi;
00483       double dx;
00484       size_t index;
00485       
00486       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00487   
00488       /* evaluate */
00489       x_hi = x_array[index + 1];
00490       x_lo = x_array[index];
00491       dx = x_hi - x_lo;
00492       if (dx > 0.0)
00493         {
00494           const double y_lo = y_array[index];
00495           const double y_hi = y_array[index + 1];
00496           const double dy = y_hi - y_lo;
00497           double delx = x - x_lo;
00498           double b_i, c_i, d_i; 
00499           coeff_calc( c, dy, dx, index,  &b_i, &c_i, &d_i);
00500           y = y_lo + delx * (b_i + delx * (c_i + delx * d_i));
00501           return gsl_success;
00502         }
00503       else
00504         {
00505           y = 0.0;
00506           return GSL_EINVAL;
00507         }
00508     }
00509 
00510 
00511     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
00512     virtual int deriv(const vec_t &x_array, const vec_t &y_array, 
00513                       size_t size, double x, double &dydx) {
00514 
00515       double x_lo, x_hi;
00516       double dx;
00517       size_t index;
00518   
00519       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00520   
00521       /* evaluate */
00522       x_hi = x_array[index + 1];
00523       x_lo = x_array[index];
00524       dx = x_hi - x_lo;
00525       if (dx > 0.0)
00526         {
00527           const double y_lo = y_array[index];
00528           const double y_hi = y_array[index + 1];
00529           const double dy = y_hi - y_lo;
00530           double delx = x - x_lo;
00531           double b_i, c_i, d_i; 
00532           coeff_calc( c, dy, dx, index,  &b_i, &c_i, &d_i);
00533           dydx = b_i + delx * (2.0 * c_i + 3.0 * d_i * delx);
00534           return gsl_success;
00535         }
00536       else
00537         {
00538           dydx = 0.0;
00539           return GSL_FAILURE;
00540         }
00541     }
00542 
00543     /** \brief Give the value of the second derivative  
00544         \f$ y^{\prime \prime}(x=x_0) \f$ .
00545     */
00546     virtual int deriv2(const vec_t &x_array, const vec_t &y_array, 
00547                        size_t size, double x, double &d2ydx2) {
00548 
00549       double x_lo, x_hi;
00550       double dx;
00551       size_t index;
00552   
00553       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00554   
00555       /* evaluate */
00556       x_hi = x_array[index + 1];
00557       x_lo = x_array[index];
00558       dx = x_hi - x_lo;
00559       if (dx > 0.0)
00560         {
00561           const double y_lo = y_array[index];
00562           const double y_hi = y_array[index + 1];
00563           const double dy = y_hi - y_lo;
00564           double delx = x - x_lo;
00565           double b_i, c_i, d_i;
00566           coeff_calc( c, dy, dx, index,  &b_i, &c_i, &d_i);
00567           d2ydx2 = 2.0 * c_i + 6.0 * d_i * delx;
00568           return gsl_success;
00569         }
00570       else
00571         {
00572           d2ydx2 = 0.0;
00573           return GSL_FAILURE;
00574         }
00575     }
00576 
00577 #ifndef DOXYGEN_INTERNAL
00578 
00579   protected:
00580 
00581     /// An internal function to assist in computing the integral
00582     double integ_eval(double ai, double bi, double ci, double di, double xi, 
00583                       double a, double b) {
00584       const double r1 = a - xi;
00585       const double r2 = b - xi;
00586       const double r12 = r1 + r2;
00587       const double bterm = 0.5 * bi * r12;
00588       const double cterm = (1.0 / 3.0) * ci * (r1 * r1 + r2 * r2 + r1 * r2);
00589       const double dterm = 0.25 * di * r12 * (r1 * r1 + r2 * r2);
00590       
00591       return (b - a) * (ai + bterm + cterm + dterm);
00592     }
00593 
00594 #endif
00595     
00596   public:
00597     
00598     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
00599     virtual int integ(const vec_t &x_array, const vec_t &y_array, 
00600                       size_t size, double a, double b, 
00601                       double &result) {
00602 
00603       size_t i, index_a, index_b;
00604   
00605       index_a=this->sv.bsearch_inc(a,x_array,0,size-1);
00606       index_b=this->sv.bsearch_inc(b,x_array,0,size-1);
00607 
00608       result = 0.0;
00609   
00610       /* interior intervals */
00611       for(i=index_a; i<=index_b; i++) {
00612         const double x_hi = x_array[i + 1];
00613         const double x_lo = x_array[i];
00614         const double y_lo = y_array[i];
00615         const double y_hi = y_array[i + 1];
00616         const double dx = x_hi - x_lo;
00617         const double dy = y_hi - y_lo;
00618         if(dx != 0.0) {
00619           double b_i, c_i, d_i; 
00620           coeff_calc( c, dy, dx, i,  &b_i, &c_i, &d_i);
00621       
00622           if (i == index_a || i == index_b)
00623             {
00624               double x1 = (i == index_a) ? a : x_lo;
00625               double x2 = (i == index_b) ? b : x_hi;
00626               result += integ_eval(y_lo, b_i, c_i, d_i, x_lo, x1, x2);
00627             }
00628           else
00629             {
00630               result += dx * (y_lo + dx*(0.5*b_i + 
00631                                          dx*(c_i/3.0 + 0.25*d_i*dx)));
00632             }
00633         }
00634         else {
00635           result = 0.0;
00636           return GSL_FAILURE;
00637         }
00638       }
00639   
00640       return gsl_success;
00641     }
00642 
00643   };
00644   
00645   /** 
00646       \brief Interpolation class
00647       
00648       This interpolation class is intended to be sufficiently general
00649       to handle any vector type. Interpolation of \ref ovector like
00650       objects is performed with the default template parameters, and
00651       \ref array_interp is provided for simple interpolation on
00652       C-style arrays.
00653 
00654       The class automatically handles decreasing arrays by converting
00655       from an object of type \c vec_t to an object ov type \c rvec_t.
00656 
00657       While \c vec_t may be any vector type which allows indexing via
00658       \c [], \c rvec_t must be a vector type which allows indexing
00659       and has a constructor with one of the two forms
00660       \verbatim
00661       rvec_t::rvec_t(vec_t &v);
00662       rvec_t::rvec_t(vec_t v);
00663       \endverbatim
00664       so that o2scl_interp can automatically "reverse" a vector if 
00665       necessary.
00666       
00667       It is important that different instances of o2scl_interp_vec
00668       and o2scl_interp not be given the same interpolation objects,
00669       as they will clash.
00670       
00671   */
00672   template<class vec_t=ovector_view, class rvec_t=ovector_const_reverse> 
00673     class o2scl_interp {
00674       
00675 #ifndef DOXYGEN_INTERNAL
00676       
00677     protected:
00678       
00679     /// Pointer to base interpolation object
00680     base_interp<vec_t> *itp;
00681       
00682     /// Pointer to base interpolation object for reversed vectors
00683     base_interp<rvec_t> *ritp;
00684       
00685 #endif
00686       
00687     public:
00688       
00689     /// Create with base interpolation objects \c it and \c rit
00690     o2scl_interp(base_interp<vec_t> &it, base_interp<rvec_t> &rit) {
00691       itp=&it;
00692       ritp=&rit;
00693     }
00694 
00695     /** \brief Create with base interpolation object \c it and use
00696         \ref def_ritp for reverse interpolation if necessary
00697     */
00698     o2scl_interp(base_interp<vec_t> &it) {
00699       itp=&it;
00700       ritp=&def_ritp;
00701     }
00702     
00703     /** \brief Create an interpolator using \ref def_itp and 
00704         \ref def_ritp
00705     */
00706     o2scl_interp() {
00707       itp=&def_itp;
00708       ritp=&def_ritp;
00709     }
00710 
00711     virtual ~o2scl_interp() {}
00712 
00713     /// Give the value of the function \f$ y(x=x_0) \f$ .
00714     virtual double interp(const double x0, size_t n, const vec_t &x, 
00715                           const vec_t &y) 
00716     {
00717         
00718       double ret=0.0;
00719       
00720       err_hnd->reset();
00721       
00722       if (n<itp->min_size) {
00723         set_err("Vector size too small in interp().",gsl_edom);
00724       }
00725       
00726       if (x[0]>x[n-1]) {
00727 
00728         // Decreasing case
00729 
00730         rvec_t ax(x);
00731         rvec_t ay(y);
00732 
00733         int r1=ritp->alloc(n);
00734         if (r1==0) {
00735           int r2=ritp->init(ax,ay,n);
00736           if (r2==0) {
00737             ritp->interp(ax,ay,n,x0,ret);
00738             ritp->free();
00739           }
00740         }
00741 
00742       } else {
00743 
00744         // Increasing case
00745           
00746         int r1=itp->alloc(n);
00747         if (r1==0) {
00748           int r2=itp->init(x,y,n);
00749           if (r2==0) {
00750             itp->interp(x,y,n,x0,ret);
00751             itp->free();
00752           }
00753         }
00754           
00755       }
00756       
00757       return ret;
00758     }                 
00759     
00760     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
00761     virtual double deriv(const double x0, size_t n, const vec_t &x, 
00762                          const vec_t &y) 
00763     {
00764       
00765       double ret=0.0;
00766       
00767       err_hnd->reset();
00768       
00769       if (n<itp->min_size) {
00770         set_err("Vector size too small in interp().",gsl_edom);
00771       }
00772       
00773       if (x[0]>x[n-1]) {
00774           
00775         rvec_t ax(x);
00776         rvec_t ay(y);
00777 
00778         int r1=ritp->alloc(n);
00779         if (r1==0) {
00780           int r2=ritp->init(ax,ay,n);
00781           if (r2==0) {
00782             ritp->deriv(ax,ay,n,x0,ret);
00783             ritp->free();
00784           }
00785         }
00786           
00787       } else {
00788 
00789         // Increasing case
00790           
00791         int r1=itp->alloc(n);
00792         if (r1==0) {
00793           int r2=itp->init(x,y,n);
00794           if (r2==0) {
00795             itp->deriv(x,y,n,x0,ret);
00796             itp->free();
00797           }
00798         }
00799       }
00800       
00801       return ret;
00802     }                 
00803     
00804     /** \brief Give the value of the second derivative  
00805         \f$ y^{\prime \prime}(x=x_0) \f$ .
00806     */
00807     virtual double deriv2(const double x0, size_t n, const vec_t &x, 
00808                           const vec_t &y) 
00809     {
00810       
00811       double ret=0.0;
00812       
00813       err_hnd->reset();
00814       
00815       if (n<itp->min_size) {
00816 
00817         set_err("Vector size too small in deriv2().",gsl_edom);
00818 
00819       }
00820       
00821       if (x[0]>x[n-1]) {
00822           
00823         // Decreasing case
00824             
00825         rvec_t ax(x);
00826         rvec_t ay(y);
00827 
00828         int r1=ritp->alloc(n);
00829         if (r1==0) {
00830           int r2=ritp->init(ax,ay,n);
00831           if (r2==0) {
00832             ritp->deriv2(ax,ay,n,x0,ret);
00833             ritp->free();
00834           }
00835         }
00836 
00837       } else {
00838 
00839         // Increasing case
00840 
00841         int r1=itp->alloc(n);
00842         if (r1==0) {
00843           int r2=itp->init(x,y,n);
00844           if (r2==0) {
00845             itp->deriv2(x,y,n,x0,ret);
00846             itp->free();
00847           }
00848         }
00849 
00850       }
00851       
00852       return ret;
00853     }                 
00854     
00855     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
00856     virtual double integ(const double x1, const double x2, size_t n, 
00857                          const vec_t &x, const vec_t &y) 
00858     {
00859       
00860       double ret=0.0;
00861       
00862       err_hnd->reset();
00863       
00864       if (n<itp->min_size) {
00865 
00866         set_err("Vector size too small in integ().",gsl_edom);
00867 
00868       }
00869       
00870       if (x[0]>x[n-1]) {
00871           
00872         // Decreasing case
00873             
00874         rvec_t ax(x);
00875         rvec_t ay(y);
00876 
00877         int r1=ritp->alloc(n);
00878         if (r1==0) {
00879           int r2=ritp->init(ax,ay,n);
00880           if (r2==0) {
00881             ritp->integ(ax,ay,n,x1,x2,ret);
00882             ritp->free();
00883           }
00884         }
00885 
00886       } else {
00887 
00888         // Increasing case
00889             
00890         int r1=itp->alloc(n);
00891         if (r1==0) {
00892           int r2=itp->init(x,y,n);
00893           if (r2==0) {
00894             itp->integ(x,y,n,x1,x2,ret);
00895             itp->free();
00896           }
00897         }
00898 
00899       }
00900       
00901       return ret;
00902     }                 
00903     
00904     /// Set base interpolation object
00905     int set_type(base_interp<vec_t> &it, base_interp<rvec_t> &rit) {
00906       itp=&it;
00907       ritp=&rit;
00908       return 0;
00909     }
00910       
00911     /** \brief Default base interpolation object (cubic spline with natural 
00912         boundary conditions)
00913     */
00914     cspline_interp<vec_t> def_itp;
00915 
00916     /** \brief Default base interpolation object for reversed vectors
00917         (cubic spline with natural boundary conditions)
00918     */
00919     cspline_interp<rvec_t> def_ritp;
00920       
00921   };
00922   
00923   /** 
00924       \brief Interpolation class for pre-specified vector
00925       
00926       This interpolation class is intended to be sufficiently general
00927       to handle any vector type. Interpolation of \ref ovector like
00928       objects is performed with the default template parameters, and
00929       \ref array_interp_vec is provided for simple interpolation on
00930       C-style arrays.
00931 
00932       The class automatically handles decreasing arrays by copying the
00933       old array to a reversed version. For several interpolations on
00934       the same data, copying the initial array can be faster than
00935       overloading operator[].
00936       
00937       It is important that different instances of o2scl_interp_vec
00938       and o2scl_interp not be given the same interpolation objects,
00939       as they will clash.
00940 
00941       \todo Need to fix constructor to behave properly if init() fails.
00942       It should free the memory and set \c ln to zero.
00943   */
00944   template<class vec_t=ovector_view, class alloc_vec_t=ovector, 
00945     class alloc_t=ovector_alloc> class o2scl_interp_vec {
00946       
00947 #ifndef DOXYGEN_INTERNAL
00948       
00949     protected:
00950       
00951     /// Memory allocator for objects of type \c alloc_vec_t
00952     alloc_t ao;
00953 
00954     /// Pointer to base interpolation object
00955     base_interp<vec_t> *itp;
00956       
00957     /// True if the original array was increasing
00958     bool inc;
00959 
00960     /// Pointer to the user-specified x vector
00961     const vec_t *lx;
00962 
00963     /// Pointer to the user-specified x vector
00964     const vec_t *ly;
00965 
00966     /// Reversed version of x
00967     alloc_vec_t lrx;
00968 
00969     /// Reversed version of y
00970     alloc_vec_t lry;
00971 
00972     /// Size of user-specified vectors
00973     size_t ln;
00974 
00975 #endif
00976       
00977     public:
00978     
00979     /// Create with base interpolation object \c it
00980     o2scl_interp_vec(base_interp<vec_t> &it, 
00981                    size_t n, const vec_t &x, const vec_t &y) {
00982       itp=&it;
00983       ln=0;
00984         
00985       if (n<itp->min_size) {
00986         set_err("Vector size too small in o2scl_interp_vec().",gsl_edom);
00987       } else {
00988         size_t iend=n-1;
00989           
00990         if (x[0]>x[iend]) {
00991             
00992           ao.allocate(lrx,n);
00993           ao.allocate(lry,n);
00994 
00995           for(size_t i=0;i<n;i++) {
00996             lrx[iend-i]=x[i];
00997             lry[iend-i]=y[i];
00998           }
00999             
01000           int r1=itp->alloc(n);
01001           if (r1==0) {
01002             itp->init(lrx,lry,n);
01003             ln=n;
01004             inc=false;
01005           }
01006             
01007         } else {
01008             
01009           int r1=itp->alloc(n);
01010           if (r1==0) {
01011             itp->init(x,y,n);
01012             ln=n;
01013             lx=&x;
01014             ly=&y;
01015             inc=true;
01016           }
01017             
01018         }
01019 
01020       } 
01021         
01022     }
01023       
01024     virtual ~o2scl_interp_vec() {
01025       if (ln>0) {
01026         itp->free();
01027       }
01028       if (inc==false) {
01029         ao.free(lrx);
01030         ao.free(lry);
01031       }
01032     }
01033     
01034     /// Give the value of the function \f$ y(x=x_0) \f$ .
01035     virtual double interp(const double x0) {
01036       double ret=0.0;
01037       if (ln>0) {
01038         if (inc) itp->interp(*lx,*ly,ln,x0,ret);
01039         else itp->interp(lrx,lry,ln,x0,ret);
01040       }
01041       return ret;
01042     }                 
01043     
01044     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
01045     virtual double deriv(const double x0) {
01046       double ret=0.0;
01047       if (ln>0) {
01048         if (inc) itp->deriv(*lx,*ly,ln,x0,ret);
01049         else itp->deriv(lrx,lry,ln,x0,ret);
01050       }
01051       return ret;
01052     }                 
01053     
01054     /** \brief Give the value of the second derivative  
01055         \f$ y^{\prime \prime}(x=x_0) \f$ .
01056     */
01057     virtual double deriv2(const double x0) {
01058       double ret=0.0;
01059       if (ln>0) {
01060         if (inc) itp->deriv2(*lx,*ly,ln,x0,ret);
01061         else itp->deriv2(lrx,lry,ln,x0,ret);
01062       }
01063       return ret;
01064     }                 
01065     
01066     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
01067     virtual double integ(const double x1, const double x2) {
01068       double ret=0.0;
01069       if (ln>0) {
01070         if (inc) itp->integ(*lx,*ly,ln,x1,x2,ret);
01071         else itp->integ(lrx,lry,ln,x1,x2,ret);
01072       }
01073       return ret;
01074     }                 
01075     
01076     /// Set base interpolation object
01077     int set_type(base_interp<vec_t> &it) {
01078       itp=&it;
01079       return 0;
01080     }
01081       
01082     /** \brief Default base interpolation object (cubic spline with natural 
01083         boundary conditions)
01084     */
01085     cspline_interp<vec_t> def_itp;
01086       
01087   };
01088   
01089   /** \brief A specialization of o2scl_interp for C-style double arrays
01090    */
01091   template<size_t n> class array_interp :
01092   public o2scl_interp<double[n],array_const_reverse<n> > {
01093     
01094   public:
01095     
01096     /// Create with base interpolation objects \c it and \c rit
01097     array_interp(base_interp<double[n]> &it, 
01098                  base_interp<array_const_reverse<n> > &rit)
01099       : o2scl_interp<double[n],array_const_reverse<n> >(it,rit) {
01100     }
01101     
01102       /** \brief Create with base interpolation object \c it and use
01103           the default base interpolation object for reversed arrays
01104       */
01105       array_interp(base_interp<double[n]> &it) : 
01106         o2scl_interp<double[n],array_const_reverse<n> >(it) {
01107       }
01108       
01109         /** \brief Create an interpolator using the default base
01110             interpolation objects
01111         */
01112         array_interp() : 
01113           o2scl_interp<double[n],array_const_reverse<n> >() {
01114         }
01115 
01116   };
01117 
01118   /** 
01119       \brief A specialization of o2scl_interp_vec for C-style arrays
01120   */
01121   template<class arr_t> class array_interp_vec : 
01122   public o2scl_interp_vec<arr_t,arr_t,array_alloc<arr_t> > {
01123     
01124   public:
01125     
01126     /// Create with base interpolation object \c it
01127     array_interp_vec(base_interp<arr_t> &it, 
01128                      size_t nv, const arr_t &x, const arr_t &y) :
01129       o2scl_interp_vec<arr_t,arr_t,array_alloc<arr_t> >(it,nv,x,y) {
01130     }
01131     
01132   };
01133   
01134 #ifndef DOXYGENP
01135 }
01136 #endif
01137 
01138 #endif

Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.

Project hosting provided by SourceForge.net Logo, O2scl Sourceforge Project Page