00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, 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_EQI_DERIV_H 00024 #define O2SCL_EQI_DERIV_H 00025 00026 #include <cmath> 00027 #include <o2scl/deriv.h> 00028 #include <o2scl/collection.h> 00029 00030 #ifndef DOXYGENP 00031 namespace o2scl { 00032 #endif 00033 00034 /** 00035 \brief Derivatives for equally-spaced abscissas 00036 00037 This is an implementation of the formulas for equally-spaced 00038 abscissas as indicated below. The level of approximation is 00039 specified in set_npoints(). The value of \f$ p \times h \f$ 00040 can be specified in \c xoff (default is zero). 00041 00042 \note Uncertainties are not computed and the code 00043 for second and third derivatives is unfinished. 00044 00045 \note The derivatives given, for example, from the 00046 five-point formula can sometimes be more accurate 00047 than computing the derivative from the interpolation class. 00048 This is especially true near the boundaries of the interpolated 00049 region. 00050 00051 \future Finish the second and third derivative formulas. 00052 00053 Two-point formula (note that this is independent of p). 00054 \f[ 00055 f^{\prime}(x_0+p h)=\frac{1}{h}\left[ 00056 f_{1}-f_{0} \right] 00057 \f] 00058 Three-point formula from Abramowitz and Stegun 00059 \f[ 00060 f^{\prime}(x_0+p h)=\frac{1}{h}\left[ 00061 \frac{2p-1}{2}f_{-1}-2 p f_{0}+\frac{2p+1}{2}f_{1}\right] 00062 \f] 00063 Four-point formula from Abramowitz and Stegun 00064 \f[ 00065 f^{\prime}(x_0+p h)=\frac{1}{h}\left[ 00066 -\frac{3 p^2-6 p+2}{6}f_{-1} 00067 +\frac{3 p^2-4 p -1}{2}f_{0} 00068 -\frac{3 p^2-2 p-2}{2}f_{1} 00069 +\frac{3 p^2-1}{6}f_{2} 00070 \right] 00071 \f] 00072 Five-point formula from Abramowitz and Stegun 00073 \f{eqnarray*} 00074 f^{\prime}(x_0+p h)&=&\frac{1}{h}\left[ 00075 \frac{2 p^3-3 p^2-p+1}{12}f_{-2} 00076 -\frac{4 p^3-3p^2-8p+4}{6}f_{-1} 00077 \right. \\ && \left. 00078 +\frac{2p^3-5p}{2}f_{0} 00079 -\frac{4p^3+3p^2-8p-4}{6}f_{1} 00080 \right. \\ && \left. 00081 +\frac{2p^3+3p^2-p-1}{12}f_{2} 00082 \right] 00083 \f} 00084 00085 The relations above can be confined to give formulas 00086 for second derivative formulas: 00087 Three-point formula 00088 \f[ 00089 f^{\prime}(x_0+p h)=\frac{1}{h^2} 00090 \left[f_{-1}-2 f_0+f_1\right] 00091 \f] 00092 Four-point formula: 00093 \f[ 00094 f^{\prime}(x_0+p h)=\frac{1}{2 h^2} 00095 \left[\left(1-2p\right)f_{-1}-\left(1-6p\right)f_0 00096 -\left(1+6p\right)f_1+\left(1+2p\right)f_2\right] 00097 \f] 00098 Five-point formula: 00099 \f[ 00100 f^{\prime}(x_0+p h)=\frac{1}{4 h^2} 00101 \left[\left(1-2p\right)^2f_{-2} 00102 +\left(8p-16 p^2\right)f_{-1} 00103 -\left(2-24 p^2\right)f_0 00104 -\left(8p+16p^2\right)f_1 00105 +\left(1+2p\right)^2 f_2\right] 00106 \f] 00107 Six-point formula: 00108 \f{eqnarray*} 00109 f^{\prime}(x_0+p h)&=&\frac{1}{12 h^2}\left[ 00110 \left(2-10p+15 p^2-6p^3\right)f_{-2} 00111 +\left(3+14p-57p^2+30p^3\right)f_{-1} 00112 \right. \\ && \left. 00113 +\left(-8+20p+78 p^2-60p^3\right)f_0 00114 +\left(-2-44p-42p^2+60p^3\right)f_1 00115 \right. \\ && \left. 00116 +\left(6+22p+3p^2-30p^3\right)f_2 00117 +\left(-1-2p+3p^2+6p^3\right)f_3 00118 \right] 00119 \f} 00120 Seven-point formula: 00121 \f{eqnarray*} 00122 f^{\prime}(x_0+p h)&=&\frac{1}{36 h^2}\left[ 00123 \left(4-24p+48p^2-36p^3+9p^4\right)f_{-3} 00124 +\left(12+12p-162p^2+180p^3-54p^4\right)f_{-2} 00125 \right. \\ && \left. 00126 +\left(-15+120p+162p^2-360p^3+135p^4\right)f_{-1} 00127 -4\left(8+48p-3p^2-90p^3+45p^4\right)f_0 00128 \right. \\ && \left. 00129 +3\left(14+32p-36p^2-60p^3+45p^4\right)f_1 00130 +\left(-12-12p+54p^2+36p^3-54p^4\right)f_2 00131 \right. \\ && \left. 00132 +\left(1-6p^2+9p^4\right)f_3 00133 \right] 00134 \f} 00135 */ 00136 template<class param_t, class func_t, class vec_t=ovector_base> 00137 class eqi_deriv : public deriv<param_t, func_t> { 00138 public: 00139 00140 eqi_deriv() { 00141 h=1.0e-4; 00142 xoff=0.0; 00143 cap=&eqi_deriv::calc_vector5; 00144 cp=&eqi_deriv::calcp5; 00145 c2p=&eqi_deriv::calc2p5; 00146 } 00147 00148 /// Stepsize (Default \f$ 10^{-4} \f$ ). 00149 double h; 00150 00151 /// Offset (default 0.0) 00152 double xoff; 00153 00154 /** \brief Set the number of points to use for first derivatives 00155 (default 5) 00156 00157 Acceptable values are 2-5 (see above). 00158 */ 00159 int set_npoints(int npoints) { 00160 if (npoints==2) { 00161 cap=&eqi_deriv::calc_vector3; 00162 cp=&eqi_deriv::calcp2; 00163 } else if (npoints==3) { 00164 cap=&eqi_deriv::calc_vector3; 00165 cp=&eqi_deriv::calcp3; 00166 } else if (npoints==4) { 00167 cap=&eqi_deriv::calc_vector4; 00168 cp=&eqi_deriv::calcp4; 00169 } else { 00170 cap=&eqi_deriv::calc_vector5; 00171 cp=&eqi_deriv::calcp5; 00172 } 00173 if (npoints<=1 || npoints>5) { 00174 O2SCL_ERR_RET("Invalid # of points in set_npoints(). Using default", 00175 gsl_einval); 00176 } 00177 return 0; 00178 } 00179 00180 /** \brief Set the number of points to use for second derivatives 00181 (default 5) 00182 00183 Acceptable values are 3-5 (see above). 00184 */ 00185 int set_npoints2(int npoints) { 00186 if (npoints==3) { 00187 c2p=&eqi_deriv::calc2p3; 00188 } else if (npoints==4) { 00189 c2p=&eqi_deriv::calc2p4; 00190 } else { 00191 c2p=&eqi_deriv::calc2p5; 00192 } 00193 if (npoints<=2 || npoints>5) { 00194 O2SCL_ERR_RET("Invalid # of points in set_npoints2(). Using default", 00195 gsl_einval); 00196 } 00197 return 0; 00198 } 00199 00200 /** 00201 \brief Calculate the first derivative of \c func w.r.t. x 00202 */ 00203 virtual int calc_err(double x, param_t &pa, func_t &func, 00204 double &dfdx, double &err) { 00205 double p=xoff/h; 00206 dfdx=(this->*cp)(x,p,pa,func)/h; 00207 err=0.0; 00208 return gsl_success; 00209 } 00210 00211 /** \brief Calculate the second derivative of \c func w.r.t. x 00212 */ 00213 virtual int calc2_err(double x, param_t &pa, func_t &func, 00214 double &dfdx, double &err) { 00215 double p=xoff/h; 00216 dfdx=(this->*c2p)(x,p,pa,func)/h/h; 00217 err=0.0; 00218 return gsl_success;; 00219 } 00220 00221 /** \brief Calculate the third derivative of \c func w.r.t. x 00222 */ 00223 virtual int calc3_err(double x, param_t &pa, func_t &func, 00224 double &dfdx, double &err) { 00225 double p=xoff/h; 00226 dfdx=(this->*c3p)(x,h,p,pa,func)/h; 00227 err=0.0; 00228 return gsl_success;; 00229 } 00230 00231 /** 00232 \brief Calculate the derivative at \c x given an array 00233 00234 This calculates the derivative at \c x given a function 00235 specified in an array \c y of size \c nx with equally spaced 00236 abscissas. The first abscissa should be given as \c x0 00237 and the distance between adjacent abscissas should be 00238 given as \c dx. The value \c x need not be one of the 00239 abscissas (i.e. it can lie in between an interval). The 00240 appropriate offset is calculated automatically. 00241 */ 00242 double calc_vector(double x, double x0, double dx, 00243 size_t nx, const vec_t &y) { 00244 size_t ix=(size_t)((x-x0)/dx); 00245 return (this->*cap)(x,x0,dx,nx,y,ix)/dx; 00246 } 00247 00248 /** 00249 \brief Calculate the second derivative at \c x given an array 00250 00251 This calculates the second derivative at \c x given a function 00252 specified in an array \c y of size \c nx with equally spaced 00253 abscissas. The first abscissa should be given as \c x0 00254 and the distance between adjacent abscissas should be 00255 given as \c dx. The value \c x need not be one of the 00256 abscissas (i.e. it can lie in between an interval). The 00257 appropriate offset is calculated automatically. 00258 */ 00259 double calc2_vector(double x, double x0, double dx, 00260 size_t nx, const vec_t &y) 00261 { 00262 size_t ix=(size_t)((x-x0)/dx); 00263 return (this->*c2ap)(x,x0,dx,nx,y,ix)/dx; 00264 } 00265 00266 /** 00267 \brief Calculate the third derivative at \c x given an array 00268 00269 This calculates the third derivative at \c x given a function 00270 specified in an array \c y of size \c nx with equally spaced 00271 abscissas. The first abscissa should be given as \c x0 and the 00272 distance between adjacent abscissas should be given as \c 00273 dx. The value \c x need not be one of the abscissas (i.e. it 00274 can lie in between an interval). The appropriate offset is 00275 calculated automatically. 00276 */ 00277 double calc3_vector(double x, double x0, double dx, 00278 size_t nx, const vec_t &y) 00279 { 00280 size_t ix=(size_t)((x-x0)/dx); 00281 return (this->*c3ap)(x,x0,dx,nx,y,ix)/dx; 00282 } 00283 00284 /** 00285 \brief Calculate the derivative of an entire array 00286 00287 Right now this uses np=5. 00288 00289 \todo generalize to other values of npoints. 00290 */ 00291 int deriv_vector(size_t nv, double dx, const vec_t &y, 00292 vec_t &dydx) 00293 { 00294 dydx[0]=(-25.0/12.0*y[0]+4.0*y[1]-3.0*y[2]+4.0/3.0*y[3]-0.25*y[4])/dx; 00295 dydx[1]=(-0.25*y[0]-5.0/6.0*y[1]+1.5*y[2]-0.5*y[3]+1.0/12.0*y[4])/dx; 00296 for(size_t i=2;i<nv-2;i++) { 00297 dydx[i]=(1.0/12.0*y[i-2]-2.0/3.0*y[i-1]+2.0/3.0*y[i+1]- 00298 1.0/12.0*y[i+2])/dx; 00299 } 00300 dydx[nv-2]=(-1.0/12.0*y[nv-5]+0.5*y[nv-4]-1.5*y[nv-3]+ 00301 5.0/6.0*y[nv-2]+0.25*y[nv-1])/dx; 00302 dydx[nv-1]=(0.25*y[nv-5]-4.0/3.0*y[nv-4]+3.0*y[nv-3]- 00303 4.0*y[nv-2]+25.0/12.0*y[nv-1])/dx; 00304 return 0; 00305 } 00306 00307 /// Return string denoting type ("eqi_deriv") 00308 virtual const char *type() { return "eqi_deriv"; } 00309 00310 #ifndef DOXYGENP 00311 00312 protected: 00313 00314 /** \brief Calculate the first derivative of \c func w.r.t. x and the 00315 uncertainty 00316 00317 This function doesn't do anything, and isn't required for 00318 this class since it computes higher-order derivatives directly. 00319 */ 00320 virtual int calc_err_int 00321 (double x, typename deriv<param_t,func_t>::dpars &pa, 00322 funct<typename deriv<param_t,func_t>::dpars> &func, 00323 double &dfdx, double &err) { 00324 return gsl_success; 00325 } 00326 00327 /// Two-point first derivative 00328 double calcp2(double x, double p, param_t &pa, 00329 func_t &func) { 00330 return (func(x+h,pa)-func(x,pa)); 00331 } 00332 00333 /// Three-point first derivative 00334 double calcp3(double x, double p, param_t &pa, 00335 func_t &func) { 00336 if (p==0.0) { 00337 return ((-0.5)*func(x-h,pa)+(0.5)*func(x+h,pa)); 00338 } 00339 return ((p-0.5)*func(x-h,pa)-2.0*p*func(x,pa)+(p+0.5)*func(x+h,pa)); 00340 } 00341 00342 /// Four-point first derivative 00343 double calcp4(double x, double p, param_t &pa, 00344 func_t &func) { 00345 double p2=p*p; 00346 return (-(3.0*p2-6.0*p+2.0)/6.0*func(x-h,pa)+ 00347 (1.5*p2-2.0*p-0.5)*func(x,pa)- 00348 (1.5*p2-p-1.0)*func(x+h,pa)+ 00349 (3.0*p2-1.0)/6.0*func(x+2.0*h,pa)); 00350 } 00351 00352 /// Five-point first derivative 00353 double calcp5(double x, double p, param_t &pa, 00354 func_t &func) { 00355 double p2=p*p, p3=p*p*p; 00356 if (p==0.0) { 00357 return ((1.0)/12.0*func(x-2.0*h,pa)- 00358 (4.0)/6.0*func(x-h,pa)- 00359 (-4.0)/6.0*func(x+h,pa)+ 00360 (-1.0)/12.0*func(x+2.0*h,pa)); 00361 } 00362 return ((2.0*p3-3.0*p2-p+1.0)/12.0*func(x-2.0*h,pa)- 00363 (4.0*p3-3.0*p2-8.0*p+4.0)/6.0*func(x-h,pa)+ 00364 (p3-2.5*p)*func(x,pa)- 00365 (4.0*p3+3.0*p2-8.0*p-4.0)/6.0*func(x+h,pa)+ 00366 (2.0*p3+3.0*p2-p-1.0)/12.0*func(x+2.0*h,pa)); 00367 } 00368 00369 00370 /// Three-point first derivative for arrays 00371 double calc_vector3(double x, double x0, double dx, size_t nx, 00372 const vec_t &y, size_t ix) { 00373 double p; 00374 if (ix>0 && ix<nx-1) { 00375 p=x-(x0+ix*dx); 00376 return ((p-0.5)*y[ix-1]-2.0*p*y[ix]+(p+0.5)*y[ix+1]); 00377 } else if (ix==0) { 00378 p=x-(x0+dx); 00379 return ((p-0.5)*y[0]-2.0*p*y[1]+(p+0.5)*y[2]); 00380 } 00381 p=x-(x0+(nx-2)*dx); 00382 return ((p-0.5)*y[nx-3]-2.0*p*y[nx-2]+(p+0.5)*y[nx-1]); 00383 } 00384 00385 /// Four-point first derivative for arrays 00386 double calc_vector4(double x, double x0, double dx, size_t nx, 00387 const vec_t &y, size_t ix) { 00388 double p, p2; 00389 if (ix>0 && ix<nx-2) { 00390 p=x-(x0+ix*dx); 00391 p2=p*p; 00392 return (-(3.0*p2-6.0*p+2.0)/6.0*y[ix-1]+ 00393 (1.5*p2-2.0*p-0.5)*y[ix]- 00394 (1.5*p2-p-1.0)*y[ix+1]+ 00395 (3.0*p2-1.0)/6.0*y[ix+2]); 00396 } else if (ix==0) { 00397 p=x-(x0+dx); 00398 p2=p*p; 00399 return (-(3.0*p2-6.0*p+2.0)/6.0*y[0]+ 00400 (1.5*p2-2.0*p-0.5)*y[1]- 00401 (1.5*p2-p-1.0)*y[2]+ 00402 (3.0*p2-1.0)/6.0*y[3]); 00403 } 00404 p=x-(x0+(nx-3)*dx); 00405 p2=p*p; 00406 return (-(3.0*p2-6.0*p+2.0)/6.0*y[nx-4]+ 00407 (1.5*p2-2.0*p-0.5)*y[nx-3]- 00408 (1.5*p2-p-1.0)*y[nx-2]+ 00409 (3.0*p2-1.0)/6.0*y[nx-1]); 00410 } 00411 00412 /// Five-point first derivative for arrays 00413 double calc_vector5(double x, double x0, 00414 double dx, size_t nx, 00415 const vec_t &y, size_t ix) { 00416 double p, p2, p3; 00417 if (ix>1 && ix<nx-2) { 00418 p=x-(x0+ix*dx); 00419 p2=p*p, p3=p*p*p; 00420 return ((2.0*p3-3.0*p2-p+1.0)/12.0*y[ix-2]- 00421 (4.0*p3-3.0*p2-8.0*p+4.0)/6.0*y[ix-1]+ 00422 (p3-2.5*p)*y[ix]- 00423 (4.0*p3+3.0*p2-8.0*p-4.0)/6.0*y[ix+1]+ 00424 (2.0*p3+3.0*p2-p-1.0)/12.0*y[ix+2]); 00425 } else if (ix<=1) { 00426 p=x-(x0+2*dx); 00427 p2=p*p, p3=p*p*p; 00428 return ((2.0*p3-3.0*p2-p+1.0)/12.0*y[0]- 00429 (4.0*p3-3.0*p2-8.0*p+4.0)/6.0*y[1]+ 00430 (p3-2.5*p)*y[2]- 00431 (4.0*p3+3.0*p2-8.0*p-4.0)/6.0*y[3]+ 00432 (2.0*p3+3.0*p2-p-1.0)/12.0*y[4]); 00433 } 00434 p=x-(x0+(nx-3)*dx); 00435 p2=p*p, p3=p*p*p; 00436 return ((2.0*p3-3.0*p2-p+1.0)/12.0*y[nx-5]- 00437 (4.0*p3-3.0*p2-8.0*p+4.0)/6.0*y[nx-4]+ 00438 (p3-2.5*p)*y[nx-3]- 00439 (4.0*p3+3.0*p2-8.0*p-4.0)/6.0*y[nx-2]+ 00440 (2.0*p3+3.0*p2-p-1.0)/12.0*y[nx-1]); 00441 } 00442 00443 /// Three-point second derivative 00444 double calc2p3(double x, double p, param_t &pa, 00445 func_t &func) { 00446 return (func(x-h,pa)-2*func(x,pa)+func(x+h,pa)); 00447 } 00448 00449 /// Four-point second derivative 00450 double calc2p4(double x, double p, param_t &pa, 00451 func_t &func) { 00452 return ((1.0-2.0*p)*func(x-h,pa)-(1.0-6.0*p)*func(x,pa) 00453 -(1.0-6.0*p)*func(x+h,pa)+(1.0+2.0*p)*func(x+2.0*h,pa))/2.0; 00454 } 00455 00456 /// Five-point second derivative 00457 double calc2p5(double x, double p, param_t &pa, 00458 func_t &func) { 00459 return ((1.0-2.0*p)*(1.0-2.0*p)*func(x-2.0*h,pa) 00460 +(8.0*p-16.0*p*p)*func(x-h,pa) 00461 -(2.0-24.0*p*p)*func(x,pa) 00462 -(8.0*p+16.0*p*p)*func(x+h,pa) 00463 +(1.0+2.0*p)*(1.0+2.0*p)*func(x+2.0*h,pa))/4.0; 00464 } 00465 00466 /// Pointer to the first derivative function 00467 double (eqi_deriv::*cp)(double x, double p, param_t &pa, 00468 func_t &func); 00469 00470 /// Pointer to the first derivative for arrays function 00471 double (eqi_deriv::*cap)(double x, double x0, 00472 double dx, size_t nx, 00473 const vec_t &y, size_t ix); 00474 00475 /// Pointer to the second derivative function 00476 double (eqi_deriv::*c2p)(double x, double p, param_t & pa, 00477 func_t &func); 00478 00479 /// Pointer to the second derivative for arrays function 00480 double (eqi_deriv::*c2ap)(double x, double x0, 00481 double dx, size_t nx, 00482 const vec_t &y, size_t ix); 00483 00484 /// Pointer to the third derivative function 00485 double (eqi_deriv::*c3p)(double x, double h, double p, param_t & pa, 00486 func_t &func); 00487 00488 /// Pointer to the third derivative for arrays function 00489 double (eqi_deriv::*c3ap)(double x, double x0, 00490 double dx, size_t nx, 00491 const vec_t &y, size_t ix); 00492 00493 #endif 00494 00495 }; 00496 00497 #ifndef DOXYGENP 00498 } 00499 #endif 00500 00501 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page