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 /* ode-initval/bsimp.c 00024 * 00025 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004 Gerard Jungman 00026 * 00027 * This program is free software; you can redistribute it and/or modify 00028 * it under the terms of the GNU General Public License as published by 00029 * the Free Software Foundation; either version 3 of the License, or (at 00030 * your option) any later version. 00031 * 00032 * This program is distributed in the hope that it will be useful, but 00033 * WITHOUT ANY WARRANTY; without even the implied warranty of 00034 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00035 * General Public License for more details. 00036 * 00037 * You should have received a copy of the GNU General Public License 00038 * along with this program; if not, write to the Free Software 00039 * Foundation, Inc., 51 Franklin Street, Fifth Floor, 00040 * Boston, MA 02110-1301, USA. 00041 */ 00042 00043 #ifndef O2SCL_GSL_BSIMP_H 00044 #define O2SCL_GSL_BSIMP_H 00045 00046 #include <gsl/gsl_math.h> 00047 00048 #include <o2scl/err_hnd.h> 00049 #include <o2scl/odestep.h> 00050 #include <o2scl/ode_jac_funct.h> 00051 #include <o2scl/vec_arith.h> 00052 #include <o2scl/lu.h> 00053 #include <o2scl/permutation.h> 00054 00055 #ifndef DOXYGENP 00056 namespace o2scl { 00057 #endif 00058 00059 /** 00060 \brief Bulirsch-Stoer implicit ODE stepper (GSL) 00061 00062 Bader-Deuflhard implicit extrapolative stepper (\ref Bader83). 00063 00064 \note The variable <tt>h_next</tt> was defined in the original 00065 GSL version has been removed here, as it was unused by the 00066 stepper routine. 00067 00068 \note At the moment, this class retains the GSL approach to 00069 handling non-integer return values in the user-specified 00070 derivative function. If the user-specified function returns \c 00071 gsl_efailed, then the stepper attempts to decrease the stepsize 00072 to continue. If the user-specified function returns a non-zero 00073 value other than \c gsl_efailed, or if the Jacobian evaluation 00074 returns any non-zero value, then the stepper aborts and returns 00075 the error value without calling the error handler. This behavior 00076 may change in the future. 00077 00078 \future Create an example with a stiff diff eq. which requires 00079 this kind of stepper 00080 00081 \future More detailed documentation about the algorithm 00082 00083 \future Figure out if this should be a child of odestep or 00084 adapt_step. The function step_local() is actually its own ODE 00085 stepper and could be reimplemented as an object of type odestep. 00086 00087 \future I don't like setting yerr to GSL_POSINF, there should 00088 be a better way to force an adaptive stepper which is calling 00089 this stepper to readjust the stepsize. 00090 00091 \future The functions deuf_kchoice() and compute_weights() can 00092 be moved out of this header file 00093 00094 \future Rework internal arrays as uvectors? 00095 */ 00096 template<class param_t, class func_t, class jac_func_t, 00097 class vec_t=ovector_base, class alloc_vec_t=ovector, 00098 class alloc_t=ovector_alloc, class mat_t=omatrix_base, 00099 class alloc_mat_t=omatrix, class mat_alloc_t=omatrix_alloc> 00100 class gsl_bsimp { 00101 00102 #ifndef DOXYGEN_INTERAL 00103 00104 protected: 00105 00106 /// Size of allocated vectors 00107 size_t dim; 00108 00109 /// Memory allocator for objects of type \c alloc_vec_t 00110 alloc_t ao; 00111 00112 /// Memory allocator for objects of type \c alloc_mat_t 00113 mat_alloc_t mo; 00114 00115 /// Function specifying derivatives 00116 func_t *funcp; 00117 00118 /// Jacobian 00119 jac_func_t *jfuncp; 00120 00121 /// User-specified paramter 00122 param_t *pap; 00123 00124 /// Number of entries in the Bader-Deuflhard extrapolation sequence 00125 static const int sequence_count=8; 00126 00127 /// Largest index of the Bader-Deuflhard extrapolation sequence 00128 static const int sequence_max=7; 00129 00130 /// Workspace for extrapolation 00131 gsl_matrix *d; 00132 00133 /// Workspace for linear system matrix 00134 gsl_matrix *a_mat; 00135 00136 /** \brief Workspace for extrapolation 00137 00138 (This state variable was named 'x' in GSL.) 00139 */ 00140 double ex_wk[sequence_max]; 00141 00142 /// \name State info 00143 //@{ 00144 size_t k_current; 00145 size_t k_choice; 00146 double eps; 00147 //@} 00148 00149 /// \name Workspace for extrapolation step 00150 //@{ 00151 double *yp; 00152 double *y_save; 00153 double *yerr_save; 00154 double *y_extrap_save; 00155 alloc_vec_t y_extrap_sequence; 00156 double *extrap_work; 00157 alloc_vec_t dfdt; 00158 alloc_vec_t y_temp; 00159 alloc_vec_t delta_temp; 00160 double *weight; 00161 alloc_mat_t dfdy; 00162 //@} 00163 00164 /// \name Workspace for the basic stepper 00165 //@{ 00166 alloc_vec_t rhs_temp; 00167 double *delta; 00168 //@} 00169 00170 /// Order of last step 00171 size_t order; 00172 00173 /// Compute weights 00174 int compute_weights(const double y[], double w[], size_t n) { 00175 size_t i; 00176 // w[i] is 1 if y[i] is zero and the absolute value of y[i] 00177 // otherwise 00178 for (i = 0; i < n; i++) { 00179 double u = fabs(y[i]); 00180 w[i] = (u > 0.0) ? u : 1.0; 00181 } 00182 return 0; 00183 } 00184 00185 /** 00186 \brief Calculate a choice for the "order" of the method, using the 00187 Deuflhard criteria. 00188 00189 Used in the allocate() function. 00190 */ 00191 size_t deuf_kchoice(double eps2, size_t dimension) { 00192 00193 const double safety_f = 0.25; 00194 const double small_eps = safety_f * eps2; 00195 00196 double a_work[sequence_count]; 00197 double alpha[sequence_max][sequence_max]; 00198 00199 /* Bader-Deuflhard extrapolation sequence */ 00200 static const int bd_sequence[sequence_count] = 00201 { 2, 6, 10, 14, 22, 34, 50, 70 }; 00202 00203 int i, k; 00204 00205 a_work[0] = bd_sequence[0] + 1.0; 00206 00207 for (k = 0; k < sequence_max; k++) { 00208 a_work[k + 1] = a_work[k] + bd_sequence[k + 1]; 00209 } 00210 00211 for (i = 0; i < sequence_max; i++) { 00212 alpha[i][i] = 1.0; 00213 for (k = 0; k < i; k++) { 00214 const double tmp1 = a_work[k + 1] - a_work[i + 1]; 00215 const double tmp2 = (a_work[i + 1] - a_work[0] + 1.0) * 00216 (2 * k + 1); 00217 alpha[k][i] = pow (small_eps, tmp1 / tmp2); 00218 } 00219 } 00220 00221 a_work[0] += dimension; 00222 00223 for (k = 0; k < sequence_max; k++) { 00224 a_work[k + 1] = a_work[k] + bd_sequence[k + 1]; 00225 } 00226 00227 for (k = 0; k < sequence_max - 1; k++) { 00228 if (a_work[k + 2] > a_work[k + 1] * alpha[k][k + 1]) { 00229 break; 00230 } 00231 } 00232 00233 return k; 00234 } 00235 00236 /** 00237 \brief Polynomial extrapolation 00238 00239 Compute the step of index <tt>i_step</tt> using polynomial 00240 extrapolation to evaulate functions by fitting a polynomial 00241 to estimates <tt>(x_i,y_i)</tt> and output the result to 00242 <tt>y_0</tt> and <tt>y_0_err</tt>. 00243 00244 The index <tt>i_step</tt> begins with zero. 00245 00246 */ 00247 int poly_extrap(gsl_matrix *dloc, const double x[], 00248 const unsigned int i_step, const double x_i, 00249 const vec_t &y_i, vec_t &y_0, vec_t &y_0_err, 00250 double work[]) { 00251 size_t j, k; 00252 00253 o2scl::vector_copy(dim,y_i,y_0_err); 00254 o2scl::vector_copy(dim,y_i,y_0); 00255 00256 if (i_step == 0) { 00257 00258 for (j = 0; j < dim; j++) { 00259 gsl_matrix_set (dloc, 0, j, y_i[j]); 00260 } 00261 00262 } else { 00263 00264 o2scl::vector_copy(dim,y_i,work); 00265 00266 for (k = 0; k < i_step; k++) { 00267 double deltaloc = 1.0 / (x[i_step - k - 1] - x_i); 00268 const double f1 = deltaloc * x_i; 00269 const double f2 = deltaloc * x[i_step - k - 1]; 00270 00271 for (j = 0; j < dim; j++) { 00272 const double q_kj = gsl_matrix_get (dloc, k, j); 00273 gsl_matrix_set (dloc, k, j, y_0_err[j]); 00274 deltaloc = work[j] - q_kj; 00275 y_0_err[j] = f1 * deltaloc; 00276 work[j] = f2 * deltaloc; 00277 y_0[j] += y_0_err[j]; 00278 } 00279 } 00280 00281 for (j = 0; j < dim; j++) { 00282 gsl_matrix_set(dloc, i_step, j, y_0_err[j]); 00283 } 00284 } 00285 return 0; 00286 } 00287 00288 /** 00289 \brief Basic implicit Bulirsch-Stoer step 00290 00291 Divide the step <tt>h_total</tt> into <tt>n_step</tt> smaller 00292 steps and do the Bader-Deuflhard semi-implicit iteration. This 00293 function starts at <tt>t0</tt> with function values 00294 <tt>y</tt>, derivatives <tt>yp_loc</tt>, and information from 00295 the Jacobian to compute the final value <tt>y_out</tt>. 00296 */ 00297 int step_local(const double t0, const double h_total, 00298 const unsigned int n_step, const double y[], 00299 const double yp_loc[], const vec_t &dfdt_loc, 00300 const mat_t &dfdy_loc, vec_t &y_out) { 00301 00302 double *const w=weight; 00303 00304 const double h = h_total / n_step; 00305 double t = t0 + h; 00306 00307 double sum; 00308 00309 /* This is the factor sigma referred to in equation 3.4 of the 00310 paper. A relative change in y exceeding sigma indicates a 00311 runaway behavior. According to the authors suitable values 00312 for sigma are >>1. I have chosen a value of 100*dim. BJG 00313 */ 00314 const double max_sum = 100.0 * dim; 00315 00316 int signum, status; 00317 size_t i, j; 00318 size_t n_inter; 00319 00320 /* Calculate the matrix for the linear system. */ 00321 for (i = 0; i < dim; i++) { 00322 for (j = 0; j < dim; j++) { 00323 gsl_matrix_set(a_mat,i,j,-h*dfdy_loc[i][j]); 00324 } 00325 gsl_matrix_set(a_mat,i,i,gsl_matrix_get(a_mat,i,i)+1.0); 00326 } 00327 00328 /* LU decomposition for the linear system. */ 00329 00330 o2scl::permutation p_vec(dim); 00331 omatrix *om=(omatrix *)a_mat; 00332 o2scl_linalg::LU_decomp(dim,*om,p_vec,signum); 00333 00334 /* Compute weighting factors */ 00335 compute_weights(y,w,dim); 00336 00337 /* Initial step. */ 00338 00339 for (i = 0; i < dim; i++) { 00340 y_temp[i]=h*(yp_loc[i]+h*dfdt_loc[i]); 00341 } 00342 00343 o2scl_linalg::LU_solve(dim,*om,p_vec,y_temp,delta_temp); 00344 00345 sum = 0.0; 00346 for (i = 0; i < dim; i++) { 00347 const double di = delta_temp[i]; 00348 delta[i] = di; 00349 y_temp[i] = y[i] + di; 00350 sum += fabs(di) / w[i]; 00351 } 00352 if (sum > max_sum) { 00353 return gsl_efailed; 00354 } 00355 00356 /* Intermediate steps. */ 00357 00358 status=(*funcp)(t,dim,y_temp,y_out,*pap); 00359 if (status) { 00360 return status; 00361 } 00362 00363 for (n_inter = 1; n_inter < n_step; n_inter++) { 00364 00365 for (i = 0; i < dim; i++) { 00366 rhs_temp[i] = h*y_out[i]-delta[i]; 00367 } 00368 00369 o2scl_linalg::LU_solve(dim,*om,p_vec,rhs_temp,delta_temp); 00370 00371 sum = 0.0; 00372 for (i = 0; i < dim; i++) { 00373 delta[i] += 2.0 * delta_temp[i]; 00374 y_temp[i] += delta[i]; 00375 sum += fabs(delta[i]) / w[i]; 00376 } 00377 if (sum > max_sum) { 00378 return gsl_efailed ; 00379 } 00380 00381 t += h; 00382 00383 status=(*funcp)(t,dim,y_temp,y_out,*pap); 00384 if (status) { 00385 return status; 00386 } 00387 } 00388 00389 00390 /* Final step. */ 00391 00392 for (i = 0; i < dim; i++) { 00393 rhs_temp[i]=h*y_out[i]-delta[i]; 00394 } 00395 00396 o2scl_linalg::LU_solve(dim,*om,p_vec,rhs_temp,delta_temp); 00397 00398 sum = 0.0; 00399 for (i = 0; i < dim; i++) { 00400 y_out[i] = y_temp[i]+delta_temp[i]; 00401 sum += fabs(delta_temp[i])/w[i]; 00402 } 00403 00404 if (sum > max_sum) { 00405 return gsl_efailed; 00406 } 00407 00408 return gsl_success; 00409 } 00410 00411 /// Allocate memory for a system of size \c n 00412 int allocate(size_t n) { 00413 00414 if (dim>0) free(); 00415 00416 dim=n; 00417 00418 d=gsl_matrix_alloc(sequence_max,n); 00419 a_mat=gsl_matrix_alloc(n,n); 00420 00421 yp=(double *)malloc(n*sizeof(double)); 00422 00423 // AWS, 12/2/08 - This was added to ensure memory reallocation 00424 // resets the stepper just like reset() does 00425 for(size_t i=0;i<n;i++) yp[i]=0.0; 00426 00427 y_save=(double *)malloc(n*sizeof(double)); 00428 yerr_save=(double *)malloc(n*sizeof(double)); 00429 y_extrap_save=(double *)malloc(n*sizeof(double)); 00430 extrap_work=(double *)malloc(n*sizeof(double)); 00431 weight=(double *)malloc(n*sizeof(double)); 00432 00433 mo.allocate(dfdy,n,n); 00434 ao.allocate(dfdt,n); 00435 ao.allocate(y_extrap_sequence,n); 00436 ao.allocate(y_temp,n); 00437 ao.allocate(rhs_temp,n); 00438 ao.allocate(delta_temp,n); 00439 00440 delta=(double *)malloc(n*sizeof(double)); 00441 00442 // This choice of epsilon is not necessarily optimal, it has 00443 // a "FIXME" comment in the original GSL code 00444 size_t k_choice_loc=deuf_kchoice(GSL_SQRT_DBL_EPSILON,n); 00445 k_choice=k_choice_loc; 00446 k_current=k_choice_loc; 00447 order=2*k_choice_loc; 00448 00449 return 0; 00450 } 00451 00452 /// Free allocated memory 00453 void free() { 00454 if (dim>0) { 00455 std::free(delta); 00456 00457 mo.free(dfdy,dim); 00458 ao.free(rhs_temp); 00459 ao.free(dfdt); 00460 ao.free(y_temp); 00461 ao.free(y_extrap_sequence); 00462 ao.free(delta_temp); 00463 00464 std::free(weight); 00465 std::free(extrap_work); 00466 std::free(y_extrap_save); 00467 std::free(y_save); 00468 std::free(yerr_save); 00469 std::free(yp); 00470 00471 gsl_matrix_free(a_mat); 00472 gsl_matrix_free(d); 00473 dim=0; 00474 } 00475 } 00476 00477 #endif 00478 00479 public: 00480 00481 gsl_bsimp() { 00482 dim=0; 00483 } 00484 00485 virtual ~gsl_bsimp() { 00486 if (dim>0) free(); 00487 } 00488 00489 /// Reset stepper 00490 int reset() { 00491 for(size_t i=0;i<dim;i++) yp[i]=0.0; 00492 return gsl_success; 00493 } 00494 00495 /** 00496 \brief Perform an integration step 00497 00498 Given initial value of the n-dimensional function in \c y and 00499 the derivative in \c dydx (which must generally be computed 00500 beforehand) at the point \c x, take a step of size \c h giving 00501 the result in \c yout, the uncertainty in \c yerr, and the new 00502 derivative in \c dydx_out (at \c x+h) using function \c derivs 00503 to calculate derivatives. Implementations which do not 00504 calculate \c yerr and/or \c dydx_out do not reference these 00505 variables so that a blank \c vec_t can be given. All of the 00506 implementations allow \c yout=y and \c dydx_out=dydx if 00507 necessary. 00508 */ 00509 virtual int step(double x, double h, size_t n, vec_t &y, vec_t &dydx, 00510 vec_t &yout, vec_t &yerr, vec_t &dydx_out, param_t &pa, 00511 func_t &derivs, jac_func_t &jac) { 00512 00513 int ret; 00514 00515 funcp=&derivs; 00516 jfuncp=&jac; 00517 pap=&pa; 00518 00519 if (n!=dim) allocate(n); 00520 00521 /* Bader-Deuflhard extrapolation sequence */ 00522 static const int bd_sequence[sequence_count] = 00523 { 2, 6, 10, 14, 22, 34, 50, 70 }; 00524 00525 double t_local=x; 00526 00527 size_t i, k; 00528 00529 if (h + t_local == t_local) { 00530 // This section is labeled with a "FIXME" comment in the 00531 // original GSL code. I'm not sure why, but an error is 00532 // sensible here. 00533 O2SCL_ERR_RET("Stepsize underflow in gsl_bsimp::step().", 00534 gsl_eundrflw); 00535 } 00536 00537 /* Save inputs */ 00538 o2scl::vector_copy(dim,y,y_extrap_save); 00539 o2scl::vector_copy(dim,y,y_save); 00540 o2scl::vector_copy(dim,yerr,yerr_save); 00541 00542 // Copy derivative 00543 o2scl::vector_copy(dim,dydx,yp); 00544 00545 // Evaluate the Jacobian for the system. */ 00546 ret=jac(t_local,dim,y,dfdy,dfdt,pa); 00547 if (ret != gsl_success) { 00548 return ret; 00549 } 00550 00551 /* Make a series of refined extrapolations, up to the specified 00552 maximum order, which was calculated based on the Deuflhard 00553 criterion in the deuf_kchoice() function (which is called by 00554 allocate() ). 00555 */ 00556 for (k = 0; k <= k_current; k++) { 00557 00558 const unsigned int N = bd_sequence[k]; 00559 const double r = (h / N); 00560 const double x_k = r * r; 00561 00562 // Each step computes a value of y_extrap_sequence, 00563 // using the number of sub-steps dictated by 00564 // the BD sequence 00565 int status=step_local(t_local,h,N,y_extrap_save,yp, 00566 dfdt,dfdy,y_extrap_sequence); 00567 00568 if (status == gsl_efailed) { 00569 /* If the local step fails, set the error to infinity in 00570 order to force a reduction in the step size 00571 */ 00572 for (i = 0; i < dim; i++) { 00573 yerr[i] = GSL_POSINF; 00574 } 00575 break; 00576 } else if (status != gsl_success) { 00577 return status; 00578 } 00579 00580 // Use the information in y_extrap_sequence to compute 00581 // the new value of y and yerr . 00582 ex_wk[k] = x_k; 00583 poly_extrap(d,ex_wk,k,x_k,y_extrap_sequence,y,yerr,extrap_work); 00584 } 00585 00586 /* Evaluate dydt_out[]. */ 00587 00588 ret=derivs(t_local+h,dim,y,dydx_out,pa); 00589 00590 // If we failed, copy the old values back to y and yerr 00591 if (ret != gsl_success) { 00592 o2scl::vector_copy(dim,y_save,y); 00593 o2scl::vector_copy(dim,yerr_save,yerr); 00594 return ret; 00595 } 00596 00597 return gsl_success; 00598 } 00599 00600 }; 00601 00602 #ifndef DOXYGENP 00603 } 00604 #endif 00605 00606 #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