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 allocate() 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 allocate(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 allocate(), 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   template<class vec_t> class cspline_interp : public base_interp<vec_t> {
00254     
00255 #ifndef DOXYGEN_INTERNAL
00256     
00257   protected:
00258 
00259     /// \name Storage for cubic spline interpolation
00260     //@{
00261     double *c;
00262     double *g;
00263     double *diag;
00264     double *offdiag;
00265     //@}
00266 
00267     /// True for periodic boundary conditions
00268     bool peri;
00269 
00270     /// Compute coefficients for cubic spline interpolation
00271     void coeff_calc(const double c_array[], double dy, double dx, 
00272                     size_t index,  double * b, double * c2, double * d) {
00273       const double c_i = c_array[index];
00274       const double c_ip1 = c_array[index + 1];
00275       *b = (dy / dx) - dx * (c_ip1 + 2.0 * c_i) / 3.0;
00276       *c2 = c_i;
00277       *d = (c_ip1 - c_i) / (3.0 * dx);
00278     }
00279 
00280 #endif
00281     
00282   public:
00283 
00284     /** \brief Create a base interpolation object with natural or
00285         periodic boundary conditions
00286     */
00287     cspline_interp(bool periodic=false) {
00288       
00289       peri=periodic;
00290       if (peri) this->min_size=2;
00291       else this->min_size=3;
00292       
00293     }
00294 
00295     virtual ~cspline_interp() {}
00296     
00297     /// Allocate memory, assuming x and y have size \c size
00298     virtual int allocate(size_t size) {
00299       
00300       c = (double *) malloc (size * sizeof (double));
00301       
00302       if ( c == NULL)
00303         {
00304           GSL_ERROR_NULL("failed to allocate space for c", GSL_ENOMEM);
00305         }
00306       
00307       g = (double *) malloc (size * sizeof (double));
00308       
00309       if ( g == NULL)
00310         {
00311           std::free(c);
00312           GSL_ERROR_NULL("failed to allocate space for g", GSL_ENOMEM);
00313         }
00314       
00315       diag = (double *) malloc (size * sizeof (double));
00316       
00317       if ( diag == NULL)
00318         {
00319           std::free(g);
00320           std::free(c);
00321           GSL_ERROR_NULL("failed to allocate space for diag", GSL_ENOMEM);
00322         }
00323       
00324       offdiag = (double *) malloc (size * sizeof (double));
00325       
00326       if ( offdiag == NULL)
00327         {
00328           std::free(diag);
00329           std::free(g);
00330           std::free(c);
00331           GSL_ERROR_NULL("failed to allocate space for offdiag", GSL_ENOMEM);
00332         }
00333       
00334       return 0;
00335     }
00336     
00337     /// Initialize interpolation routine
00338     virtual int init(const vec_t &xa, const vec_t &ya, size_t size) {
00339 
00340       if (peri) {
00341 
00342         /// Periodic boundary conditions
00343          
00344         size_t i;
00345         size_t num_points = size;
00346         /* Engeln-Mullges + Uhlig "n" */
00347         size_t max_index = num_points - 1;  
00348         /* linear system is sys_size x sys_size */
00349         size_t sys_size = max_index;    
00350         
00351         if (sys_size == 2) {
00352           /* solve 2x2 system */
00353           
00354           const double h0 = xa[1] - xa[0];
00355           const double h1 = xa[2] - xa[1];
00356           const double h2 = xa[3] - xa[2];
00357           const double A = 2.0*(h0 + h1);
00358           const double B = h0 + h1;
00359           double gx[2];
00360           double det;
00361           
00362           gx[0] = 3.0 * ((ya[2] - ya[1]) / h1 - (ya[1] - ya[0]) / h0);
00363           gx[1] = 3.0 * ((ya[1] - ya[2]) / h2 - (ya[2] - ya[1]) / h1);
00364           
00365           det = 3.0 * (h0 + h1) * (h0 + h1);
00366           c[1] = ( A * gx[0] - B * gx[1])/det;
00367           c[2] = (-B * gx[0] + A * gx[1])/det;
00368           c[0] =  c[2];
00369           
00370           return gsl_success;
00371 
00372         } else {
00373           
00374           for (i = 0; i < sys_size-1; i++) {
00375             const double h_i       = xa[i + 1] - xa[i];
00376             const double h_ip1     = xa[i + 2] - xa[i + 1];
00377             const double ydiff_i   = ya[i + 1] - ya[i];
00378             const double ydiff_ip1 = ya[i + 2] - ya[i + 1];
00379             const double g_i = (h_i != 0.0) ? 1.0 / h_i : 0.0;
00380             const double g_ip1 = (h_ip1 != 0.0) ? 1.0 / h_ip1 : 0.0;
00381             offdiag[i] = h_ip1;
00382             diag[i] = 2.0 * (h_ip1 + h_i);
00383             g[i] = 3.0 * (ydiff_ip1 * g_ip1 - ydiff_i * g_i);
00384           }
00385           
00386           i = sys_size - 1;
00387           {
00388             const double h_i       = xa[i + 1] - xa[i];
00389             const double h_ip1     = xa[1] - xa[0];
00390             const double ydiff_i   = ya[i + 1] - ya[i];
00391             const double ydiff_ip1 = ya[1] - ya[0];
00392             const double g_i = (h_i != 0.0) ? 1.0 / h_i : 0.0;
00393             const double g_ip1 = (h_ip1 != 0.0) ? 1.0 / h_ip1 : 0.0;
00394             offdiag[i] = h_ip1;
00395             diag[i] = 2.0 * (h_ip1 + h_i);
00396             g[i] = 3.0 * (ydiff_ip1 * g_ip1 - ydiff_i * g_i);
00397           }
00398           
00399           {
00400             gsl_vector_view g_vec = gsl_vector_view_array( g, sys_size);
00401             gsl_vector_view diag_vec = 
00402               gsl_vector_view_array( diag,sys_size);
00403             gsl_vector_view offdiag_vec = 
00404               gsl_vector_view_array( offdiag, sys_size);
00405             gsl_vector_view solution_vec = 
00406               gsl_vector_view_array (( c) + 1, sys_size);
00407             int status = 
00408               gsl_linalg_solve_symm_cyc_tridiag(&diag_vec.vector,
00409                                                 &offdiag_vec.vector,
00410                                                 &g_vec.vector,
00411                                                 &solution_vec.vector);
00412             c[0] =  c[max_index];
00413             
00414             return status;
00415           }
00416         }
00417 
00418       } else {
00419 
00420         /// Natural boundary conditions
00421 
00422         size_t i;
00423         size_t num_points = size;
00424         size_t max_index = num_points - 1;
00425         size_t sys_size = max_index - 1;
00426 
00427         c[0] = 0.0;
00428         c[max_index] = 0.0;
00429 
00430         for (i = 0; i < sys_size; i++)
00431           {
00432             const double h_i   = xa[i + 1] - xa[i];
00433             const double h_ip1 = xa[i + 2] - xa[i + 1];
00434             const double ydiff_i   = ya[i + 1] - ya[i];
00435             const double ydiff_ip1 = ya[i + 2] - ya[i + 1];
00436             const double g_i = (h_i != 0.0) ? 1.0 / h_i : 0.0;
00437             const double g_ip1 = (h_ip1 != 0.0) ? 1.0 / h_ip1 : 0.0;
00438             offdiag[i] = h_ip1;
00439             diag[i] = 2.0 * (h_ip1 + h_i);
00440             g[i] = 3.0 * (ydiff_ip1 * g_ip1 -  ydiff_i * g_i);
00441           }
00442 
00443         if (sys_size == 1)
00444           {
00445             c[1] =  g[0] /  diag[0];
00446             return gsl_success;
00447           }
00448         else
00449           {
00450             gsl_vector_view g_vec = gsl_vector_view_array( g, sys_size);
00451             gsl_vector_view diag_vec = gsl_vector_view_array( diag, sys_size);
00452             gsl_vector_view offdiag_vec = 
00453               gsl_vector_view_array( offdiag, sys_size - 1);
00454             gsl_vector_view solution_vec = 
00455               gsl_vector_view_array (( c) + 1, sys_size);
00456             
00457             int status = gsl_linalg_solve_symm_tridiag(&diag_vec.vector, 
00458                                                        &offdiag_vec.vector, 
00459                                                        &g_vec.vector, 
00460                                                        &solution_vec.vector);
00461             return status;
00462           }
00463       }
00464     }
00465       
00466     
00467     /// Free allocated memory
00468     virtual int free() {
00469       std::free(c);
00470       std::free(g);
00471       std::free(diag);
00472       std::free(offdiag);
00473 
00474       return 0;
00475     }
00476 
00477     /// Give the value of the function \f$ y(x=x_0) \f$ .
00478     virtual int interp(const vec_t &x_array, const vec_t &y_array, 
00479                        size_t size, double x, double &y) {
00480       
00481       double x_lo, x_hi;
00482       double dx;
00483       size_t index;
00484       
00485       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00486   
00487       /* evaluate */
00488       x_hi = x_array[index + 1];
00489       x_lo = x_array[index];
00490       dx = x_hi - x_lo;
00491       if (dx > 0.0)
00492         {
00493           const double y_lo = y_array[index];
00494           const double y_hi = y_array[index + 1];
00495           const double dy = y_hi - y_lo;
00496           double delx = x - x_lo;
00497           double b_i, c_i, d_i; 
00498           coeff_calc( c, dy, dx, index,  &b_i, &c_i, &d_i);
00499           y = y_lo + delx * (b_i + delx * (c_i + delx * d_i));
00500           return gsl_success;
00501         }
00502       else
00503         {
00504           y = 0.0;
00505           return GSL_EINVAL;
00506         }
00507     }
00508 
00509 
00510     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
00511     virtual int deriv(const vec_t &x_array, const vec_t &y_array, 
00512                       size_t size, double x, double &dydx) {
00513 
00514       double x_lo, x_hi;
00515       double dx;
00516       size_t index;
00517   
00518       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00519   
00520       /* evaluate */
00521       x_hi = x_array[index + 1];
00522       x_lo = x_array[index];
00523       dx = x_hi - x_lo;
00524       if (dx > 0.0)
00525         {
00526           const double y_lo = y_array[index];
00527           const double y_hi = y_array[index + 1];
00528           const double dy = y_hi - y_lo;
00529           double delx = x - x_lo;
00530           double b_i, c_i, d_i; 
00531           coeff_calc( c, dy, dx, index,  &b_i, &c_i, &d_i);
00532           dydx = b_i + delx * (2.0 * c_i + 3.0 * d_i * delx);
00533           return gsl_success;
00534         }
00535       else
00536         {
00537           dydx = 0.0;
00538           return GSL_FAILURE;
00539         }
00540     }
00541 
00542     /** \brief Give the value of the second derivative  
00543         \f$ y^{\prime \prime}(x=x_0) \f$ .
00544     */
00545     virtual int deriv2(const vec_t &x_array, const vec_t &y_array, 
00546                        size_t size, double x, double &d2ydx2) {
00547 
00548       double x_lo, x_hi;
00549       double dx;
00550       size_t index;
00551   
00552       index=this->sv.bsearch_inc(x,x_array,0,size-1);
00553   
00554       /* evaluate */
00555       x_hi = x_array[index + 1];
00556       x_lo = x_array[index];
00557       dx = x_hi - x_lo;
00558       if (dx > 0.0)
00559         {
00560           const double y_lo = y_array[index];
00561           const double y_hi = y_array[index + 1];
00562           const double dy = y_hi - y_lo;
00563           double delx = x - x_lo;
00564           double b_i, c_i, d_i;
00565           coeff_calc( c, dy, dx, index,  &b_i, &c_i, &d_i);
00566           d2ydx2 = 2.0 * c_i + 6.0 * d_i * delx;
00567           return gsl_success;
00568         }
00569       else
00570         {
00571           d2ydx2 = 0.0;
00572           return GSL_FAILURE;
00573         }
00574     }
00575 
00576 #ifndef DOXYGEN_INTERNAL
00577 
00578   protected:
00579 
00580     /// An internal function to assist in computing the integral
00581     double integ_eval(double ai, double bi, double ci, double di, double xi, 
00582                       double a, double b) {
00583       const double r1 = a - xi;
00584       const double r2 = b - xi;
00585       const double r12 = r1 + r2;
00586       const double bterm = 0.5 * bi * r12;
00587       const double cterm = (1.0 / 3.0) * ci * (r1 * r1 + r2 * r2 + r1 * r2);
00588       const double dterm = 0.25 * di * r12 * (r1 * r1 + r2 * r2);
00589       
00590       return (b - a) * (ai + bterm + cterm + dterm);
00591     }
00592 
00593 #endif
00594     
00595   public:
00596     
00597     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
00598     virtual int integ(const vec_t &x_array, const vec_t &y_array, 
00599                       size_t size, double a, double b, 
00600                       double &result) {
00601 
00602       size_t i, index_a, index_b;
00603   
00604       index_a=this->sv.bsearch_inc(a,x_array,0,size-1);
00605       index_b=this->sv.bsearch_inc(b,x_array,0,size-1);
00606 
00607       result = 0.0;
00608   
00609       /* interior intervals */
00610       for(i=index_a; i<=index_b; i++) {
00611         const double x_hi = x_array[i + 1];
00612         const double x_lo = x_array[i];
00613         const double y_lo = y_array[i];
00614         const double y_hi = y_array[i + 1];
00615         const double dx = x_hi - x_lo;
00616         const double dy = y_hi - y_lo;
00617         if(dx != 0.0) {
00618           double b_i, c_i, d_i; 
00619           coeff_calc( c, dy, dx, i,  &b_i, &c_i, &d_i);
00620       
00621           if (i == index_a || i == index_b)
00622             {
00623               double x1 = (i == index_a) ? a : x_lo;
00624               double x2 = (i == index_b) ? b : x_hi;
00625               result += integ_eval(y_lo, b_i, c_i, d_i, x_lo, x1, x2);
00626             }
00627           else
00628             {
00629               result += dx * (y_lo + dx*(0.5*b_i + 
00630                                          dx*(c_i/3.0 + 0.25*d_i*dx)));
00631             }
00632         }
00633         else {
00634           result = 0.0;
00635           return GSL_FAILURE;
00636         }
00637       }
00638   
00639       return gsl_success;
00640     }
00641 
00642   };
00643   
00644   /** 
00645       \brief Cubic spline interpolation with periodic 
00646       boundary conditions (GSL)
00647 
00648       This is convenient to allow interpolation objects to
00649       be supplied as template parameters
00650   */
00651   template<class vec_t> class cspline_peri_interp : 
00652   public cspline_interp<vec_t> {
00653 
00654   public:
00655     cspline_peri_interp() : cspline_interp<vec_t>(true) {
00656     }
00657 
00658   };
00659 
00660   /** 
00661       \brief Interpolation class
00662       
00663       This interpolation class is intended to be sufficiently general
00664       to handle any vector type. Interpolation of \ref ovector like
00665       objects is performed with the default template parameters, and
00666       \ref array_interp is provided for simple interpolation on
00667       C-style arrays.
00668 
00669       The class automatically handles decreasing arrays by converting
00670       from an object of type \c vec_t to an object ov type \c rvec_t.
00671 
00672       While \c vec_t may be any vector type which allows indexing via
00673       \c [], \c rvec_t must be a vector type which allows indexing
00674       and has a constructor with one of the two forms
00675       \verbatim
00676       rvec_t::rvec_t(vec_t &v);
00677       rvec_t::rvec_t(vec_t v);
00678       \endverbatim
00679       so that o2scl_interp can automatically "reverse" a vector if 
00680       necessary.
00681       
00682       It is important that different instances of o2scl_interp_vec
00683       and o2scl_interp not be given the same interpolation objects,
00684       as they will clash.
00685       
00686   */
00687   template<class vec_t=ovector_view, class rvec_t=ovector_const_reverse> 
00688     class o2scl_interp {
00689       
00690 #ifndef DOXYGEN_INTERNAL
00691       
00692     protected:
00693       
00694     /// Pointer to base interpolation object
00695     base_interp<vec_t> *itp;
00696       
00697     /// Pointer to base interpolation object for reversed vectors
00698     base_interp<rvec_t> *ritp;
00699       
00700 #endif
00701       
00702     public:
00703       
00704     /// Create with base interpolation objects \c it and \c rit
00705     o2scl_interp(base_interp<vec_t> &it, base_interp<rvec_t> &rit) {
00706       itp=&it;
00707       ritp=&rit;
00708     }
00709 
00710     /** \brief Create with base interpolation object \c it and use
00711         \ref def_ritp for reverse interpolation if necessary
00712     */
00713     o2scl_interp(base_interp<vec_t> &it) {
00714       itp=&it;
00715       ritp=&def_ritp;
00716     }
00717     
00718     /** \brief Create an interpolator using \ref def_itp and 
00719         \ref def_ritp
00720     */
00721     o2scl_interp() {
00722       itp=&def_itp;
00723       ritp=&def_ritp;
00724     }
00725 
00726     virtual ~o2scl_interp() {}
00727 
00728     /// Give the value of the function \f$ y(x=x_0) \f$ .
00729     virtual double interp(const double x0, size_t n, const vec_t &x, 
00730                           const vec_t &y) 
00731     {
00732         
00733       double ret=0.0;
00734       
00735       err_hnd->reset();
00736       
00737       if (n<itp->min_size) {
00738         set_err("Vector size too small in interp().",gsl_edom);
00739       }
00740       
00741       if (x[0]>x[n-1]) {
00742 
00743         // Decreasing case
00744 
00745         rvec_t ax(x);
00746         rvec_t ay(y);
00747 
00748         int r1=ritp->allocate(n);
00749         if (r1==0) {
00750           int r2=ritp->init(ax,ay,n);
00751           if (r2==0) {
00752             ritp->interp(ax,ay,n,x0,ret);
00753             ritp->free();
00754           }
00755         }
00756 
00757       } else {
00758 
00759         // Increasing case
00760           
00761         int r1=itp->allocate(n);
00762         if (r1==0) {
00763           int r2=itp->init(x,y,n);
00764           if (r2==0) {
00765             itp->interp(x,y,n,x0,ret);
00766             itp->free();
00767           }
00768         }
00769           
00770       }
00771       
00772       return ret;
00773     }                 
00774     
00775     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
00776     virtual double deriv(const double x0, size_t n, const vec_t &x, 
00777                          const vec_t &y) 
00778     {
00779       
00780       double ret=0.0;
00781       
00782       err_hnd->reset();
00783       
00784       if (n<itp->min_size) {
00785         set_err("Vector size too small in interp().",gsl_edom);
00786       }
00787       
00788       if (x[0]>x[n-1]) {
00789           
00790         rvec_t ax(x);
00791         rvec_t ay(y);
00792 
00793         int r1=ritp->allocate(n);
00794         if (r1==0) {
00795           int r2=ritp->init(ax,ay,n);
00796           if (r2==0) {
00797             ritp->deriv(ax,ay,n,x0,ret);
00798             ritp->free();
00799           }
00800         }
00801           
00802       } else {
00803 
00804         // Increasing case
00805           
00806         int r1=itp->allocate(n);
00807         if (r1==0) {
00808           int r2=itp->init(x,y,n);
00809           if (r2==0) {
00810             itp->deriv(x,y,n,x0,ret);
00811             itp->free();
00812           }
00813         }
00814       }
00815       
00816       return ret;
00817     }                 
00818     
00819     /** \brief Give the value of the second derivative  
00820         \f$ y^{\prime \prime}(x=x_0) \f$ .
00821     */
00822     virtual double deriv2(const double x0, size_t n, const vec_t &x, 
00823                           const vec_t &y) 
00824     {
00825       
00826       double ret=0.0;
00827       
00828       err_hnd->reset();
00829       
00830       if (n<itp->min_size) {
00831 
00832         set_err("Vector size too small in deriv2().",gsl_edom);
00833 
00834       }
00835       
00836       if (x[0]>x[n-1]) {
00837           
00838         // Decreasing case
00839             
00840         rvec_t ax(x);
00841         rvec_t ay(y);
00842 
00843         int r1=ritp->allocate(n);
00844         if (r1==0) {
00845           int r2=ritp->init(ax,ay,n);
00846           if (r2==0) {
00847             ritp->deriv2(ax,ay,n,x0,ret);
00848             ritp->free();
00849           }
00850         }
00851 
00852       } else {
00853 
00854         // Increasing case
00855 
00856         int r1=itp->allocate(n);
00857         if (r1==0) {
00858           int r2=itp->init(x,y,n);
00859           if (r2==0) {
00860             itp->deriv2(x,y,n,x0,ret);
00861             itp->free();
00862           }
00863         }
00864 
00865       }
00866       
00867       return ret;
00868     }                 
00869     
00870     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
00871     virtual double integ(const double x1, const double x2, size_t n, 
00872                          const vec_t &x, const vec_t &y) 
00873     {
00874       
00875       double ret=0.0;
00876       
00877       err_hnd->reset();
00878       
00879       if (n<itp->min_size) {
00880 
00881         set_err("Vector size too small in integ().",gsl_edom);
00882 
00883       }
00884       
00885       if (x[0]>x[n-1]) {
00886           
00887         // Decreasing case
00888             
00889         rvec_t ax(x);
00890         rvec_t ay(y);
00891 
00892         int r1=ritp->allocate(n);
00893         if (r1==0) {
00894           int r2=ritp->init(ax,ay,n);
00895           if (r2==0) {
00896             ritp->integ(ax,ay,n,x1,x2,ret);
00897             ritp->free();
00898           }
00899         }
00900 
00901       } else {
00902 
00903         // Increasing case
00904             
00905         int r1=itp->allocate(n);
00906         if (r1==0) {
00907           int r2=itp->init(x,y,n);
00908           if (r2==0) {
00909             itp->integ(x,y,n,x1,x2,ret);
00910             itp->free();
00911           }
00912         }
00913 
00914       }
00915       
00916       return ret;
00917     }                 
00918     
00919     /// Set base interpolation object
00920     int set_type(base_interp<vec_t> &it, base_interp<rvec_t> &rit) {
00921       itp=&it;
00922       ritp=&rit;
00923       return 0;
00924     }
00925       
00926     /** \brief Default base interpolation object (cubic spline with natural 
00927         boundary conditions)
00928     */
00929     cspline_interp<vec_t> def_itp;
00930 
00931     /** \brief Default base interpolation object for reversed vectors
00932         (cubic spline with natural boundary conditions)
00933     */
00934     cspline_interp<rvec_t> def_ritp;
00935       
00936   };
00937   
00938   /** 
00939       \brief Interpolation class for pre-specified vector
00940       
00941       This interpolation class is intended to be sufficiently general
00942       to handle any vector type. Interpolation of \ref ovector like
00943       objects is performed with the default template parameters, and
00944       \ref array_interp_vec is provided for simple interpolation on
00945       C-style arrays.
00946 
00947       The class automatically handles decreasing arrays by copying the
00948       old array to a reversed version. For several interpolations on
00949       the same data, copying the initial array can be faster than
00950       overloading operator[].
00951       
00952       It is important that different instances of o2scl_interp_vec
00953       and o2scl_interp not be given the same interpolation objects,
00954       as they will clash.
00955 
00956       \todo Need to fix constructor to behave properly if init() fails.
00957       It should free the memory and set \c ln to zero.
00958   */
00959   template<class vec_t=ovector_view, class alloc_vec_t=ovector, 
00960     class alloc_t=ovector_alloc> class o2scl_interp_vec {
00961       
00962 #ifndef DOXYGEN_INTERNAL
00963       
00964     protected:
00965       
00966     /// Memory allocator for objects of type \c alloc_vec_t
00967     alloc_t ao;
00968 
00969     /// Pointer to base interpolation object
00970     base_interp<vec_t> *itp;
00971       
00972     /// True if the original array was increasing
00973     bool inc;
00974 
00975     /// Pointer to the user-specified x vector
00976     const vec_t *lx;
00977 
00978     /// Pointer to the user-specified x vector
00979     const vec_t *ly;
00980 
00981     /// Reversed version of x
00982     alloc_vec_t lrx;
00983 
00984     /// Reversed version of y
00985     alloc_vec_t lry;
00986 
00987     /// Size of user-specified vectors
00988     size_t ln;
00989 
00990 #endif
00991       
00992     public:
00993     
00994     /// Create with base interpolation object \c it
00995     o2scl_interp_vec(base_interp<vec_t> &it, 
00996                    size_t n, const vec_t &x, const vec_t &y) {
00997       itp=&it;
00998       ln=0;
00999         
01000       if (n<itp->min_size) {
01001         set_err("Vector size too small in o2scl_interp_vec().",gsl_edom);
01002       } else {
01003         size_t iend=n-1;
01004           
01005         if (x[0]>x[iend]) {
01006             
01007           ao.allocate(lrx,n);
01008           ao.allocate(lry,n);
01009 
01010           for(size_t i=0;i<n;i++) {
01011             lrx[iend-i]=x[i];
01012             lry[iend-i]=y[i];
01013           }
01014             
01015           int r1=itp->allocate(n);
01016           if (r1==0) {
01017             itp->init(lrx,lry,n);
01018             ln=n;
01019             inc=false;
01020           }
01021             
01022         } else {
01023             
01024           int r1=itp->allocate(n);
01025           if (r1==0) {
01026             itp->init(x,y,n);
01027             ln=n;
01028             lx=&x;
01029             ly=&y;
01030             inc=true;
01031           }
01032             
01033         }
01034 
01035       } 
01036         
01037     }
01038       
01039     virtual ~o2scl_interp_vec() {
01040       if (ln>0) {
01041         itp->free();
01042       }
01043       if (inc==false) {
01044         ao.free(lrx);
01045         ao.free(lry);
01046       }
01047     }
01048     
01049     /// Give the value of the function \f$ y(x=x_0) \f$ .
01050     virtual double interp(const double x0) {
01051       double ret=0.0;
01052       if (ln>0) {
01053         if (inc) itp->interp(*lx,*ly,ln,x0,ret);
01054         else itp->interp(lrx,lry,ln,x0,ret);
01055       }
01056       return ret;
01057     }                 
01058     
01059     /// Give the value of the derivative \f$ y^{\prime}(x=x_0) \f$ .
01060     virtual double deriv(const double x0) {
01061       double ret=0.0;
01062       if (ln>0) {
01063         if (inc) itp->deriv(*lx,*ly,ln,x0,ret);
01064         else itp->deriv(lrx,lry,ln,x0,ret);
01065       }
01066       return ret;
01067     }                 
01068     
01069     /** \brief Give the value of the second derivative  
01070         \f$ y^{\prime \prime}(x=x_0) \f$ .
01071     */
01072     virtual double deriv2(const double x0) {
01073       double ret=0.0;
01074       if (ln>0) {
01075         if (inc) itp->deriv2(*lx,*ly,ln,x0,ret);
01076         else itp->deriv2(lrx,lry,ln,x0,ret);
01077       }
01078       return ret;
01079     }                 
01080     
01081     /// Give the value of the integral \f$ \int_a^{b}y(x)~dx \f$ .
01082     virtual double integ(const double x1, const double x2) {
01083       double ret=0.0;
01084       if (ln>0) {
01085         if (inc) itp->integ(*lx,*ly,ln,x1,x2,ret);
01086         else itp->integ(lrx,lry,ln,x1,x2,ret);
01087       }
01088       return ret;
01089     }                 
01090     
01091     /// Set base interpolation object
01092     int set_type(base_interp<vec_t> &it) {
01093       itp=&it;
01094       return 0;
01095     }
01096       
01097     /** \brief Default base interpolation object (cubic spline with natural 
01098         boundary conditions)
01099     */
01100     cspline_interp<vec_t> def_itp;
01101       
01102   };
01103   
01104   /** \brief A specialization of o2scl_interp for C-style double arrays
01105    */
01106   template<size_t n> class array_interp :
01107   public o2scl_interp<double[n],array_const_reverse<n> > {
01108     
01109   public:
01110     
01111     /// Create with base interpolation objects \c it and \c rit
01112     array_interp(base_interp<double[n]> &it, 
01113                  base_interp<array_const_reverse<n> > &rit)
01114       : o2scl_interp<double[n],array_const_reverse<n> >(it,rit) {
01115     }
01116     
01117       /** \brief Create with base interpolation object \c it and use
01118           the default base interpolation object for reversed arrays
01119       */
01120       array_interp(base_interp<double[n]> &it) : 
01121         o2scl_interp<double[n],array_const_reverse<n> >(it) {
01122       }
01123       
01124         /** \brief Create an interpolator using the default base
01125             interpolation objects
01126         */
01127         array_interp() : 
01128           o2scl_interp<double[n],array_const_reverse<n> >() {
01129         }
01130 
01131   };
01132 
01133   /** 
01134       \brief A specialization of o2scl_interp_vec for C-style arrays
01135   */
01136   template<class arr_t> class array_interp_vec : 
01137   public o2scl_interp_vec<arr_t,arr_t,array_alloc<arr_t> > {
01138     
01139   public:
01140     
01141     /// Create with base interpolation object \c it
01142     array_interp_vec(base_interp<arr_t> &it, 
01143                      size_t nv, const arr_t &x, const arr_t &y) :
01144       o2scl_interp_vec<arr_t,arr_t,array_alloc<arr_t> >(it,nv,x,y) {
01145     }
01146     
01147   };
01148   
01149 #ifndef DOXYGENP
01150 }
01151 #endif
01152 
01153 #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