00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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_GSL_DERIV_H 00024 #define O2SCL_GSL_DERIV_H 00025 00026 #include <iostream> 00027 #include <o2scl/deriv.h> 00028 #include <gsl/gsl_deriv.h> 00029 #include <gsl/gsl_errno.h> 00030 00031 #ifndef DOXYGENP 00032 namespace o2scl { 00033 #endif 00034 00035 /** 00036 \brief Numerical differentiation (GSL) 00037 00038 This class computes the numerical derivative of a function. The 00039 stepsize \ref h should be specified before use. If similar 00040 functions are being differentiated in succession, the user may 00041 be able to increase the speed of later derivatives by setting 00042 the new stepsize equal to the optimized stepsize from the 00043 previous differentiation, by setting \ref h to \ref h_opt. 00044 00045 \note Second and third derivatives are computed by naive nested 00046 applications of the formula for the first derivative and the 00047 uncertainty for these will likely be underestimated. 00048 00049 \future Include the forward and backward GSL derivatives 00050 00051 */ 00052 template<class param_t, class func_t> 00053 class gsl_deriv : public deriv<param_t,func_t> { 00054 public: 00055 00056 gsl_deriv() { 00057 h=0.0; 00058 h_opt=0.0; 00059 } 00060 00061 virtual ~gsl_deriv() {} 00062 00063 /** 00064 \brief Initial stepsize 00065 00066 This should be specified before a call to calc() or 00067 calc_err(). If it is zero, then \f$ x 10^{-4} \f$ will used, or 00068 if \c x is zero, then \f$ 10^{-4} \f$ will be used. 00069 */ 00070 double h; 00071 00072 /** 00073 \brief The last value of the optimized stepsize 00074 00075 This is initialized to zero in the constructor and 00076 set by calc_err() to the most recent value of the 00077 optimized stepsize. 00078 */ 00079 double h_opt; 00080 00081 /** \brief Calculate the first derivative of \c func w.r.t. x and 00082 uncertainty 00083 */ 00084 virtual int calc_err(double x, param_t &pa, func_t &func, 00085 double &dfdx, double &err) { 00086 00087 double hh; 00088 if (h==0.0) { 00089 if (x==0.0) hh=1.0e-4; 00090 else hh=1.0e-4*x; 00091 } else { 00092 hh=h; 00093 } 00094 00095 double r_0, round, trunc, error; 00096 00097 central_deriv(x,h,r_0,round,trunc,func,pa); 00098 error = round + trunc; 00099 00100 if (round < trunc && (round > 0 && trunc > 0)) { 00101 double r_opt, round_opt, trunc_opt, error_opt; 00102 00103 /* Compute an optimised stepsize to minimize the total error, 00104 using the scaling of the truncation error (O(h^2)) and 00105 rounding error (O(1/h)). */ 00106 00107 h_opt = h * pow (round / (2.0 * trunc), 1.0 / 3.0); 00108 central_deriv(x,h_opt,r_opt,round_opt,trunc_opt,func,pa); 00109 error_opt = round_opt + trunc_opt; 00110 00111 /* Check that the new error is smaller, and that the new derivative 00112 is consistent with the error bounds of the original estimate. */ 00113 00114 if (error_opt < error && fabs (r_opt - r_0) < 4.0 * error) { 00115 r_0 = r_opt; 00116 error = error_opt; 00117 } 00118 } 00119 00120 dfdx=r_0; 00121 err=error; 00122 00123 return 0; 00124 } 00125 00126 /// Return string denoting type ("gsl_deriv") 00127 virtual const char *type() { return "gsl_deriv"; } 00128 00129 #ifndef DOXYGEN_INTERNAL 00130 00131 protected: 00132 00133 /** \brief Internal version of calc_err() for second 00134 and third derivatives 00135 */ 00136 virtual int calc_err_int 00137 (double x, typename deriv<param_t,func_t>::dpars &pa, 00138 funct<typename deriv<param_t,func_t>::dpars> &func, 00139 double &dfdx, double &err) { 00140 00141 err_hnd->reset(); 00142 00143 double hh; 00144 if (h==0.0) { 00145 if (x==0.0) hh=1.0e-4; 00146 else hh=1.0e-4*x; 00147 } else { 00148 hh=h; 00149 } 00150 00151 double r_0, round, trunc, error; 00152 00153 central_deriv(x,h,r_0,round,trunc,func,pa); 00154 error = round + trunc; 00155 00156 if (round < trunc && (round > 0 && trunc > 0)) { 00157 double r_opt, round_opt, trunc_opt, error_opt; 00158 00159 /* Compute an optimised stepsize to minimize the total error, 00160 using the scaling of the truncation error (O(h^2)) and 00161 rounding error (O(1/h)). 00162 */ 00163 00164 // This value 'hopt' is different from 'h_opt' above so that 00165 // the two derivative functions don't get confused about which 00166 // 'hopt' they are using. 00167 double hopt = h * pow (round / (2.0 * trunc), 1.0 / 3.0); 00168 central_deriv(x,hopt,r_opt,round_opt,trunc_opt,func,pa); 00169 error_opt = round_opt + trunc_opt; 00170 00171 /* Check that the new error is smaller, and that the new derivative 00172 is consistent with the error bounds of the original estimate. 00173 */ 00174 00175 if (error_opt < error && fabs (r_opt - r_0) < 4.0 * error) { 00176 r_0 = r_opt; 00177 error = error_opt; 00178 } 00179 } 00180 00181 dfdx=r_0; 00182 err=error; 00183 00184 return 0; 00185 } 00186 00187 /** 00188 \brief Compute derivative using 5-point rule 00189 00190 Compute the derivative using the 5-point rule (x-h, x-h/2, x, 00191 x+h/2, x+h) and the error using the difference between the 00192 5-point and the 3-point rule (x-h,x,x+h). Note that the 00193 central point is not used for either. 00194 00195 This must be a class template because it is used by 00196 both calc_err() and calc_err_int(). 00197 */ 00198 template<class func2_t, class param2_t> 00199 int central_deriv(double x, double hh, double &result, 00200 double &abserr_round, double &abserr_trunc, 00201 func2_t &func, param2_t &pa) { 00202 00203 double fm1, fp1, fmh, fph; 00204 00205 func(x-hh,fm1,pa); 00206 func(x+hh,fp1,pa); 00207 00208 func(x-hh/2,fmh,pa); 00209 func(x+hh/2,fph,pa); 00210 00211 double r3 = 0.5 * (fp1 - fm1); 00212 double r5 = (4.0 / 3.0) * (fph - fmh) - (1.0 / 3.0) * r3; 00213 00214 double e3 = (fabs (fp1) + fabs (fm1)) * GSL_DBL_EPSILON; 00215 double e5 = 2.0 * (fabs (fph) + fabs (fmh)) * GSL_DBL_EPSILON + e3; 00216 00217 /* The next term is due to finite precision in x+h = O (eps * x) */ 00218 00219 double dy=GSL_MAX(fabs(r3/hh),fabs(r5/hh))*fabs(x/hh)*GSL_DBL_EPSILON; 00220 00221 /* The truncation error in the r5 approximation itself is O(h^4). 00222 However, for safety, we estimate the error from r5-r3, which is 00223 O(h^2). By scaling h we will minimise this estimated error, not 00224 the actual truncation error in r5. 00225 */ 00226 00227 result = r5 / hh; 00228 /* Estimated truncation error O(h^2) */ 00229 abserr_trunc = fabs ((r5 - r3) / hh); 00230 /* Rounding error (cancellations) */ 00231 abserr_round = fabs (e5 / hh) + dy; 00232 00233 if (this->verbose>0) { 00234 std::cout << "gsl_deriv: " << std::endl; 00235 std::cout << "step: " << hh << std::endl; 00236 std::cout << "abcissas: " << x-hh/2 << " " << x-hh << " " 00237 << x+hh/2 << " " << x+hh << std::endl; 00238 std::cout << "ordinates: " << fm1 << " " << fmh << " " << fph << " " 00239 << fp1 << std::endl; 00240 std::cout << "res: " << result << " trc: " << abserr_trunc 00241 << " rnd: " << abserr_round << std::endl; 00242 } 00243 00244 return 0; 00245 } 00246 00247 #endif 00248 00249 }; 00250 00251 #ifndef DOXYGENP 00252 } 00253 #endif 00254 00255 #endif 00256 00257 00258
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